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

#include <bilininteg.hpp>

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

Public Member Functions

 VectorMassIntegrator ()
 Construct an integrator with coefficient 1.0.
 
 VectorMassIntegrator (Coefficient &q, int qo=0)
 
 VectorMassIntegrator (Coefficient &q, const IntegrationRule *ir)
 
 VectorMassIntegrator (VectorCoefficient &q, int qo=0)
 Construct an integrator with diagonal coefficient q.
 
 VectorMassIntegrator (MatrixCoefficient &q, int qo=0)
 Construct an integrator with matrix coefficient q.
 
int GetVDim () const
 
void SetVDim (int vdim_)
 
virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 Given a particular Finite Element computes the element matrix elmat.
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly.
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly.
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action.
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssembleNURBSPA (const FiniteElementSpace &fes)
 Method defining partial assembly on NURBS patches.
 
virtual void AssemblePABoundary (const FiniteElementSpace &fes)
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag)
 Assemble diagonal of \(A D A^T\) ( \(A\) is this integrator) and add it to diag.
 
virtual void AddMultNURBSPA (const Vector &x, Vector &y) const
 Method for partially assembled action on NURBS patches.
 
virtual void AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action.
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true)
 Method defining element assembly.
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEABoundaryFaces (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
virtual void AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
virtual void AssembleFaceMatrix (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleFaceMatrix (const FiniteElement &trial_face_fe, const FiniteElement &test_fe1, const FiniteElement &test_fe2, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleTraceFaceMatrix (int elem, const FiniteElement &trial_face_fe, const FiniteElement &test_fe, FaceElementTransformations &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator resulting from a face integral term. Note that the default implementation in the base class is general but not efficient.
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix.
 
virtual void AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term.
 
virtual void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual real_t ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
virtual bool RequiresFaceNormalDerivatives () const
 For bilinear forms on element faces, specifies if the normal derivatives are needed on the faces or just the face restriction.
 
virtual void AddMultPAFaceNormalDerivatives (const Vector &x, const Vector &dxdn, Vector &y, Vector &dydn) const
 Method for partially assembled action.
 
virtual ~BilinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL).
 
void SetIntegrationMode (Mode m)
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 For patchwise integration, SetNURBSPatchIntRule must be called.
 
bool HasNURBSPatchIntRule () const
 
bool Patchwise () const
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use.
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL).
 
virtual real_t GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy.
 
virtual void AssembleGradPA (const Vector &x, const FiniteElementSpace &fes)
 Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.
 
virtual real_t GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly.
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action.
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly.
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Attributes

CoefficientQ
 
VectorCoefficientVQ
 
MatrixCoefficientMQ
 
Vector pa_data
 
const DofToQuadmaps
 Not owned.
 
const GeometricFactorsgeom
 Not owned.
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int quad1D
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 
- Protected Member Functions inherited from mfem::BilinearFormIntegrator
 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Detailed Description

Class for integrating the bilinear form \(a(u,v) := (Q u, v)\), where \(u=(u_1,\dots,u_n)\) and \(v=(v_1,\dots,v_n)\), \(u_i\) and \(v_i\) are defined by scalar FE through standard transformation.

Definition at line 2456 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ VectorMassIntegrator() [1/5]

mfem::VectorMassIntegrator::VectorMassIntegrator ( )
inline

Construct an integrator with coefficient 1.0.

Definition at line 2477 of file bilininteg.hpp.

◆ VectorMassIntegrator() [2/5]

mfem::VectorMassIntegrator::VectorMassIntegrator ( Coefficient & q,
int qo = 0 )
inline

Construct an integrator with scalar coefficient q. If possible, save memory by using a scalar integrator since the resulting matrix is block diagonal with the same diagonal block repeated.

Definition at line 2482 of file bilininteg.hpp.

◆ VectorMassIntegrator() [3/5]

mfem::VectorMassIntegrator::VectorMassIntegrator ( Coefficient & q,
const IntegrationRule * ir )
inline

Definition at line 2484 of file bilininteg.hpp.

◆ VectorMassIntegrator() [4/5]

mfem::VectorMassIntegrator::VectorMassIntegrator ( VectorCoefficient & q,
int qo = 0 )
inline

Construct an integrator with diagonal coefficient q.

Definition at line 2488 of file bilininteg.hpp.

◆ VectorMassIntegrator() [5/5]

mfem::VectorMassIntegrator::VectorMassIntegrator ( MatrixCoefficient & q,
int qo = 0 )
inline

Construct an integrator with matrix coefficient q.

Definition at line 2491 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultMF()

virtual void mfem::VectorMassIntegrator::AddMultMF ( const Vector & x,
Vector & y ) const
virtual

Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.

This method can be called only after the method AssembleMF() has been called.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AddMultPA()

virtual void mfem::VectorMassIntegrator::AddMultPA ( const Vector & x,
Vector & y ) const
virtual

Method for partially assembled action.

Perform the action of integrator on the input x and add the result to the output y. Both x and y are E-vectors, i.e. they represent the element-wise discontinuous version of the FE space.

This method can be called only after the method AssemblePA() has been called.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleDiagonalMF()

virtual void mfem::VectorMassIntegrator::AssembleDiagonalMF ( Vector & diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleDiagonalPA()

virtual void mfem::VectorMassIntegrator::AssembleDiagonalPA ( Vector & diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleElementMatrix()

void mfem::VectorMassIntegrator::AssembleElementMatrix ( const FiniteElement & el,
ElementTransformation & Trans,
DenseMatrix & elmat )
virtual

Given a particular Finite Element computes the element matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1542 of file bilininteg.cpp.

◆ AssembleElementMatrix2()

void mfem::VectorMassIntegrator::AssembleElementMatrix2 ( const FiniteElement & trial_fe,
const FiniteElement & test_fe,
ElementTransformation & Trans,
DenseMatrix & elmat )
virtual

Compute the local matrix representation of a bilinear form \(a(u,v)\) defined on different trial (given by \(u\)) and test (given by \(v\)) spaces. The rows in the local matrix correspond to the test dofs and the columns – to the trial dofs.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1624 of file bilininteg.cpp.

◆ AssembleMF()

virtual void mfem::VectorMassIntegrator::AssembleMF ( const FiniteElementSpace & fes)
virtual

Method defining matrix-free assembly.

Used with BilinearFormIntegrators that have different spaces. The result of fully matrix-free assembly is stored internally so that it can be used later in the methods AddMultMF() and AddMultTransposeMF().

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssemblePA() [1/2]

virtual void mfem::VectorMassIntegrator::AssemblePA ( const FiniteElementSpace & fes)
virtual

Method defining partial assembly.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA() and AddMultTransposePA().

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssemblePA() [2/2]

void mfem::BilinearFormIntegrator::AssemblePA ( const FiniteElementSpace & trial_fes,
const FiniteElementSpace & test_fes )
virtual

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 50 of file bilininteg.cpp.

◆ GetVDim()

int mfem::VectorMassIntegrator::GetVDim ( ) const
inline

Definition at line 2494 of file bilininteg.hpp.

◆ SetVDim()

void mfem::VectorMassIntegrator::SetVDim ( int vdim_)
inline

Definition at line 2495 of file bilininteg.hpp.

◆ SupportsCeed()

bool mfem::VectorMassIntegrator::SupportsCeed ( ) const
inlinevirtual

Indicates whether this integrator can use a Ceed backend.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2511 of file bilininteg.hpp.

Member Data Documentation

◆ dim

int mfem::VectorMassIntegrator::dim
protected

Definition at line 2473 of file bilininteg.hpp.

◆ dofs1D

int mfem::VectorMassIntegrator::dofs1D
protected

Definition at line 2473 of file bilininteg.hpp.

◆ geom

const GeometricFactors* mfem::VectorMassIntegrator::geom
protected

Not owned.

Definition at line 2472 of file bilininteg.hpp.

◆ maps

const DofToQuad* mfem::VectorMassIntegrator::maps
protected

Not owned.

Definition at line 2471 of file bilininteg.hpp.

◆ MQ

MatrixCoefficient* mfem::VectorMassIntegrator::MQ
protected

Definition at line 2468 of file bilininteg.hpp.

◆ ne

int mfem::VectorMassIntegrator::ne
protected

Definition at line 2473 of file bilininteg.hpp.

◆ nq

int mfem::VectorMassIntegrator::nq
protected

Definition at line 2473 of file bilininteg.hpp.

◆ pa_data

Vector mfem::VectorMassIntegrator::pa_data
protected

Definition at line 2470 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::VectorMassIntegrator::Q
protected

Definition at line 2466 of file bilininteg.hpp.

◆ quad1D

int mfem::VectorMassIntegrator::quad1D
protected

Definition at line 2473 of file bilininteg.hpp.

◆ VQ

VectorCoefficient* mfem::VectorMassIntegrator::VQ
protected

Definition at line 2467 of file bilininteg.hpp.


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