67 #include "../common/mfem-common.hpp" 94 int main(
int argc,
char *argv[])
97 const char *mesh_file =
"../../data/inline-quad.mesh";
101 bool visualization =
true;
104 bool static_cond =
false;
108 args.
AddOption(&mesh_file,
"-m",
"--mesh",
109 "Mesh file to use.");
111 "Finite element order (polynomial degree).");
112 args.
AddOption(&delta_order,
"-do",
"--delta-order",
113 "Order enrichment for DPG test space.");
115 "Epsilon coefficient");
116 args.
AddOption(&ref,
"-ref",
"--num-refinements",
117 "Number of uniform refinements");
118 args.
AddOption(&theta,
"-theta",
"--theta",
119 "Theta parameter for AMR");
120 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case" 121 " 0: manufactured, 1: Erickson-Johnson ");
123 "Vector Coefficient beta");
124 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
125 "--no-static-condensation",
"Enable static condensation.");
126 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
127 "--no-visualization",
128 "Enable or disable GLVis visualization.");
136 if (iprob > 1) { iprob = 1; }
141 mesh_file =
"../../data/inline-quad.mesh";
144 Mesh mesh(mesh_file, 1, 1);
146 MFEM_VERIFY(
dim > 1,
"Dimension = 1 is not supported in this example");
187 int test_order = order+delta_order;
211 trial_fes.
Append(sigma_fes);
212 trial_fes.
Append(hatu_fes);
213 trial_fes.
Append(hatf_fes);
224 a->StoreMatrices(
true);
228 TrialSpace::u_space, TestSpace::v_space);
232 TrialSpace::sigma_space, TestSpace::v_space);
236 TrialSpace::u_space, TestSpace::tau_space);
240 TrialSpace::sigma_space, TestSpace::tau_space);
244 TrialSpace::hatu_space, TestSpace::tau_space);
248 TrialSpace::hatf_space, TestSpace::v_space);
257 TestSpace::v_space, TestSpace::v_space);
260 TestSpace::v_space, TestSpace::v_space);
263 TestSpace::v_space, TestSpace::v_space);
266 TestSpace::tau_space, TestSpace::tau_space);
269 TestSpace::tau_space, TestSpace::tau_space);
288 std::cout <<
"\n Ref |" 293 <<
" Rate |" << endl;
294 std::cout << std::string(64,
'-') << endl;
296 if (static_cond) {
a->EnableStaticCondensation(); }
297 for (
int it = 0; it<=ref; it++)
309 ess_bdr_uhat = 1; ess_bdr_fhat = 0;
322 int n = ess_tdof_list_uhat.
Size();
323 int m = ess_tdof_list_fhat.
Size();
325 for (
int j = 0; j < n; j++)
327 ess_tdof_list[j] = ess_tdof_list_uhat[j]
331 for (
int j = 0; j < m; j++)
333 ess_tdof_list[j+n] = ess_tdof_list_fhat[j]
342 for (
int i = 0; i<trial_fes.
Size(); i++)
344 offsets[i+1] = trial_fes[i]->GetVSize();
345 dofs += trial_fes[i]->GetTrueVSize();
357 a->FormLinearSystem(ess_tdof_list,x,Ah,X,B);
376 a->RecoverFEMSolution(X,x);
384 double L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
386 Vector & residuals =
a->ComputeResidual(x);
387 double residual = residuals.
Norml2();
389 double rate_err = (it) ?
dim*log(err0/L2Error)/log((
double)dof0/dofs) : 0.0;
390 double rate_res = (it) ?
dim*log(res0/residual)/log((
double)dof0/dofs) : 0.0;
396 std::ios oldState(
nullptr);
397 oldState.copyfmt(std::cout);
398 std::cout << std::right << std::setw(5) << it <<
" | " 399 << std::setw(10) << dof0 <<
" | " 400 << std::setprecision(3)
401 << std::setw(10) << std::scientific << err0 <<
" | " 402 << std::setprecision(2)
403 << std::setw(6) << std::fixed << rate_err <<
" | " 404 << std::setprecision(3)
405 << std::setw(10) << std::scientific << res0 <<
" | " 406 << std::setprecision(2)
407 << std::setw(6) << std::fixed << rate_res <<
" | " 409 std::cout.copyfmt(oldState);
413 const char * keys = (it == 0 &&
dim == 2) ?
"jRcm\n" :
nullptr;
417 "Numerical u", 0,0, 500, 500, keys);
419 "Numerical flux", 501,0,500, 500, keys);
428 double max_resid = residuals.
Max();
429 for (
int iel = 0; iel<mesh.
GetNE(); iel++)
431 if (residuals[iel] > theta * max_resid)
433 elements_to_refine.
Append(iel);
438 for (
int i =0; i<trial_fes.
Size(); i++)
440 trial_fes[i]->Update(
false);
472 if (X.
Size() == 3) { z = X[2]; }
480 double denom = exp(-r2) - exp(-r1);
482 double g1 = exp(r2*(x-1.));
483 double g2 = exp(r1*(x-1.));
486 return g * cos(M_PI * y)/denom;
491 double alpha = M_PI * (x + y + z);
503 if (X.
Size() == 3) { z = X[2]; }
513 double denom = exp(-r2) - exp(-r1);
515 double g1 = exp(r2*(x-1.));
517 double g2 = exp(r1*(x-1.));
520 double g_x = g1_x - g2_x;
522 du[0] = g_x * cos(M_PI * y)/denom;
523 du[1] = -M_PI * g * sin(M_PI*y)/denom;
528 double alpha = M_PI * (x + y + z);
530 for (
int i = 0; i<du.
Size(); i++)
532 du[i] = M_PI * cos(
alpha);
544 if (X.
Size() == 3) { z = X[2]; }
553 double denom = exp(-r2) - exp(-r1);
555 double g1 = exp(r2*(x-1.));
557 double g1_xx = r2*g1_x;
558 double g2 = exp(r1*(x-1.));
560 double g2_xx = r1*g2_x;
562 double g_xx = g1_xx - g2_xx;
564 double u = g * cos(M_PI * y)/denom;
565 double u_xx = g_xx * cos(M_PI * y)/denom;
566 double u_yy = -M_PI * M_PI *
u;
572 double alpha = M_PI * (x + y + z);
574 return -M_PI*M_PI *
u * X.
Size();
598 for (
int i = 0; i<hatf.
Size(); i++)
612 for (
int i = 0; i<du.
Size(); i++)
614 s +=
beta[i] * du[i];
624 for (
int i = 0; i < mesh->
GetNE(); i++)
627 double c1 = min(
epsilon/volume, 1.);
628 double c2 = min(1./
epsilon, 1./volume);
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
A matrix coefficient that is constant in space and time.
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Class for domain integration L(v) := (f, v)
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Data type for scaled Jacobi-type smoother of sparse matrix.
A class to handle Vectors in a block fashion.
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
int Dimension() const
Dimension of the reference space used within the elements.
double exact_u(const Vector &X)
void SetSize(int s)
Resize the vector to size s.
void PrintUsage(std::ostream &out) const
Print the usage message.
Vector coefficient that is constant in space and time.
Pointer to an Operator of a specified type.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Data type dense matrix using column-major storage.
double exact_laplacian_u(const Vector &X)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int main(int argc, char *argv[])
bool Good() const
Return true if the command line options were parsed successfully.
(Q div u, div v) for RT elements
void exact_gradu(const Vector &X, Vector &du)
void exact_hatf(const Vector &X, Vector &hatf)
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
int NumRowBlocks() const
Return the number of row blocks.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
void SetMaxIter(int max_it)
double f_exact(const Vector &X)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Array< int > & RowOffsets()
Return the row offsets for block starts.
void VisualizeField(socketstream &sock, const char *vishost, int visport, ParGridFunction &gf, const char *title, int x=0, int y=0, int w=400, int h=400, bool vec=false)
FiniteElementSpace * FESpace()
void setup_test_norm_coeffs(GridFunction &c1_gf, GridFunction &c2_gf)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
void exact_sigma(const Vector &X, Vector &sigma)
A general vector function coefficient.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Mesh * GetMesh() const
Returns the mesh.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
double Max() const
Returns the maximal element of the vector.
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
int GetNE() const
Returns number of elements.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
double Norml2() const
Returns the l2 norm of the vector.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
int Size() const
Return the logical size of the array.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
A general function coefficient.
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
double sigma(const Vector &x)
double u(const Vector &xvec)
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
double GetElementVolume(int i)
double exact_hatu(const Vector &X)
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
void Neg()
(*this) = -(*this)