55int main(
int argc,
char *argv[])
58 const char *mesh_file =
"../../data/cube-nurbs.mesh";
63 bool static_cond =
false;
65 const char *device_config =
"cpu";
67 bool visualization = 1;
70 args.
AddOption(&mesh_file,
"-m",
"--mesh",
72 args.
AddOption(&ref_levels,
"-r",
"--refine",
73 "Number of times to refine the mesh uniformly, -1 for auto.");
75 "Finite element order (polynomial degree).");
76 args.
AddOption(&NURBS,
"-n",
"--nurbs",
"-nn",
"--no-nurbs",
79 "Choose between 0: grad, 1: curl, 2: div");
80 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
81 "--no-static-condensation",
"Enable static condensation.");
82 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
83 "--no-partial-assembly",
"Enable Partial Assembly.");
84 args.
AddOption(&device_config,
"-d",
"--device",
85 "Device configuration string, see Device::Configure().");
86 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
88 "Enable or disable GLVis visualization.");
90 args.
AddOption(&visport,
"-p",
"--send-port",
"Socket for GLVis.");
102 Device device(device_config);
108 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
112 MFEM_ABORT(
"The curl problem is only defined in 3D.");
123 ref_levels = (int)floor(log(50000./mesh->
GetNE())/log(2.)/
dim);
125 for (
int l = 0; l < ref_levels; l++)
154 mfem::out<<
"Create NURBS fec and ext"<<std::endl;
183 cout <<
"Number of Nedelec finite element unknowns: " << test_size << endl;
184 cout <<
"Number of H1 finite element unknowns: " << trial_size << endl;
188 cout <<
"Number of Nedelec finite element unknowns: " << trial_size << endl;
189 cout <<
"Number of Raviart-Thomas finite element unknowns: " << test_size <<
194 cout <<
"Number of Raviart-Thomas finite element unknowns: "
195 << trial_size << endl;
196 cout <<
"Number of L2 finite element unknowns: " << test_size << endl;
232 a.SetAssemblyLevel(AssemblyLevel::PARTIAL);
256 if (static_cond) {
a.EnableStaticCondensation(); }
259 if (!pa) {
a.Finalize(); }
266 a_mixed.
Mult(gftrial, x);
271 mixed.
Mult(gftrial, x);
328 cout <<
"\n Solution of (E_h,v) = (grad p_h,v) for E_h and v in H(curl): "
329 "|| E_h - grad p ||_{L_2} = " << errSol <<
'\n' << endl;
330 cout <<
" Projection E_h of exact grad p in H(curl): || E_h - grad p "
331 "||_{L_2} = " << errProj <<
'\n' << endl;
338 cout <<
"\n Solution of (E_h,w) = (curl v_h,w) for E_h and w in H(div): "
339 "|| E_h - curl v ||_{L_2} = " << errSol <<
'\n' << endl;
340 cout <<
" Projection E_h of exact curl v in H(div): || E_h - curl v "
341 "||_{L_2} = " << errProj <<
'\n' << endl;
345 int order_quad = max(3, 2*order+1);
355 cout <<
"\n Solution of (f_h,q) = (div v_h,q) for f_h and q in L_2: "
356 "|| f_h - div v ||_{L_2} = " << errSol <<
'\n' << endl;
358 cout <<
" Projection f_h of exact div v in L_2: || f_h - div v "
359 "||_{L_2} = " << errProj <<
'\n' << endl;
364 ofstream mesh_ofs(
"refined.mesh");
365 mesh_ofs.precision(8);
366 mesh->
Print(mesh_ofs);
367 ofstream sol_ofs(
"sol.gf");
368 sol_ofs.precision(8);
376 sol_sock.precision(8);
377 sol_sock <<
"solution\n" << *mesh << x << flush;
392 return sin(x(0)) * sin(x(1)) * sin(x(2));
396 return sin(x(0)) * sin(x(1));
406 f(0) = cos(x(0)) * sin(x(1)) * sin(x(2));
407 f(1) = sin(x(0)) * cos(x(1)) * sin(x(2));
408 f(2) = sin(x(0)) * sin(x(1)) * cos(x(2));
412 f(0) = cos(x(0)) * sin(x(1));
413 f(1) = sin(x(0)) * cos(x(1));
414 if (x.
Size() == 3) {
f(2) = 0.0; }
422 return -3.0 * sin(x(0)) * sin(x(1)) * sin(x(2));
426 return -2.0 * sin(x(0)) * sin(x(1));
436 v(0) = sin(
kappa * x(1));
437 v(1) = sin(
kappa * x(2));
438 v(2) = sin(
kappa * x(0));
442 v(0) = sin(
kappa * x(1));
443 v(1) = sin(
kappa * x(0));
444 if (x.
Size() == 3) { v(2) = 0.0; }
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
A coefficient that is constant across space and time.
Data type for scaled Jacobi-type smoother of sparse matrix.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
NURBSExtension * StealNURBSext()
A general function coefficient.
Class for grid function - Vector with associated FE space.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||exsol - u_h||_L2 for scalar or vector H1 or L2 elements.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Arbitrary order H1-conforming (continuous) finite elements.
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Arbitrary order "L2-conforming" discontinuous finite elements.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Arbitrary order H(curl)-conforming Nedelec finite elements.
NURBSExtension generally contains multiple NURBSPatch objects spanning an entire Mesh....
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Arbitrary order H(curl) NURBS finite elements.
Arbitrary order H(div) NURBS finite elements.
Jacobi smoothing for a given bilinear form (no matrix necessary).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
A general vector function coefficient.
int Size() const
Returns the size of the vector.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
std::function< real_t(const Vector &)> f(real_t mass_coeff)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p_exact(const Vector &x)
void curlv_exact(const Vector &x, Vector &cv)
void v_exact(const Vector &x, Vector &v)
void gradp_exact(const Vector &, Vector &)
real_t div_gradp_exact(const Vector &x)