MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::MomentFittingIntRules Class Reference

Class for subdomain IntegrationRules by means of moment-fitting. More...

#include <intrules_cut.hpp>

Inheritance diagram for mfem::MomentFittingIntRules:
[legend]
Collaboration diagram for mfem::MomentFittingIntRules:
[legend]

Public Member Functions

 MomentFittingIntRules (int order, Coefficient &lvlset, int lsO)
 Constructor to set up the generated cut IntegrationRules.
 
void SetOrder (int order) override
 Change the order of the constructed IntegrationRule.
 
void GetSurfaceIntegrationRule (ElementTransformation &Tr, IntegrationRule &result) override
 Construct a cut-surface IntegrationRule.
 
void GetVolumeIntegrationRule (ElementTransformation &Tr, IntegrationRule &result, const IntegrationRule *sir=nullptr) override
 Construct a cut-volume IntegrationRule.
 
void GetSurfaceWeights (ElementTransformation &Tr, const IntegrationRule &sir, Vector &weights) override
 Compute transformation quadrature weights for surface integration.
 
 ~MomentFittingIntRules () override
 Destructor of MomentFittingIntRules.
 
virtual void SetLevelSetCoefficient (Coefficient &ls)
 Change the Coefficient whose zero level set specifies the cut.
 
virtual void SetLevelSetProjectionOrder (int order)
 
- Public Member Functions inherited from mfem::CutIntegrationRules
virtual ~CutIntegrationRules ()
 Destructor of CutIntegrationRules.
 

Protected Member Functions

void InitSurface (int order, Coefficient &levelset, int lsO, ElementTransformation &Tr)
 Initialization for surface IntegrationRule.
 
void InitVolume (int order, Coefficient &levelset, int lsO, ElementTransformation &Tr)
 Initialization for volume IntegrationRule.
 
void Init (int order, Coefficient &levelset, int lsO)
 Initialize the MomentFittingIntRules.
 
void Clear ()
 Clear stored data of the MomentFittingIntRules.
 
void ComputeFaceWeights (ElementTransformation &Tr)
 Compute the IntegrationRules on the faces.
 
void ComputeSurfaceWeights1D (ElementTransformation &Tr)
 Compute 1D quadrature weights.
 
void ComputeVolumeWeights1D (ElementTransformation &Tr, const IntegrationRule *sir)
 Compute the 1D quadrature weights.
 
void ComputeSurfaceWeights2D (ElementTransformation &Tr)
 Compute 2D quadrature weights.
 
void ComputeVolumeWeights2D (ElementTransformation &Tr, const IntegrationRule *sir)
 Compute the 2D quadrature weights.
 
void ComputeSurfaceWeights3D (ElementTransformation &Tr)
 Compute 2D quadrature weights.
 
void ComputeVolumeWeights3D (ElementTransformation &Tr, const IntegrationRule *sir)
 Compute the 3D quadrature weights.
 
void DivFreeBasis2D (const IntegrationPoint &ip, DenseMatrix &shape)
 A divergence free basis on the element [-1,1]^2.
 
void OrthoBasis2D (const IntegrationPoint &ip, DenseMatrix &shape)
 A orthogonalized divergence free basis on the element [-1,1]^2.
 
void OrthoBasis3D (const IntegrationPoint &ip, DenseMatrix &shape)
 A orthogonalized divergence free basis on the element [-1,1]^3.
 
void mGSStep (DenseMatrix &shape, DenseTensor &shapeMFN, int step)
 A step of the modified Gram-Schmidt algorithm.
 
void Basis2D (const IntegrationPoint &ip, Vector &shape)
 Monomial basis on the element [-1,1]^2.
 
void BasisAD2D (const IntegrationPoint &ip, DenseMatrix &shape)
 Antiderivatives of the monomial basis on the element [-1,1]^2.
 
void Basis3D (const IntegrationPoint &ip, Vector &shape)
 Monomial basis on the element [-1,1]^3.
 
void BasisAD3D (const IntegrationPoint &ip, DenseMatrix &shape)
 Antiderivatives of the monomial basis on the element [-1,1]^3.
 
- Protected Member Functions inherited from mfem::CutIntegrationRules
 CutIntegrationRules (int order, Coefficient &lvlset, int lsO=2)
 Constructor to set up the generated cut IntegrationRules.
 

Protected Attributes

int dim
 Space Dimension of the element.
 
int nBasis
 Number of divergence-free basis functions for surface integration.
 
int nBasisVolume
 Number of basis functions for volume integration.
 
IntegrationRule ir
 IntegrationRule representing the reused IntegrationPoints.
 
