56 return 1. - (pow(X(0), 2.) + pow(X(1), 2.));
58 return 1. - (pow(X(0) / 1.5, 2.) + pow(X(1) / .75, 2.));
60 return 1. - (pow(X(0), 2.) + pow(X(1), 2.) + pow(X(2), 2.));
62 return 1. - (pow(X(0) / 1.5, 2.) + pow(X(1) / .75, 2.) + pow(X(2) / .5, 2.));
76 return 3. * pow(X(0), 2.) - pow(X(1), 2.);
80 return 4. - 3. * pow(X(0), 2.) + 2. * pow(X(1), 2.) - pow(X(2), 2.);
98 return 7.26633616541076;
100 return 40. / 3. * M_PI;
102 return 9.90182151329315;
114 return pow(.55, 3.) / 3.;
118 return 9. / 8. * M_PI;
122 return 3. / 4. * M_PI;
138 int method, ir_order, ls_order;
158 SIntegrationRule(
int method_,
int Order,
160 : method(method_), ir_order(Order), ls_order(lsOrder),
161 level_set(LvlSet),
dim(mesh->Dimension())
164 if (method == 1) {
return; }
166#ifdef MFEM_USE_LAPACK
172 mf_ir.GetSurfaceIntegrationRule(Tr, ir);
183 mf_ir.GetSurfaceWeights(Tr, ir, w);
184 SurfaceWeights.
SetCol(0, w);
206 for (
int elem = 1; elem < mesh->
GetNE(); elem++)
209 mf_ir.GetSurfaceIntegrationRule(Tr, ir);
210 mf_ir.GetSurfaceWeights(Tr, ir, w);
211 SurfaceWeights.
SetCol(elem, w);
227 MFEM_ABORT(
"Moment-fitting requires MFEM to be built with LAPACK!");
239#ifdef MFEM_USE_ALGOIM
241 a_ir.GetSurfaceIntegrationRule(Tr, *
this);
243 a_ir.GetSurfaceWeights(Tr, *
this, w);
250 MFEM_ABORT(
"MFEM is not built with Algoim support!");
280 int method, ir_order, ls_order;
298 CIntegrationRule(
int method_,
int Order,
300 : method(method_), ir_order(Order), ls_order(lsOrder),
301 level_set(LvlSet),
dim(mesh->Dimension())
304 if (method == 1) {
return; }
306#ifdef MFEM_USE_LAPACK
312 mf_ir.GetVolumeIntegrationRule(Tr, ir);
341 for (
int elem = 1; elem < mesh->
GetNE(); elem++)
344 mf_ir.GetVolumeIntegrationRule(Tr, ir);
354 Weights(2 * ip, elem) = ir.
IntPoint(ip).
x;
360 MFEM_ABORT(
"Moment-fitting requires MFEM to be built with LAPACK!");
369#ifdef MFEM_USE_ALGOIM
371 a_ir.GetVolumeIntegrationRule(Tr, *
this);
374 MFEM_ABORT(
"MFEM is not built with Algoim support!");
405 SIntegrationRule* SIntRule;
425 SIntegrationRule* ir)
448 SIntRule->SetElementAndSurfaceWeight(Tr);
450 for (
int ip = 0; ip < SIntRule->
GetNPoints(); ip++)
475 CIntegrationRule* CIntRule;
495 CIntegrationRule* ir)
518 CIntRule->SetElement(Tr);
520 for (
int ip = 0; ip < CIntRule->
GetNPoints(); ip++)
533int main(
int argc,
char *argv[])
535#if defined(MFEM_USE_LAPACK) || defined(MFEM_USE_ALGOIM)
540 const char *inttype =
"surface2d";
541 bool visualization =
true;
545 args.
AddOption(&order,
"-o",
"--order",
"Order of quadrature rule");
546 args.
AddOption(&ref_levels,
"-r",
"--refine",
"Number of meh refinements");
547 args.
AddOption(&method,
"-m",
"--method",
548 "Cut integration method: 0 for moments-based, 1 for Algoim.");
549 args.
AddOption(&inttype,
"-i",
"--integration-type",
550 "IntegrationType to demonstrate");
551 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
552 "--no-visualization",
553 "Enable or disable GLVis visualization.");
556 if (strcmp(inttype,
"volumetric1d") == 0
557 || strcmp(inttype,
"Volumetric1D") == 0)
561 else if (strcmp(inttype,
"surface2d") == 0
562 || strcmp(inttype,
"Surface2D") == 0)
566 else if (strcmp(inttype,
"volumetric2d") == 0
567 || strcmp(inttype,
"Volumetric2D") == 0)
571 else if (strcmp(inttype,
"surface3d") == 0
572 || strcmp(inttype,
"Surface3d") == 0)
576 else if (strcmp(inttype,
"volumetric3d") == 0
577 || strcmp(inttype,
"Volumetric3d") == 0)
583 Mesh *mesh =
nullptr;
586 mesh =
new Mesh(
"../data/inline-segment.mesh");
591 mesh =
new Mesh(2, 4, 1, 0, 2);
602 mesh =
new Mesh(3, 8, 1, 0, 3);
611 mesh->
AddHex(0,1,2,3,4,5,6,7);
615 for (
int lev = 0; lev < ref_levels; lev++)
631 SIntegrationRule* sir =
new SIntegrationRule(method, order,
633 CIntegrationRule* cir = NULL;
638 cir =
new CIntegrationRule(method, order, levelset, 2, mesh);
658 int nbasis = 2 * (order + 1) + (
int)(order * (order + 1) / 2);
665 cout <<
"============================================" << endl;
666 cout <<
"Mesh size dx: ";
669 cout << 3.2 / pow(2., (
real_t)ref_levels) << endl;
673 cout << .25 / pow(2., (
real_t)ref_levels) << endl;
678 cout <<
"Number of div free basis functions: " << nbasis << endl;
679 cout <<
"Number of quadrature points: " << ir.
GetNPoints() << endl;
681 cout << scientific << setprecision(10);
682 cout <<
"============================================" << endl;
683 cout <<
"Computed value of surface integral: " << surface.
Sum() << endl;
684 cout <<
"True value of surface integral: " <<
Surface() << endl;
685 cout <<
"Absolute Error (Surface): ";
686 cout << abs(surface.
Sum() - Surface()) << endl;
687 cout <<
"Relative Error (Surface): ";
688 cout << abs(surface.
Sum() - Surface()) /
Surface() << endl;
693 cout <<
"--------------------------------------------" << endl;
694 cout <<
"Computed value of volume integral: " << volume.
Sum() << endl;
695 cout <<
"True value of volume integral: " <<
Volume() << endl;
696 cout <<
"Absolute Error (Volume): ";
697 cout << abs(volume.
Sum() -
Volume()) << endl;
698 cout <<
"Relative Error (Volume): ";
701 cout <<
"============================================" << endl;
714 sol_sock.precision(8);
715 sol_sock <<
"solution\n" << *mesh << lgf << flush;
716 sol_sock <<
"keys pppppppppppppppppppppppppppcmmlRj\n";
717 sol_sock <<
"levellines " << 0. <<
" " << 0. <<
" " << 1 <<
"\n" << flush;
726 cout <<
"MFEM must be built with LAPACK or ALGOIM for this example." << endl;
727 return MFEM_SKIP_RETURN_VALUE;
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
Data type dense matrix using column-major storage.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
void SetCol(int c, const real_t *col)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Abstract class for all finite elements.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
A general function coefficient.
Class for grid function - Vector with associated FE space.
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 integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Container class for integration rules.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Adds a quadrilateral to the mesh given by 4 vertices v1 through v4.
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Adds a hexahedron to the mesh given by 8 vertices v1 through v8.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Class for subdomain IntegrationRules by means of moment-fitting.
void ParseCheck(std::ostream &out=mfem::out)
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...
real_t Sum() const
Return the sum of the vector entries.
void SetSize(int s)
Resize the vector to size s.
real_t Surface()
Analytic surface integral.
real_t integrand(const Vector &X)
Function that should be integrated.
real_t Volume()
Analytic volume integral over subdomain with positive level-set.
real_t lvlset(const Vector &X)
Level-set function defining the implicit interface.
IntegrationType
Integration rule the example should demonstrate.
real_t u(const Vector &xvec)
void add(const Vector &v1, const Vector &v2, Vector &v)