MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | List of all members
mfem::BilinearFormIntegrator Class Reference

Abstract base class BilinearFormIntegrator. More...

#include <bilininteg.hpp>

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

Public Member Functions

virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
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 AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. More...
 
virtual void AddMultTransposePA (const Vector &x, Vector &y) const
 Method for partially assembled transposed action. 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 AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 Given a particular Finite Element computes the element matrix elmat. More...
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
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 ()
 
- 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

 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::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Detailed Description

Abstract base class BilinearFormIntegrator.

Definition at line 35 of file bilininteg.hpp.

Constructor & Destructor Documentation

mfem::BilinearFormIntegrator::BilinearFormIntegrator ( const IntegrationRule ir = NULL)
inlineprotected

Definition at line 38 of file bilininteg.hpp.

virtual mfem::BilinearFormIntegrator::~BilinearFormIntegrator ( )
inlinevirtual

Definition at line 254 of file bilininteg.hpp.

Member Function Documentation

void mfem::BilinearFormIntegrator::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::NonlinearFormIntegrator.

Reimplemented in mfem::VectorDiffusionIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, and mfem::SumIntegrator.

Definition at line 105 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::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::NonlinearFormIntegrator.

Reimplemented in mfem::IdentityInterpolator, mfem::GradientInterpolator, mfem::DGTraceIntegrator, mfem::VectorDiffusionIntegrator, mfem::DivDivIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEMassIntegrator, mfem::CurlCurlIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::GradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::MixedScalarCurlIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 87 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AddMultTransposeMF ( const Vector x,
Vector y 
) const
virtual

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 in mfem::SumIntegrator.

Definition at line 111 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AddMultTransposePA ( const Vector x,
Vector y 
) const
virtual

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 in mfem::IdentityInterpolator, mfem::GradientInterpolator, mfem::DGTraceIntegrator, mfem::VectorDivergenceIntegrator, mfem::VectorFEMassIntegrator, mfem::VectorFEDivergenceIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::GradientIntegrator, mfem::MixedVectorWeakCurlIntegrator, mfem::MixedVectorCurlIntegrator, mfem::MixedVectorGradientIntegrator, mfem::MixedScalarCurlIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 93 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleDiagonalMF ( Vector diag)
virtual
void mfem::BilinearFormIntegrator::AssembleDiagonalPA ( Vector diag)
virtual
void mfem::BilinearFormIntegrator::AssembleDiagonalPA_ADAt ( const Vector D,
Vector diag 
)
virtual

Assemble diagonal of ADA^T (A is this integrator) and add it to diag.

Reimplemented in mfem::VectorFEDivergenceIntegrator.

Definition at line 81 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleEA ( const FiniteElementSpace fes,
Vector emat,
const bool  add = true 
)
virtual

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 in mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 54 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleEABoundaryFaces ( const FiniteElementSpace fes,
Vector ea_data_bdr,
const bool  add = true 
)
virtual

Reimplemented in mfem::DGTraceIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 72 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleEAInteriorFaces ( const FiniteElementSpace fes,
Vector ea_data_int,
Vector ea_data_ext,
const bool  add = true 
)
virtual

Reimplemented in mfem::DGTraceIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 62 of file bilininteg.cpp.

virtual void mfem::BilinearFormIntegrator::AssembleElementGrad ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
DenseMatrix elmat 
)
inlinevirtual

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 178 of file bilininteg.hpp.

void mfem::BilinearFormIntegrator::AssembleElementMatrix ( const FiniteElement el,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual
void mfem::BilinearFormIntegrator::AssembleElementMatrix2 ( const FiniteElement trial_fe,
const FiniteElement test_fe,
ElementTransformation Trans,
DenseMatrix elmat 
)
virtual
void mfem::BilinearFormIntegrator::AssembleElementVector ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
Vector elvect 
)
virtual

Perform the local action of the BilinearFormIntegrator. Note that the default implementation in the base class is general but not efficient.

Reimplemented from mfem::NonlinearFormIntegrator.

Reimplemented in mfem::VectorDiffusionIntegrator, and mfem::DiffusionIntegrator.

Definition at line 156 of file bilininteg.cpp.

virtual void mfem::BilinearFormIntegrator::AssembleFaceGrad ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Tr,
const Vector elfun,
DenseMatrix elmat 
)
inlinevirtual

Assemble the local action of the gradient of the NonlinearFormIntegrator resulting from a face integral term.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 183 of file bilininteg.hpp.