DenseMatrixSVDVolumeSVD
 SVD of the matrix for volumetric IntegrationRules.
 
Array< IntegrationPointFaceIP
 Array of face integration points.
 
DenseMatrix FaceWeights
 Column-wise Matrix for the face quadrature weights.
 
Vector FaceWeightsComp
 Indicates the already computed face IntegrationRules.
 
- Protected Attributes inherited from mfem::CutIntegrationRules
int Order
 Order of the IntegrationRule.
 
CoefficientLvlSet
 The zero level set of this Coefficient defines the cut surface.
 
int lsOrder
 Space order for the LS projection.
 

Detailed Description

Class for subdomain IntegrationRules by means of moment-fitting.

Class for subdomain (surface and subdomain) IntegrationRules by means of moment-fitting. The class provides different functions to construct the IntegrationRules. The construction is done as described in Mueller et al. in "Highly accurate surface and volume integration on implicit domains by means of moment-fitting" (2013, see https://onlinelibrary.wiley.com/doi/full/10.1002/nme.4569).

Definition at line 120 of file intrules_cut.hpp.

Constructor & Destructor Documentation

◆ MomentFittingIntRules()

mfem::MomentFittingIntRules::MomentFittingIntRules ( int order,
Coefficient & lvlset,
int lsO )
inline

Constructor to set up the generated cut IntegrationRules.

Parameters
[in]orderOrder of the constructed IntegrationRule.
[in]lvlsetCoefficient whose zero level set specifies the cut.
[in]lsOPolynomial degree for projecting the level-set Coefficient to a GridFunction, which is used to compute gradients and normals.

Definition at line 289 of file intrules_cut.hpp.

◆ ~MomentFittingIntRules()

mfem::MomentFittingIntRules::~MomentFittingIntRules ( )
inlineoverride

Destructor of MomentFittingIntRules.

Definition at line 343 of file intrules_cut.hpp.

Member Function Documentation

◆ Basis2D()

void mfem::MomentFittingIntRules::Basis2D ( const IntegrationPoint & ip,
Vector & shape )
protected

Monomial basis on the element [-1,1]^2.

Definition at line 1320 of file intrules_cut.cpp.

◆ Basis3D()

void mfem::MomentFittingIntRules::Basis3D ( const IntegrationPoint & ip,
Vector & shape )
protected

Monomial basis on the element [-1,1]^3.

Definition at line 1365 of file intrules_cut.cpp.

◆ BasisAD2D()

void mfem::MomentFittingIntRules::BasisAD2D ( const IntegrationPoint & ip,
DenseMatrix & shape )
protected

Antiderivatives of the monomial basis on the element [-1,1]^2.

Definition at line 1340 of file intrules_cut.cpp.

◆ BasisAD3D()

void mfem::MomentFittingIntRules::BasisAD3D ( const IntegrationPoint & ip,
DenseMatrix & shape )
protected

Antiderivatives of the monomial basis on the element [-1,1]^3.

Definition at line 1386 of file intrules_cut.cpp.

◆ Clear()

void mfem::MomentFittingIntRules::Clear ( )
protected

Clear stored data of the MomentFittingIntRules.

Definition at line 1417 of file intrules_cut.cpp.

◆ ComputeFaceWeights()

void mfem::MomentFittingIntRules::ComputeFaceWeights ( ElementTransformation & Tr)
protected

Compute the IntegrationRules on the faces.

Compute the IntegrationRules on the (cut) faces of an element. These will be saved and reused if possible.

Parameters
[in]TrElementTransformation of the element the IntegrationRules on the faces are computed

Definition at line 128 of file intrules_cut.cpp.

◆ ComputeSurfaceWeights1D()

void mfem::MomentFittingIntRules::ComputeSurfaceWeights1D ( ElementTransformation & Tr)
protected

Compute 1D quadrature weights.

Compute the quadrature weights for the 1D surface quadrature rule.

Parameters
[in]TrElementTransformation of the current element

Definition at line 210 of file intrules_cut.cpp.

◆ ComputeSurfaceWeights2D()

void mfem::MomentFittingIntRules::ComputeSurfaceWeights2D ( ElementTransformation & Tr)
protected

Compute 2D quadrature weights.

Compute the quadrature weights for the 2D surface quadrature rule by means of moment-fitting. To construct the quadrature rule, special integrals are reduced to integrals over the edges of the subcell where the level-set is positive.

Parameters
[in]TrElementTransformation of the current element

Definition at line 309 of file intrules_cut.cpp.

◆ ComputeSurfaceWeights3D()

void mfem::MomentFittingIntRules::ComputeSurfaceWeights3D ( ElementTransformation & Tr)
protected

Compute 2D quadrature weights.

Compute the quadrature weights for the 3D surface quadrature rule by means of moment-fitting. To construct the quadrature rule, special integrals are reduced to integrals over the edges of the subcell where the level-set is positive.

Parameters
[in]TrElementTransformation of the current element

Definition at line 818 of file intrules_cut.cpp.

◆ ComputeVolumeWeights1D()

void mfem::MomentFittingIntRules::ComputeVolumeWeights1D ( ElementTransformation & Tr,
const IntegrationRule * sir )
protected

Compute the 1D quadrature weights.

Compute the 1D quadrature weights for the volumetric subdomain quadrature rule.

Parameters
[in]TrElementTransformation of the current element
[in]sircorresponding IntegrationRule on surface

Definition at line 257 of file intrules_cut.cpp.

◆ ComputeVolumeWeights2D()

void mfem::MomentFittingIntRules::ComputeVolumeWeights2D ( ElementTransformation & Tr,
const IntegrationRule * sir )
protected

Compute the 2D quadrature weights.

Compute the 2D quadrature weights for the volumetric subdomain quadrature rule by means of moment-fitting. To construct the quadrature rule, special integrals are reduced to integrals over the boundary of the subcell where the level-set is positive.

Parameters
[in]TrElementTransformation of the current element
[in]sircorresponding IntegrationRule on surface

Definition at line 560 of file intrules_cut.cpp.

◆ ComputeVolumeWeights3D()

void mfem::MomentFittingIntRules::ComputeVolumeWeights3D ( ElementTransformation & Tr,
const IntegrationRule * sir )
protected

Compute the 3D quadrature weights.

Compute the 3D quadrature weights for the volumetric subdomain quadrature rule by means of moment-fitting. To construct the quadrature rule, special integrals are reduced to integrals over the boundary of the subcell where the level-set is positive.

Parameters
[in]TrElementTransformation of the current element
[in]sircorresponding IntegrationRule on surface

Definition at line 999 of file intrules_cut.cpp.

◆ DivFreeBasis2D()

void mfem::MomentFittingIntRules::DivFreeBasis2D ( const IntegrationPoint & ip,
DenseMatrix & shape )
protected

A divergence free basis on the element [-1,1]^2.

Definition at line 1192 of file intrules_cut.cpp.

◆ GetSurfaceIntegrationRule()

void mfem::MomentFittingIntRules::GetSurfaceIntegrationRule ( ElementTransformation & Tr,
IntegrationRule & result )
overridevirtual

Construct a cut-surface IntegrationRule.

Construct an IntegrationRule to integrate on the surface given by the already specified level set function, for the element given by Tr.

Parameters
[in]TrSpecifies the IntegrationRule's associated mesh element.
[out]resultIntegrationRule on the cut-surface

Implements mfem::CutIntegrationRules.

Definition at line 1435 of file intrules_cut.cpp.

◆ GetSurfaceWeights()

void mfem::MomentFittingIntRules::GetSurfaceWeights ( ElementTransformation & Tr,
const IntegrationRule & sir,
Vector & weights )
overridevirtual

Compute transformation quadrature weights for surface integration.

Compute the transformation weights for integration over the cut-surface in reference space.

Parameters
[in]TrSpecifies the IntegrationRule's associated element.
[in]sirIntegrationRule defining the IntegrationPoints
[out]weightsVector containing the transformation weights.

Implements mfem::CutIntegrationRules.

Definition at line 1536 of file intrules_cut.cpp.

◆ GetVolumeIntegrationRule()

void mfem::MomentFittingIntRules::GetVolumeIntegrationRule ( ElementTransformation & Tr,
IntegrationRule & result,
const IntegrationRule * sir = nullptr )
overridevirtual

Construct a cut-volume IntegrationRule.

Construct an IntegrationRule to integrate in the subdomain given by the positive values of the already specified level set function, for the element given by Tr.

Parameters
[in]TrSpecifies the IntegrationRule's associated mesh element.
[out]resultIntegrationRule for the cut-volume
[in]sirCorresponding IntegrationRule for the surface, which can be used to avoid computations.

Implements mfem::CutIntegrationRules.

Definition at line 1476 of file intrules_cut.cpp.

◆ Init()

void mfem::MomentFittingIntRules::Init ( int order,
Coefficient & levelset,
int lsO )
inlineprotected

Initialize the MomentFittingIntRules.

Definition at line 171 of file intrules_cut.hpp.

◆ InitSurface()

void mfem::MomentFittingIntRules::InitSurface ( int order,
Coefficient & levelset,
int lsO,
ElementTransformation & Tr )
protected

Initialization for surface IntegrationRule.

Initialize the members for computation of surface IntegrationRule. This is called when the first IntegrationRule is computed or when Order or level-set change.

Parameters
[in]orderOrder of the IntegrationRule
[in]levelsetlevel-set function defining the implicit interface
[in]lsOpolynomial degree for approximation of level-set function
[in]TrElemenTransformation to initalize the members with

Definition at line 36 of file intrules_cut.cpp.

◆ InitVolume()

void mfem::MomentFittingIntRules::InitVolume ( int order,
Coefficient & levelset,
int lsO,
ElementTransformation & Tr )
protected

Initialization for volume IntegrationRule.

Initialize the members for computation of surface IntegrationRule. This is called when the first IntegrationRule is computed or when Order or level-set change.

Parameters
[in]orderOrder of the IntegrationRule
[in]levelsetlevel-set function defining the implicit interface
[in]lsOpolynomial degree for approximation of level-set function
[in]TrElemenTransformation to initalize the members with

Definition at line 76 of file intrules_cut.cpp.

◆ mGSStep()

void mfem::MomentFittingIntRules::mGSStep ( DenseMatrix & shape,
DenseTensor & shapeMFN,
int step )
protected

A step of the modified Gram-Schmidt algorithm.

Definition at line 1277 of file intrules_cut.cpp.

◆ OrthoBasis2D()

void mfem::MomentFittingIntRules::OrthoBasis2D ( const IntegrationPoint & ip,
DenseMatrix & shape )
protected

A orthogonalized divergence free basis on the element [-1,1]^2.

Definition at line 1235 of file intrules_cut.cpp.

◆ OrthoBasis3D()

void mfem::MomentFittingIntRules::OrthoBasis3D ( const IntegrationPoint & ip,
DenseMatrix & shape )
protected

A orthogonalized divergence free basis on the element [-1,1]^3.

Definition at line 1266 of file intrules_cut.cpp.

◆ SetLevelSetCoefficient()

virtual void mfem::CutIntegrationRules::SetLevelSetCoefficient ( Coefficient & ls)
inlinevirtual

Change the Coefficient whose zero level set specifies the cut.

Reimplemented from mfem::CutIntegrationRules.

Definition at line 56 of file intrules_cut.hpp.

◆ SetLevelSetProjectionOrder()

void mfem::CutIntegrationRules::SetLevelSetProjectionOrder ( int order)
virtual

Change the polynomial degree for projecting the level set Coefficient to a GridFunction, which is used to compute local gradients and normals.

Reimplemented from mfem::CutIntegrationRules.

Definition at line 60 of file intrules_cut.cpp.

◆ SetOrder()

void mfem::MomentFittingIntRules::SetOrder ( int order)
overridevirtual

Change the order of the constructed IntegrationRule.

Reimplemented from mfem::CutIntegrationRules.

Definition at line 1429 of file intrules_cut.cpp.

Member Data Documentation

◆ dim

int mfem::MomentFittingIntRules::dim
protected

Space Dimension of the element.

Definition at line 124 of file intrules_cut.hpp.

◆ FaceIP

Array<IntegrationPoint> mfem::MomentFittingIntRules::FaceIP
protected

Array of face integration points.

Definition at line 134 of file intrules_cut.hpp.

◆ FaceWeights

DenseMatrix mfem::MomentFittingIntRules::FaceWeights
protected

Column-wise Matrix for the face quadrature weights.

Definition at line 136 of file intrules_cut.hpp.

◆ FaceWeightsComp

Vector mfem::MomentFittingIntRules::FaceWeightsComp
protected

Indicates the already computed face IntegrationRules.

Definition at line 138 of file intrules_cut.hpp.

◆ ir

IntegrationRule mfem::MomentFittingIntRules::ir
protected

IntegrationRule representing the reused IntegrationPoints.

Definition at line 130 of file intrules_cut.hpp.

◆ nBasis

int mfem::MomentFittingIntRules::nBasis
protected

Number of divergence-free basis functions for surface integration.

Definition at line 126 of file intrules_cut.hpp.

◆ nBasisVolume

int mfem::MomentFittingIntRules::nBasisVolume
protected

Number of basis functions for volume integration.

Definition at line 128 of file intrules_cut.hpp.

◆ VolumeSVD

DenseMatrixSVD* mfem::MomentFittingIntRules::VolumeSVD
protected

SVD of the matrix for volumetric IntegrationRules.

Definition at line 132 of file intrules_cut.hpp.


The documentation for this class was generated from the following files: