// This routine simply glues together many of the routines that are already
// written in the Poisson solver library
//
// phi( 1:SubNumPhysNodes       ) is a scalar quantity.  
//
// E1 ( 1:NumElems, 1:kmax2d ) is a vector quantity.
// E2 ( 1:NumElems, 1:kmax2d ) is a vector quantity.
//
// See also: ConvertEfieldOntoDGbasis
void ComputeElectricField( const double t, const mesh& Mesh, const dTensorBC5& q,
    dTensor2& E1, dTensor2& E2)
{

    //
    const int       mx = q.getsize(1);   assert_eq(mx,dogParamsCart2.get_mx());
    const int       my = q.getsize(2);   assert_eq(my,dogParamsCart2.get_my());
    const int NumElems = q.getsize(3);
    const int     meqn = q.getsize(4);
    const int     kmax = q.getsize(5);

    const int space_order = dogParams.get_space_order();

    // unstructured parameters:
    const int kmax2d    = E2.getsize(2);
    const int NumBndNodes  = Mesh.get_NumBndNodes();
    const int NumPhysNodes = Mesh.get_NumPhysNodes();

    // Quick error check
    if( !Mesh.get_is_submesh() )
    {
        printf("ERROR: mesh needs to have subfactor set to %d\n", space_order);
        printf("Go to Unstructured mesh and remesh the problem\n");
        exit(-1);
    }
    const int SubFactor    = Mesh.get_SubFactor();

    assert_eq( NumElems, Mesh.get_NumElems() );

    // -- Step 1: Compute rho -- //
    dTensor3 rho(NumElems, 1, kmax2d );
    void ComputeDensity( const mesh& Mesh, const dTensorBC5& q, dTensor3& rho );
    ComputeDensity( Mesh, q, rho );

    // -- Step 2: Figure out how large phi needs to be
    int SubNumPhysNodes = 0;
    int SubNumBndNodes  = 0;
    switch( dogParams.get_space_order() )
    {
        case 1:
            SubNumPhysNodes = NumPhysNodes;
            SubNumBndNodes  = NumBndNodes;
            break;

        case 2:
            SubNumPhysNodes = Mesh.get_SubNumPhysNodes();
            SubNumBndNodes  = Mesh.get_SubNumBndNodes();
            if(SubFactor!=2)
            {
                printf("\n");
                printf(" Error: for space_order = %i, need SubFactor = %i\n",space_order,2);
                printf("      SubFactor = %i\n",SubFactor);
                printf("\n");
                exit(1);
            }
            break;

        case 3:
            SubNumPhysNodes = Mesh.get_SubNumPhysNodes();
            SubNumBndNodes  = Mesh.get_SubNumBndNodes();
            if(SubFactor!=3)
            {
                printf("\n");
                printf(" Error: for space_order = %i, need SubFactor = %i\n",space_order,3);
                printf("      SubFactor = %i\n",SubFactor);
                printf("\n");
                exit(1);
            }
            break;

        default:
            printf("\n");
            printf(" ERROR in RunDogpack_unst.cpp: space_order value not supported.\n");
            printf("       space_order = %i\n",space_order);
            printf("\n");
            exit(1);
    }

    // local storage:
    dTensor1 rhs(SubNumPhysNodes);
    dTensor1 phi(SubNumPhysNodes);

    // Get Cholesky factorization matrix R
    //
    // TODO - this should be saved earlier in the code rather than reading
    // from file every time we with to run a Poisson solve!
    //
    SparseCholesky R(SubNumPhysNodes);
    string outputdir = dogParams.get_outputdir();
    R.init(outputdir);
    R.read(outputdir);

    // Create right-hand side vector
    void Rhs2D_unst(const int space_order,
            const mesh& Mesh, const dTensor3& rhs_dg,
            dTensor1& rhs);
    Rhs2D_unst(space_order, Mesh, rho, rhs);

    // Call Poisson solver  
    void PoissonSolver2D_unst(const int space_order,
            const mesh& Mesh,
            const SparseCholesky& R,
            const dTensor1& rhs,
            dTensor1& phi,
            dTensor2& E1,
            dTensor2& E2);
    PoissonSolver2D_unst(space_order, Mesh, R, rhs, phi, E1, E2);

    // Compare errors with the exact Electric field:
    //
    void L2Project_Unst(
        const double time,
        const dTensor2* vel_vec,
        const int istart, 
        const int iend, 
        const int QuadOrder, 
        const int BasisOrder_qin,
        const int BasisOrder_auxin,
        const int BasisOrder_fout,
        const mesh& Mesh, 
        const dTensor3* qin, 
        const dTensor3* auxin, 
        dTensor3* fout, 
        void (*Func)(const double t, const dTensor2* vel_vec, const dTensor2&,const dTensor2&,
            const dTensor2&,dTensor2&));

    const int sorder = dogParams.get_space_order();
    dTensor3 qtmp   (NumElems, 2, kmax2d );  qtmp.setall(0.);
    dTensor3 auxtmp (NumElems, 0, kmax2d );
    dTensor3 ExactE (NumElems, 2, kmax2d );
    L2Project_Unst( t, NULL, 1, NumElems, 
        sorder, sorder, sorder, sorder, Mesh, 
        &qtmp, &auxtmp, &ExactE, 
        &ExactElectricField );

    // Compute errors on these two:
    //
    double err = 0.;
    for( int n=1; n <= NumElems; n++ )
    for( int k=1; k <= kmax2d;   k++ )
    {
        err += Mesh.get_area_prim(n)*pow( ExactE.get(n,1,k) - E1.get(n,k), 2 );
        err += Mesh.get_area_prim(n)*pow( ExactE.get(n,2,k) - E2.get(n,k), 2 );
    }
    printf("error = %2.15e\n", err );

}
示例#2
0
void ConstructA_CG2(const mesh& Mesh, FullMatrix& A)
{
  const int NumPhysElems = Mesh.get_NumPhysElems();
  const int NumBndNodes  = Mesh.get_SubNumBndNodes();
  const int Asize = Mesh.get_SubNumPhysNodes();

  assert_eq(Asize,A.get_NumRows());
  assert_eq(Asize,A.get_NumCols());
  
  dTensor1 A1(6);
  dTensor1 A2(6);
  dTensor1 A3(6);
  dTensor1 A4(6);
  dTensor1 A5(6);
  dTensor1 A6(6);

  A1.set(1, -oneninth     );
  A1.set(2,  4.0*oneninth );
  A1.set(3, -oneninth     );
  A1.set(4,  4.0*oneninth );
  A1.set(5,  4.0*oneninth );
  A1.set(6, -oneninth     );
  
  A2.set(1, -onethird     );
  A2.set(2,  0.0          );
  A2.set(3,  onethird     );
  A2.set(4, -4.0*onethird );
  A2.set(5,  4.0*onethird );
  A2.set(6,  0.0          );
  
  A3.set(1, -onethird     );
  A3.set(2, -4.0*onethird );
  A3.set(3,  0.0          );
  A3.set(4,  0.0          );
  A3.set(5,  4.0*onethird );
  A3.set(6,  onethird     );
  
  A4.set(1,  4.0          );
  A4.set(2, -4.0          );
  A4.set(3,  0.0          );
  A4.set(4, -4.0          );
  A4.set(5,  4.0          );
  A4.set(6,  0.0          );

  A5.set(1,  2.0          );
  A5.set(2, -4.0          );
  A5.set(3,  2.0          );
  A5.set(4,  0.0          );
  A5.set(5,  0.0          );
  A5.set(6,  0.0          );
  
  A6.set(1,  2.0          );
  A6.set(2,  0.0          );
  A6.set(3,  0.0          );
  A6.set(4, -4.0          );
  A6.set(5,  0.0          );
  A6.set(6,  2.0          );

  dTensor2 spts(3,2);
  spts.set(1,1,  1.0/3.0 );
  spts.set(1,2, -1.0/6.0 );
  
  spts.set(2,1, -1.0/6.0 );
  spts.set(2,2, -1.0/6.0 );
  
  spts.set(3,1, -1.0/6.0 );
  spts.set(3,2,  1.0/3.0 );
  
  dTensor1 wgts(3);
  wgts.set(1, 1.0/6.0 );
  wgts.set(2, 1.0/6.0 );
  wgts.set(3, 1.0/6.0 );
  
  // Loop over all elements in the mesh
  for (int i=1; i<=NumPhysElems; i++)
    {
      // Information for element i
      iTensor1 tt(6);
      for (int k=1; k<=6; k++)
	{  tt.set(k, Mesh.get_node_subs(i,k) );  }
      
      // Evaluate gradients of the Lagrange polynomials on Gauss quadrature points      
      dTensor2 gpx(6,3);
      dTensor2 gpy(6,3);
      
      for (int m=1; m<=3; m++)
	{
	  double  xi = spts.get(m,1);
	  double eta = spts.get(m,2);
	  
	  for (int k=1; k<=6; k++)
	    {
	      double gp_xi  = A2.get(k) + 2.0*A5.get(k)*xi + A4.get(k)*eta;
	      double gp_eta = A3.get(k) + A4.get(k)*xi + 2.0*A6.get(k)*eta;

	      gpx.set(k,m, Mesh.get_jmat(i,1,1)*gp_xi
		         + Mesh.get_jmat(i,1,2)*gp_eta );
	      gpy.set(k,m, Mesh.get_jmat(i,2,1)*gp_xi
		         + Mesh.get_jmat(i,2,2)*gp_eta );
	    }
	}

      // Entries of the stiffness matrix A
      double Area = Mesh.get_area_prim(i);
      for (int j=1; j<=6; j++)
	for (int k=1; k<=6; k++)
	  {
	    double tmp = A.get(tt.get(j),tt.get(k));
	    for (int m=1; m<=3; m++)
	      {
		tmp = tmp + 2.0*Area*wgts.get(m)*(gpx.get(j,m)*gpx.get(k,m)+gpy.get(j,m)*gpy.get(k,m));
	      }
	    A.set(tt.get(j),tt.get(k), tmp );
	  }
    }

  // Replace boundary node equations by Dirichlet boundary condition enforcement
  for (int i=1; i<=NumBndNodes; i++)
    {
      const int j=Mesh.get_sub_bnd_node(i);
      
      for (int k=1; k<=A.get_NumCols(); k++)
	{
	  A.set(j,k, 0.0 );	  
	}
      for (int k=1; k<=A.get_NumRows(); k++)
	{
	  A.set(k,j, 0.0 );
	}
      A.set(j,j, 1.0 );
    }

  // Get sparse structure representation
  A.Sparsify();
  
}