void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual
void mfem::BilinearFormIntegrator::AssembleFaceMatrix ( const FiniteElement trial_face_fe,
const FiniteElement test_fe1,
const FiniteElement test_fe2,
FaceElementTransformations Trans,
DenseMatrix elmat 
)
virtual

Abstract method used for assembling TraceFaceIntegrators in a MixedBilinearForm.

Reimplemented in mfem::NormalTraceJumpIntegrator, mfem::TraceJumpIntegrator, and mfem::SumIntegrator.

Definition at line 147 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssembleFaceVector ( const FiniteElement el1,
const FiniteElement el2,
FaceElementTransformations Tr,
const Vector elfun,
Vector elvect 
)
virtual

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.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 167 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::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::NonlinearFormIntegrator.

Reimplemented in mfem::VectorDiffusionIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, and mfem::SumIntegrator.

Definition at line 99 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::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::NonlinearFormIntegrator.

Reimplemented in mfem::VectorDiffusionIntegrator, mfem::DivDivIntegrator, mfem::VectorFEMassIntegrator, mfem::CurlCurlIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 23 of file bilininteg.cpp.

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

Reimplemented in mfem::DGTraceIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 42 of file bilininteg.cpp.

void mfem::BilinearFormIntegrator::AssemblePAInteriorFaces ( const FiniteElementSpace fes)
virtual

Reimplemented in mfem::DGTraceIntegrator, mfem::SumIntegrator, and mfem::TransposeIntegrator.

Definition at line 36 of file bilininteg.cpp.

virtual void mfem::BilinearFormIntegrator::ComputeElementFlux ( const FiniteElement el,
ElementTransformation Trans,
Vector u,
const FiniteElement fluxelem,
Vector flux,
bool  with_coef = true,
const IntegrationRule ir = NULL 
)
inlinevirtual

Virtual method required for Zienkiewicz-Zhu type error estimators.

The purpose of the method is to compute a local "flux" finite element function given a local finite element solution. The "flux" function has to be computed in terms of its coefficients (represented by the Vector flux) which multiply the basis functions defined by the FiniteElement fluxelem. Typically, the "flux" function will have more than one component and consequently flux should be store the coefficients of all components: first all coefficient for component 0, then all coefficients for component 1, etc. What the "flux" function represents depends on the specific integrator. For example, in the case of DiffusionIntegrator, the flux is the gradient of the solution multiplied by the diffusion coefficient.

Parameters
[in]elFiniteElement of the solution.
[in]TransThe ElementTransformation describing the physical position of the mesh element.
[in]uSolution coefficients representing the expansion of the solution function in the basis of el.
[in]fluxelemFiniteElement of the "flux".
[out]flux"Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem. The size of flux as a Vector has to be set by this method, e.g. using Vector::SetSize().
[in]with_coefIf zero (the default value is 1) the implementation of the method may choose not to scale the "flux" function by any coefficients describing the integrator.
[in]irIf passed (the default value is NULL), the implementation of the method will ignore the integration rule provided by the fluxelem parameter and, instead, compute the discrete flux at the points specified by the integration rule ir.

Reimplemented in mfem::ElasticityIntegrator, mfem::CurlCurlIntegrator, and mfem::DiffusionIntegrator.

Definition at line 223 of file bilininteg.hpp.

virtual double mfem::BilinearFormIntegrator::ComputeFluxEnergy ( const FiniteElement fluxelem,
ElementTransformation Trans,
Vector flux,
Vector d_energy = NULL 
)
inlinevirtual

Virtual method required for Zienkiewicz-Zhu type error estimators.

The purpose of this method is to compute a local number that measures the energy of a given "flux" function (see ComputeElementFlux() for a description of the "flux" function). Typically, the energy of a "flux" function should be equal to a_local(u,u), if the "flux" is defined from a solution u; here a_local(.,.) denotes the element-local bilinear form represented by the integrator.

Parameters
[in]fluxelemFiniteElement of the "flux".
[in]TransThe ElementTransformation describing the physical position of the mesh element.
[in]flux"Flux" coefficients representing the expansion of the "flux" function in the basis of fluxelem.
[out]d_energyIf not NULL, the given Vector should be set to represent directional energy split that can be used for anisotropic error estimation.
Returns
The computed energy.

Reimplemented in mfem::ElasticityIntegrator, mfem::CurlCurlIntegrator, and mfem::DiffusionIntegrator.

Definition at line 249 of file bilininteg.hpp.


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