MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::DiffusionIntegrator Class Reference

#include <bilininteg.hpp>

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

Classes

struct  Kernels
 

Public Types

using ApplyKernelType
 
using DiagonalKernelType
 
- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 

Public Member Functions

 MFEM_REGISTER_KERNELS (ApplyPAKernels, ApplyKernelType,(int, int, int))
 
 MFEM_REGISTER_KERNELS (DiagonalPAKernels, DiagonalKernelType,(int, int, int))
 
 DiffusionIntegrator (const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with coefficient Q = 1.
 
 DiffusionIntegrator (Coefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a scalar coefficient q.
 
 DiffusionIntegrator (VectorCoefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a vector coefficient q.
 
 DiffusionIntegrator (MatrixCoefficient &q, const IntegrationRule *ir=nullptr)
 Construct a diffusion integrator with a matrix coefficient q.
 
void AssembleElementMatrix (const FiniteElement &el, ElementTransformation &Trans, DenseMatrix &elmat) override
 
void AssembleElementMatrix2 (const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans, DenseMatrix &elmat) override
 
void AssemblePatchMatrix (const int patch, const FiniteElementSpace &fes, SparseMatrix *&smat) override
 
void AssembleNURBSPA (const FiniteElementSpace &fes) override
 Method defining partial assembly on NURBS patches.
 
void AssemblePatchPA (const int patch, const FiniteElementSpace &fes)
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) override
 Perform the local action of the BilinearFormIntegrator.
 
void ComputeElementFlux (const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true, const IntegrationRule *ir=NULL) override
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
real_t ComputeFluxEnergy (const FiniteElement &fluxelem, ElementTransformation &Trans, Vector &flux, Vector *d_energy=NULL) override
 Virtual method required for Zienkiewicz-Zhu type error estimators.
 
void AssembleMF (const FiniteElementSpace &fes) override
 Method defining matrix-free assembly.
 
void AssemblePA (const FiniteElementSpace &fes) override
 Method defining partial assembly.
 
void AssembleEA (const FiniteElementSpace &fes, Vector &emat, const bool add) override
 Method defining element assembly.
 
void AssembleDiagonalPA (Vector &diag) override
 Assemble diagonal and add it to Vector diag.
 
void AssembleDiagonalMF (Vector &diag) override
 Assemble diagonal and add it to Vector diag.
 
void AddMultMF (const Vector &, Vector &) const override
 
void AddMultPA (const Vector &, Vector &) const override
 Method for partially assembled action.
 
void AddMultTransposePA (const Vector &, Vector &) const override
 Method for partially assembled transposed action.
 
void AddMultNURBSPA (const Vector &, Vector &) const override
 Method for partially assembled action on NURBS patches.
 
void AddMultPatchPA (const int patch, const Vector &x, Vector &y) const
 
bool SupportsCeed () const override
 Indicates whether this integrator can use a Ceed backend.
 
CoefficientGetCoefficient () const
 
void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes) override
 
- Public Member Functions inherited from mfem::BilinearFormIntegrator
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 ADAT ( A is this integrator) and add it to diag.
 
virtual void AddMultTransposeMF (const Vector &x, Vector &y) const
 
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 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 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 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.
 
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.
 

Static Public Member Functions

static const IntegrationRuleGetRule (const FiniteElement &trial_fe, const FiniteElement &test_fe)
 
template<int DIM, int D1D, int Q1D>
static void AddSpecialization ()
 

Protected Member Functions

const IntegrationRuleGetDefaultIntegrationRule (const FiniteElement &trial_fe, const FiniteElement &test_fe, const ElementTransformation &trans) const override
 Subclasses should override to choose a default integration rule.
 
- 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)
 

Protected Attributes

CoefficientQ
 
VectorCoefficientVQ
 
MatrixCoefficientMQ
 
- 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
 

Detailed Description

Class for integrating the bilinear form a(u,v):=(Qu,v) where Q can be a scalar or a matrix coefficient.

Definition at line 2167 of file bilininteg.hpp.

Member Typedef Documentation

◆ ApplyKernelType

Initial value:
void(*)(const int, const bool, const Array<real_t>&,
const Array<real_t>&, const Array<real_t>&,
const Array<real_t>&,
const Vector&, const Vector&,
Vector&, const int, const int)

Definition at line 2171 of file bilininteg.hpp.

◆ DiagonalKernelType

Initial value:
void(*)(const int, const bool, const Array<real_t>&,
const Array<real_t>&, const Vector&, Vector&,
const int, const int)

Definition at line 2177 of file bilininteg.hpp.

Constructor & Destructor Documentation

◆ DiffusionIntegrator() [1/4]

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

Construct a diffusion integrator with coefficient Q = 1.

Definition at line 877 of file bilininteg.cpp.

◆ DiffusionIntegrator() [2/4]

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

Construct a diffusion integrator with a scalar coefficient q.

Definition at line 884 of file bilininteg.cpp.

◆ DiffusionIntegrator() [3/4]

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

Construct a diffusion integrator with a vector coefficient q.

Definition at line 891 of file bilininteg.cpp.

◆ DiffusionIntegrator() [4/4]

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

Construct a diffusion integrator with a matrix coefficient q.

