MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::VectorFEMassIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 VectorFEMassIntegrator ()
 
 VectorFEMassIntegrator (Coefficient *q_)
 
 VectorFEMassIntegrator (Coefficient &q)
 
 VectorFEMassIntegrator (DiagonalMatrixCoefficient *dq_)
 
 VectorFEMassIntegrator (DiagonalMatrixCoefficient &dq)
 
 VectorFEMassIntegrator (MatrixCoefficient *mq_)
 
 VectorFEMassIntegrator (MatrixCoefficient &mq)
 
void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) override
 Given a particular Finite Element computes the element matrix elmat.
 
void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) override
 
void AssemblePA (const FiniteElementSpace &fes) override
 Method defining partial assembly.
 
void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
 
void AddMultPA (const Vector &x, Vector &y) const override
 Method for partially assembled action.
 
void AddMultTransposePA (const Vector &x, Vector &y) const override
 Method for partially assembled transposed action.
 
void AssembleDiagonalPA (Vector &diag) override
 Assemble diagonal and add it to Vector diag.
 
void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add) override
 Method defining element assembly.
 
const CoefficientGetCoefficient () const
 
- 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.
 
void AssembleMF (const FiniteElementSpace &fes) override
 Method defining matrix-free assembly.
 
void AddMultMF (const Vector &x, Vector &y) const override
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag.
 
virtual void AssembleEABoundary (const FiniteElementSpace &fes, Vector &ea_data_bdr, const bool add=true)
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &fes, Vector &ea_data_int, Vector &ea_data_ext, const bool add=true)
 
virtual void AssembleEAInteriorFaces (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes, Vector &emat, const bool add=true)
 Method defining element assembly for mixed trace integrators.
 
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_fe1, const FiniteElement &test_fe1, const FiniteElement &trial_fe2, const FiniteElement &test_fe2, 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)
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) override
 Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.
 
void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) override
 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.
 
void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat) override
 Assemble the local gradient matrix.
 
void AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) override
 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
void SetIntegrationMode (Mode m)
 
bool Patchwise () const
 
void SetPAMemoryType (MemoryType mt)
 
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.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 
- Public Member Functions inherited from mfem::Integrator
 Integrator (const IntegrationRule *ir=NULL)
 Create a new Integrator, optionally providing a prescribed quadrature rule to use in assembly.
 
virtual void SetIntRule (const IntegrationRule *ir)
 Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate rule.
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. Sets the NURBS patch integration rule to null.
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 Sets an integration rule for use on NURBS patches.
 
bool HasNURBSPatchIntRule () const
 Check if a NURBS patch integration rule has been set.
 
const IntegrationRuleGetIntRule () const
 Directly return the IntRule pointer (possibly null) without checking for NURBS patch rules or falling back on a default.
 
const IntegrationRuleGetIntegrationRule () const
 Equivalent to GetIntRule, but retained for backward compatibility with applications.
 

Protected Attributes

CoefficientQ
 
DiagonalMatrixCoefficientDQ
 
MatrixCoefficientMQ
 
Vector pa_data
 
const DofToQuadmapsO
 Not owned. DOF-to-quad map, open.
 
const DofToQuadmapsC
 Not owned. DOF-to-quad map, closed.
 
const DofToQuadmapsOtest
 Not owned. DOF-to-quad map, open.
 
const DofToQuadmapsCtest
 Not owned. DOF-to-quad map, closed.
 
const GeometricFactorsgeom
 Not owned.
 
int dim
 
int ne
 
int nq
 
int dofs1D
 
int dofs1Dtest
 
int quad1D
 
int trial_fetype
 
int test_fetype
 
bool symmetric = true
 False if using a nonsymmetric matrix coefficient.
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
Mode integrationMode = Mode::ELEMENTWISE
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 
- Protected Attributes inherited from mfem::Integrator
const IntegrationRuleIntRule
 
NURBSMeshRulespatchRules = nullptr
 

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)
 
- Protected Member Functions inherited from mfem::Integrator
const IntegrationRuleGetIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state of the Integrator object.
 
const IntegrationRuleGetIntegrationRule (const FiniteElement &el, const ElementTransformation &trans) const
 Returns an integration rule based on the arguments and internal state. (Version for identical trial_fe and test_fe)
 
virtual const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const
 Subclasses should override to choose a default integration rule.
 

Detailed Description

