MFEM  v4.6.0
Finite element discretization library
Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
mfem::DiffusionIntegrator Class Reference

#include <bilininteg.hpp>

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

Public Member Functions

 DiffusionIntegrator (const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with coefficient Q = 1. More...
 
 DiffusionIntegrator (Coefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a scalar coefficient q. More...
 
 DiffusionIntegrator (VectorCoefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a vector coefficient q. More...
 
 DiffusionIntegrator (MatrixCoefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a matrix coefficient q. More...
 
virtual void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat)
 
virtual void AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat)
 
virtual void AssembleNURBSPA (const FiniteElementSpace &fes)
 Method defining partial assembly on NURBS patches. More...
 
void AssemblePatchPA (const int patch, const FiniteElementSpace &fes)
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the BilinearFormIntegrator. 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 void AssembleMF (const FiniteElementSpace &fes)
 Method defining matrix-free assembly. More...
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add)
 Method defining element assembly. More...
 
virtual void AssembleDiagonalPA (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AssembleDiagonalMF (Vector &diag)
 Assemble diagonal and add it to Vector diag. More...
 
virtual void AddMultMF (const Vector &, Vector &) const
 
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 AddMultNURBSPA (const Vector &, Vector &) const
 Method for partially assembled action on NURBS patches. More...
 
void AddMultPatchPA (const int patch, const Vector &x, Vector &y) const
 
bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend. More...
 
CoefficientGetCoefficient () const
 
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)
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
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 ADA^T (A is this integrator) and add it to diag. More...
 
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 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 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 ~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 SetIntegrationMode (Mode m)
 
void SetNURBSPatchIntRule (NURBSMeshRules *pr)
 For patchwise integration, SetNURBSPatchIntRule must be called. More...
 
bool HasNURBSPatchIntRule () const
 
bool Patchwise () const
 
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...
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Static Public Member Functions

static const IntegrationRuleGetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe)
 

Protected Attributes

CoefficientQ
 
VectorCoefficientVQ
 
MatrixCoefficientMQ
 
- 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 grad u, grad v) where Q can be a scalar or a matrix coefficient.

Definition at line 2095 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ DiffusionIntegrator() [1/4]

mfem::DiffusionIntegrator::DiffusionIntegrator ( const IntegrationRule ir = nullptr)
inline

Construct a diffusion integrator with coefficient Q = 1.

Definition at line 2173 of file bilininteg.hpp.

◆ DiffusionIntegrator() [2/4]

mfem::DiffusionIntegrator::DiffusionIntegrator ( Coefficient q,
const IntegrationRule ir = nullptr 
)
inline

Construct a diffusion integrator with a scalar coefficient q.

Definition at line 2178 of file bilininteg.hpp.

◆ DiffusionIntegrator() [3/4]

mfem::DiffusionIntegrator::DiffusionIntegrator ( VectorCoefficient q,
const IntegrationRule ir = nullptr 
)
inline

Construct a diffusion integrator with a vector coefficient q.

Definition at line 2183 of file bilininteg.hpp.

◆ DiffusionIntegrator() [4/4]

mfem::DiffusionIntegrator::DiffusionIntegrator ( MatrixCoefficient q,
const IntegrationRule ir = nullptr 
)
inline

Construct a diffusion integrator with a matrix coefficient q.

Definition at line 2189 of file bilininteg.hpp.

Member Function Documentation

◆ AddMultMF()

virtual void mfem::DiffusionIntegrator::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.

◆ AddMultNURBSPA()

virtual void mfem::DiffusionIntegrator::AddMultNURBSPA ( const Vector x,
Vector y 
) const
virtual

Method for partially assembled action on NURBS patches.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AddMultPA()

virtual void mfem::DiffusionIntegrator::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.

◆ AddMultPatchPA()

void mfem::DiffusionIntegrator::AddMultPatchPA ( const int  patch,
const Vector x,
Vector y 
) const

◆ AddMultTransposePA()

virtual void mfem::DiffusionIntegrator::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 from mfem::BilinearFormIntegrator.

◆ AssembleDiagonalMF()

virtual void mfem::DiffusionIntegrator::AssembleDiagonalMF ( Vector diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleDiagonalPA()

virtual void mfem::DiffusionIntegrator::AssembleDiagonalPA ( Vector diag)
virtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssembleEA()

virtual void mfem::DiffusionIntegrator::AssembleEA ( const FiniteElementSpace fes,
Vector emat,
const bool  add 
)
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 from mfem::BilinearFormIntegrator.

◆ AssembleElementMatrix()

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

Given a particular Finite Element computes the element stiffness matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 832 of file bilininteg.cpp.

◆ AssembleElementMatrix2()

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

Given a trial and test Finite Element computes the element stiffness matrix elmat.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 924 of file bilininteg.cpp.

◆ AssembleElementVector()

void mfem::DiffusionIntegrator::AssembleElementVector ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun,
Vector elvect 
)
virtual

Perform the local action of the BilinearFormIntegrator.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1008 of file bilininteg.cpp.

◆ AssembleMF()

virtual void mfem::DiffusionIntegrator::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.

◆ AssembleNURBSPA()

virtual void mfem::DiffusionIntegrator::AssembleNURBSPA ( const FiniteElementSpace fes)
virtual

Method defining partial assembly on NURBS patches.

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

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssemblePA() [1/5]

void mfem::BilinearFormIntegrator::AssemblePA

Used with BilinearFormIntegrators that have different spaces.

Definition at line 35 of file bilininteg.cpp.

◆ AssemblePA() [2/5]

void mfem::NonlinearFormIntegrator::AssemblePA

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() [3/5]

void mfem::BilinearFormIntegrator::AssemblePA

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() [4/5]

void mfem::NonlinearFormIntegrator::AssemblePA

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() [5/5]

virtual void mfem::DiffusionIntegrator::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.

◆ AssemblePatchMatrix()

virtual void mfem::DiffusionIntegrator::AssemblePatchMatrix ( const int  patch,
const FiniteElementSpace fes,
SparseMatrix *&  smat 
)
virtual

Given a particular NURBS patch, computes the patch matrix as a SparseMatrix smat.

Reimplemented from mfem::BilinearFormIntegrator.

◆ AssemblePatchPA()

void mfem::DiffusionIntegrator::AssemblePatchPA ( const int  patch,
const FiniteElementSpace fes 
)

◆ ComputeElementFlux()

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

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

Definition at line 1091 of file bilininteg.cpp.

◆ ComputeFluxEnergy()

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

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

Definition at line 1190 of file bilininteg.cpp.

◆ GetCoefficient()

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

Definition at line 2257 of file bilininteg.hpp.

◆ GetRule()

const IntegrationRule & mfem::DiffusionIntegrator::GetRule ( const FiniteElement trial_fe,
const FiniteElement test_fe 
)
static

Definition at line 1265 of file bilininteg.cpp.

◆ SupportsCeed()

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

Indicates whether this integrator can use a Ceed backend.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2255 of file bilininteg.hpp.

Member Data Documentation

◆ MQ

MatrixCoefficient* mfem::DiffusionIntegrator::MQ
protected

Definition at line 2100 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::DiffusionIntegrator::Q
protected

Definition at line 2098 of file bilininteg.hpp.

◆ VQ

VectorCoefficient* mfem::DiffusionIntegrator::VQ
protected

Definition at line 2099 of file bilininteg.hpp.


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