MFEM  v4.5.2
Finite element discretization library
Public Member Functions | Protected Member Functions | List of all members
mfem::MixedVectorGradientIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 MixedVectorGradientIntegrator ()
 
 MixedVectorGradientIntegrator (Coefficient &q)
 
 MixedVectorGradientIntegrator (DiagonalMatrixCoefficient &dq)
 
 MixedVectorGradientIntegrator (MatrixCoefficient &mq)
 
- Public Member Functions inherited from mfem::MixedVectorIntegrator
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementMatrix (const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &elmat)
 Support for use in BilinearForm. Can be used only when appropriate. More...
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePAInteriorFaces (const FiniteElementSpace &fes)
 
virtual void AssemblePABoundaryFaces (const FiniteElementSpace &fes)
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleDiagonalPA_ADAt (const Vector &D, Vector &diag)
 Assemble diagonal of ADA^T (A is this integrator) and add it to diag. More...
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add=true)
 Method defining element assembly. More...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
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 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 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. More...
 
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. More...
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix. More...
 
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. More...
 
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. More...
 
virtual double ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL)
 Virtual method required for Zienkiewicz-Zhu type error estimators. More...
 
virtual ~BilinearFormIntegrator ()
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- 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). More...
 
void SetIntegrationRule (const IntegrationRule &ir)
 Prescribe a fixed IntegrationRule to use. More...
 
void SetPAMemoryType (MemoryType mt)
 
const IntegrationRuleGetIntegrationRule () const
 Get the integration rule of the integrator (possibly NULL). More...
 
virtual double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. More...
 
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. More...
 
virtual double GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly. More...
 
virtual void AddMultGradPA (const Vector &x, Vector &y) const
 Method for partially assembled gradient action. More...
 
virtual void AssembleGradDiagonalPA (Vector &diag) const
 Method for computing the diagonal of the gradient with partial assembly. More...
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend. More...
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Member Functions

virtual bool VerifyFiniteElementTypes (const FiniteElement &trial_fe, const FiniteElement &test_fe) const
 
virtual const char * FiniteElementTypeFailureMessage () const
 
virtual int GetTrialVDim (const FiniteElement &trial_fe)
 
virtual void CalcTrialShape (const FiniteElement &trial_fe, ElementTransformation &Trans, DenseMatrix &shape)
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AddMultPA (const Vector &, Vector &) const
 Method for partially assembled action. More...
 
virtual void AddMultTransposePA (const Vector &, Vector &) const
 Method for partially assembled transposed action. More...
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Protected Member Functions inherited from mfem::MixedVectorIntegrator
 MixedVectorIntegrator ()
 
 MixedVectorIntegrator (Coefficient &q)
 
 MixedVectorIntegrator (VectorCoefficient &vq, bool diag=true)
 
 MixedVectorIntegrator (MatrixCoefficient &mq)
 
virtual int GetIntegrationOrder (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
 
virtual int GetTestVDim (const FiniteElement &test_fe)
 
virtual void CalcTestShape (const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &shape)
 
- Protected Member Functions inherited from mfem::BilinearFormIntegrator
 BilinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Additional Inherited Members

- Protected Attributes inherited from mfem::MixedVectorIntegrator
bool same_calc_shape
 
int space_dim
 
CoefficientQ
 
VectorCoefficientVQ
 
DiagonalMatrixCoefficientDQ
 
MatrixCoefficientMQ
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Detailed Description

Class for integrating the bilinear form a(u,v) := (Q grad u, v) in either 2D or 3D and where Q is an optional coefficient (of type scalar, matrix, or diagonal matrix) u is in H1 and v is in H(Curl) or H(Div). Partial assembly (PA) is supported but could be further optimized by using more efficient threading and shared memory.

Definition at line 1808 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ MixedVectorGradientIntegrator() [1/4]

mfem::MixedVectorGradientIntegrator::MixedVectorGradientIntegrator ( )
inline

Definition at line 1811 of file bilininteg.hpp.

◆ MixedVectorGradientIntegrator() [2/4]

mfem::MixedVectorGradientIntegrator::MixedVectorGradientIntegrator ( Coefficient q)
inline

Definition at line 1812 of file bilininteg.hpp.

◆ MixedVectorGradientIntegrator() [3/4]

mfem::MixedVectorGradientIntegrator::MixedVectorGradientIntegrator ( DiagonalMatrixCoefficient dq)
inline

Definition at line 1814 of file bilininteg.hpp.

◆ MixedVectorGradientIntegrator() [4/4]

mfem::MixedVectorGradientIntegrator::MixedVectorGradientIntegrator ( MatrixCoefficient mq)
inline

Definition at line 1816 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultPA()

void mfem::MixedVectorGradientIntegrator::AddMultPA ( const Vector x,
Vector y 
) const
protectedvirtual

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 1115 of file bilininteg_vectorfe.cpp.

◆ AddMultTransposePA()

void mfem::MixedVectorGradientIntegrator::AddMultTransposePA ( const Vector x,
Vector y 
) const
protectedvirtual

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 1129 of file bilininteg_vectorfe.cpp.

◆ AssemblePA() [1/5]

void mfem::BilinearFormIntegrator::AssemblePA
protected

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().

Definition at line 23 of file bilininteg.cpp.

◆ AssemblePA() [2/5]

void mfem::BilinearFormIntegrator::AssemblePA
protected

Used with BilinearFormIntegrators that have different spaces.

Definition at line 29 of file bilininteg.cpp.

◆ AssemblePA() [3/5]

void mfem::NonlinearFormIntegrator::AssemblePA
protected

Method defining partial assembly.

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

Definition at line 25 of file nonlininteg.cpp.

◆ AssemblePA() [4/5]

void mfem::NonlinearFormIntegrator::AssemblePA
protected

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultPA(). Used with BilinearFormIntegrators that have different spaces.

Definition at line 31 of file nonlininteg.cpp.

◆ AssemblePA() [5/5]

void mfem::MixedVectorGradientIntegrator::AssemblePA ( const FiniteElementSpace trial_fes,
const FiniteElementSpace test_fes 
)
protectedvirtual

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1054 of file bilininteg_vectorfe.cpp.

◆ CalcTrialShape()

virtual void mfem::MixedVectorGradientIntegrator::CalcTrialShape ( const FiniteElement trial_fe,
ElementTransformation Trans,
DenseMatrix shape 
)
inlineprotectedvirtual

Reimplemented from mfem::MixedVectorIntegrator.

Definition at line 1838 of file bilininteg.hpp.

◆ FiniteElementTypeFailureMessage()

virtual const char* mfem::MixedVectorGradientIntegrator::FiniteElementTypeFailureMessage ( ) const
inlineprotectedvirtual

Reimplemented from mfem::MixedVectorIntegrator.

Definition at line 1828 of file bilininteg.hpp.

◆ GetTrialVDim()

virtual int mfem::MixedVectorGradientIntegrator::GetTrialVDim ( const FiniteElement trial_fe)
inlineprotectedvirtual

Reimplemented from mfem::MixedVectorIntegrator.

Definition at line 1835 of file bilininteg.hpp.

◆ VerifyFiniteElementTypes()

virtual bool mfem::MixedVectorGradientIntegrator::VerifyFiniteElementTypes ( const FiniteElement trial_fe,
const FiniteElement test_fe 
) const
inlineprotectedvirtual

Reimplemented from mfem::MixedVectorIntegrator.

Definition at line 1820 of file bilininteg.hpp.


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