56 int main(
int argc,
char *argv[])
59 const char *mesh_file =
"../data/beam-hex.mesh";
62 bool static_cond =
false;
64 const char *device_config =
"cpu";
65 bool visualization = 1;
68 args.
AddOption(&mesh_file,
"-m",
"--mesh",
71 "Finite element order (polynomial degree).");
72 args.
AddOption(&prob,
"-p",
"--problem-type",
73 "Choose between 0: grad, 1: curl, 2: div");
74 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
75 "--no-static-condensation",
"Enable static condensation.");
76 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
77 "--no-partial-assembly",
"Enable Partial Assembly.");
78 args.
AddOption(&device_config,
"-d",
"--device",
79 "Device configuration string, see Device::Configure().");
80 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
82 "Enable or disable GLVis visualization.");
95 Device device(device_config);
101 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
110 int ref_levels = (int)floor(
log(50000./mesh->
GetNE())/
log(2.)/
dim);
111 for (
int l = 0; l < ref_levels; l++)
146 cout <<
"Number of Nedelec finite element unknowns: " << test_size << endl;
147 cout <<
"Number of H1 finite element unknowns: " << trial_size << endl;
151 cout <<
"Number of Nedelec finite element unknowns: " << trial_size << endl;
152 cout <<
"Number of Raviart-Thomas finite element unknowns: " << test_size <<
157 cout <<
"Number of Raviart-Thomas finite element unknowns: "
158 << trial_size << endl;
159 cout <<
"Number of L2 finite element unknowns: " << test_size << endl;
229 a_mixed.
Mult(gftrial, x);
234 mixed.
Mult(gftrial, x);
284 dlo.
Mult(gftrial, discreteInterpolant);
308 double errInterp = discreteInterpolant.
ComputeL2Error(gradp_coef);
311 cout <<
"\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl): "
312 "|| E_h - grad p ||_{L_2} = " << errSol <<
'\n' << endl;
313 cout <<
" Gradient interpolant E_h = grad p_h in H(curl): || E_h - grad p"
314 " ||_{L_2} = " << errInterp <<
'\n' << endl;
315 cout <<
" Projection E_h of exact grad p in H(curl): || E_h - grad p "
316 "||_{L_2} = " << errProj <<
'\n' << endl;
321 double errInterp = discreteInterpolant.
ComputeL2Error(curlv_coef);
324 cout <<
"\n Solution of (E_h,w) = (curl v_h,w) for E_h and w in H(div): "
325 "|| E_h - curl v ||_{L_2} = " << errSol <<
'\n' << endl;
326 cout <<
" Curl interpolant E_h = curl v_h in H(div): || E_h - curl v "
327 "||_{L_2} = " << errInterp <<
'\n' << endl;
328 cout <<
" Projection E_h of exact curl v in H(div): || E_h - curl v "
329 "||_{L_2} = " << errProj <<
'\n' << endl;
333 int order_quad = max(2, 2*order+1);
335 for (
int i=0; i < Geometry::NumGeom; ++i)
341 double errInterp = discreteInterpolant.
ComputeL2Error(divgradp_coef, irs);
344 cout <<
"\n Solution of (f_h,q) = (div v_h,q) for f_h and q in L_2: "
345 "|| f_h - div v ||_{L_2} = " << errSol <<
'\n' << endl;
346 cout <<
" Divergence interpolant f_h = div v_h in L_2: || f_h - div v "
347 "||_{L_2} = " << errInterp <<
'\n' << endl;
348 cout <<
" Projection f_h of exact div v in L_2: || f_h - div v "
349 "||_{L_2} = " << errProj <<
'\n' << endl;
354 ofstream mesh_ofs(
"refined.mesh");
355 mesh_ofs.precision(8);
356 mesh->
Print(mesh_ofs);
357 ofstream sol_ofs(
"sol.gf");
358 sol_ofs.precision(8);
367 sol_sock.precision(8);
368 sol_sock <<
"solution\n" << *mesh << x << flush;
383 return sin(x(0)) *
sin(x(1)) *
sin(x(2));
387 return sin(x(0)) *
sin(x(1));
405 if (x.
Size() == 3) {
f(2) = 0.0; }
413 return -3.0 *
sin(x(0)) *
sin(x(1)) *
sin(x(2));
417 return -2.0 *
sin(x(0)) *
sin(x(1));
435 if (x.
Size() == 3) { v(2) = 0.0; }
Conjugate gradient method.
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
A coefficient that is constant across space and time.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
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 ...
void curlv_exact(const Vector &x, Vector &cv)
Jacobi smoothing for a given bilinear form (no matrix necessary).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetMaxIter(int max_it)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void v_exact(const Vector &x, Vector &v)
void PrintUsage(std::ostream &out) const
Print the usage message.
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
double p_exact(const Vector &x)
A general vector function coefficient.
int SpaceDimension() const
void SetRelTol(double rtol)
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
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...
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
virtual void Print(std::ostream &os=mfem::out) const
void PrintOptions(std::ostream &out) const
Print the options.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
A general function coefficient.
Arbitrary order H(curl)-conforming Nedelec finite elements.
double div_gradp_exact(const Vector &x)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
void gradp_exact(const Vector &, Vector &)
bool Good() const
Return true if the command line options were parsed successfully.