Integrator for \((Q u, v)\), where \(Q\) is an optional coefficient (of type scalar, vector (diagonal matrix), or matrix), trial function \(u\) is in \(H(curl\) or \(H(div)\), and test function \(v\) is in \(H(curl\), \(H(div)\), or \(v=(v_1,\dots,v_n)\), where \(v_i\) are in \(H^1\).

Definition at line 2885 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ VectorFEMassIntegrator() [1/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( )
inline

Definition at line 2916 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [2/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( Coefficient * q_)
inline

Definition at line 2917 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [3/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( Coefficient & q)
inline

Definition at line 2918 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [4/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( DiagonalMatrixCoefficient * dq_)
inline

Definition at line 2919 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [5/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( DiagonalMatrixCoefficient & dq)
inline

Definition at line 2920 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [6/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( MatrixCoefficient * mq_)
inline

Definition at line 2921 of file bilininteg.hpp.

◆ VectorFEMassIntegrator() [7/7]

mfem::VectorFEMassIntegrator::VectorFEMassIntegrator ( MatrixCoefficient & mq)
inline

Definition at line 2922 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultPA()

void mfem::VectorFEMassIntegrator::AddMultPA ( const Vector & x,
Vector & y ) const
overridevirtual

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.

Definition at line 216 of file bilininteg_vectorfemass_pa.cpp.

◆ AddMultTransposePA()

void mfem::VectorFEMassIntegrator::AddMultTransposePA ( const Vector & x,
Vector & y ) const
overridevirtual

Method for partially assembled transposed action.

Perform the transpose 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.

Definition at line 316 of file bilininteg_vectorfemass_pa.cpp.

◆ AssembleDiagonalPA()

void mfem::VectorFEMassIntegrator::AssembleDiagonalPA ( Vector & diag)
overridevirtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 146 of file bilininteg_vectorfemass_pa.cpp.

◆ AssembleEA()

void mfem::VectorFEMassIntegrator::AssembleEA ( const FiniteElementSpace & fes,
Vector & emat,
const bool add )
overridevirtual

Method defining element assembly.

The result of the element assembly is added to the emat Vector if add is true. Otherwise, if add is false, we set emat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 256 of file bilininteg_hdiv_ea.cpp.

◆ AssembleElementMatrix()

void mfem::VectorFEMassIntegrator::AssembleElementMatrix ( const FiniteElement & el,
ElementTransformation & Trans,
DenseMatrix & elmat )
overridevirtual

Given a particular Finite Element computes the element matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 2602 of file bilininteg.cpp.

◆ AssembleElementMatrix2()

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

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 2668 of file bilininteg.cpp.

◆ AssemblePA() [1/2]

void mfem::VectorFEMassIntegrator::AssemblePA ( const FiniteElementSpace & fes)
overridevirtual

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.

Definition at line 23 of file bilininteg_vectorfemass_pa.cpp.

◆ AssemblePA() [2/2]

void mfem::VectorFEMassIntegrator::AssemblePA ( const FiniteElementSpace & trial_fes,
const FiniteElementSpace & test_fes )
overridevirtual

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 28 of file bilininteg_vectorfemass_pa.cpp.

◆ GetCoefficient()

const Coefficient * mfem::VectorFEMassIntegrator::GetCoefficient ( ) const
inline

Definition at line 2941 of file bilininteg.hpp.

Member Data Documentation

◆ dim

int mfem::VectorFEMassIntegrator::dim
protected

Definition at line 2912 of file bilininteg.hpp.

◆ dofs1D

int mfem::VectorFEMassIntegrator::dofs1D
protected

Definition at line 2912 of file bilininteg.hpp.

◆ dofs1Dtest

int mfem::VectorFEMassIntegrator::dofs1Dtest
protected

Definition at line 2912 of file bilininteg.hpp.

◆ DQ

DiagonalMatrixCoefficient* mfem::VectorFEMassIntegrator::DQ
protected

Definition at line 2902 of file bilininteg.hpp.

◆ geom

const GeometricFactors* mfem::VectorFEMassIntegrator::geom
protected

Not owned.

Definition at line 2911 of file bilininteg.hpp.

◆ mapsC

const DofToQuad* mfem::VectorFEMassIntegrator::mapsC
protected

Not owned. DOF-to-quad map, closed.

Definition at line 2908 of file bilininteg.hpp.

◆ mapsCtest

const DofToQuad* mfem::VectorFEMassIntegrator::mapsCtest
protected

Not owned. DOF-to-quad map, closed.

Definition at line 2910 of file bilininteg.hpp.

◆ mapsO

const DofToQuad* mfem::VectorFEMassIntegrator::mapsO
protected

Not owned. DOF-to-quad map, open.

Definition at line 2907 of file bilininteg.hpp.

◆ mapsOtest

const DofToQuad* mfem::VectorFEMassIntegrator::mapsOtest
protected

Not owned. DOF-to-quad map, open.

Definition at line 2909 of file bilininteg.hpp.

◆ MQ

MatrixCoefficient* mfem::VectorFEMassIntegrator::MQ
protected

Definition at line 2903 of file bilininteg.hpp.

◆ ne

int mfem::VectorFEMassIntegrator::ne
protected

Definition at line 2912 of file bilininteg.hpp.

◆ nq

int mfem::VectorFEMassIntegrator::nq
protected

Definition at line 2912 of file bilininteg.hpp.

◆ pa_data

Vector mfem::VectorFEMassIntegrator::pa_data
protected

Definition at line 2906 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::VectorFEMassIntegrator::Q
protected

Definition at line 2901 of file bilininteg.hpp.

◆ quad1D

int mfem::VectorFEMassIntegrator::quad1D
protected

Definition at line 2912 of file bilininteg.hpp.

◆ symmetric

bool mfem::VectorFEMassIntegrator::symmetric = true
protected

False if using a nonsymmetric matrix coefficient.

Definition at line 2913 of file bilininteg.hpp.

◆ test_fetype

int mfem::VectorFEMassIntegrator::test_fetype
protected

Definition at line 2912 of file bilininteg.hpp.

◆ trial_fetype

int mfem::VectorFEMassIntegrator::trial_fetype
protected

Definition at line 2912 of file bilininteg.hpp.


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