Definition at line 898 of file bilininteg.cpp.

Member Function Documentation

◆ AddMultMF()

void mfem::DiffusionIntegrator::AddMultMF ( const Vector & x,
Vector & y ) const
overridevirtual

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.

Definition at line 61 of file bilininteg_diffusion_mf.cpp.

◆ AddMultNURBSPA()

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

Method for partially assembled action on NURBS patches.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 365 of file bilininteg_diffusion_pa.cpp.

◆ AddMultPA()

void mfem::DiffusionIntegrator::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 40 of file bilininteg_diffusion_pa.cpp.

◆ AddMultPatchPA()

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

Definition at line 169 of file bilininteg_diffusion_pa.cpp.

◆ AddMultTransposePA()

void mfem::DiffusionIntegrator::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 76 of file bilininteg_diffusion_pa.cpp.

◆ AddSpecialization()

template<int DIM, int D1D, int Q1D>
static void mfem::DiffusionIntegrator::AddSpecialization ( )
inlinestatic

Definition at line 2337 of file bilininteg.hpp.

◆ AssembleDiagonalMF()

void mfem::DiffusionIntegrator::AssembleDiagonalMF ( Vector & diag)
overridevirtual

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 48 of file bilininteg_diffusion_mf.cpp.

◆ AssembleDiagonalPA()

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

Assemble diagonal and add it to Vector diag.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 22 of file bilininteg_diffusion_pa.cpp.

◆ AssembleEA()

void mfem::DiffusionIntegrator::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 245 of file bilininteg_diffusion_ea.cpp.

◆ AssembleElementMatrix()

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

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 905 of file bilininteg.cpp.

◆ AssembleElementMatrix2()

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

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 979 of file bilininteg.cpp.

◆ AssembleElementVector()

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

Perform the local action of the BilinearFormIntegrator.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1062 of file bilininteg.cpp.

◆ AssembleMF()

void mfem::DiffusionIntegrator::AssembleMF ( const FiniteElementSpace & fes)
overridevirtual

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.

Definition at line 19 of file bilininteg_diffusion_mf.cpp.

◆ AssembleNURBSPA()

void mfem::DiffusionIntegrator::AssembleNURBSPA ( const FiniteElementSpace & fes)
overridevirtual

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.

Definition at line 144 of file bilininteg_diffusion_pa.cpp.

◆ AssemblePA() [1/2]

void mfem::DiffusionIntegrator::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 89 of file bilininteg_diffusion_pa.cpp.

◆ AssemblePA() [2/2]

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

Used with BilinearFormIntegrators that have different spaces.

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 52 of file bilininteg.cpp.

◆ AssemblePatchMatrix()

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

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

Reimplemented from mfem::BilinearFormIntegrator.

Definition at line 1234 of file bilininteg_diffusion_patch.cpp.

◆ AssemblePatchPA()

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

Definition at line 158 of file bilininteg_diffusion_pa.cpp.

◆ 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 )
overridevirtual

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

◆ ComputeFluxEnergy()

real_t mfem::DiffusionIntegrator::ComputeFluxEnergy ( const FiniteElement & fluxelem,
ElementTransformation & Trans,
Vector & flux,
Vector * d_energy = NULL )
overridevirtual

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

◆ GetCoefficient()

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

Definition at line 2334 of file bilininteg.hpp.

◆ GetDefaultIntegrationRule()

const IntegrationRule * mfem::DiffusionIntegrator::GetDefaultIntegrationRule ( const FiniteElement & trial_fe,
const FiniteElement & test_fe,
const ElementTransformation & trans ) const
inlineoverrideprotectedvirtual

Subclasses should override to choose a default integration rule.

This method is intended to be overridden by subclasses to choose an appropriate integration rule based on the finite element spaces and/or element transformation. The trial_fe and test_fe should be equal for linear forms. The default base-class implementation returns null, which assumes that an appropriate rule is provided by another means, or that null integration rules are handled appropriately by the caller.

Reimplemented from mfem::Integrator.

Definition at line 2343 of file bilininteg.hpp.

◆ GetRule()

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

Definition at line 1318 of file bilininteg.cpp.

◆ MFEM_REGISTER_KERNELS() [1/2]

mfem::DiffusionIntegrator::MFEM_REGISTER_KERNELS ( ApplyPAKernels ,
ApplyKernelType ,
(int, int, int)  )

◆ MFEM_REGISTER_KERNELS() [2/2]

mfem::DiffusionIntegrator::MFEM_REGISTER_KERNELS ( DiagonalPAKernels ,
DiagonalKernelType ,
(int, int, int)  )

◆ SupportsCeed()

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

Indicates whether this integrator can use a Ceed backend.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 2332 of file bilininteg.hpp.

Member Data Documentation

◆ MQ

MatrixCoefficient* mfem::DiffusionIntegrator::MQ
protected

Definition at line 2188 of file bilininteg.hpp.

◆ Q

Coefficient* mfem::DiffusionIntegrator::Q
protected

Definition at line 2186 of file bilininteg.hpp.

◆ VQ

VectorCoefficient* mfem::DiffusionIntegrator::VQ
protected

Definition at line 2187 of file bilininteg.hpp.


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