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

This class is used to express the local action of a general nonlinear finite element operator. In addition it may provide the capability to assemble the local gradient operator and to compute the local energy. More...

#include <nonlininteg.hpp>

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

Public Types

enum  Mode { ELEMENTWISE = 0, PATCHWISE = 1, PATCHWISE_REDUCED = 2 }
 

Public Member Functions

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 void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator. More...
 
virtual void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator resulting from a face integral term. 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 double GetElementEnergy (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun)
 Compute the local energy. 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 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 AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. 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...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
ceed::OperatorGetCeedOp ()
 
virtual ~NonlinearFormIntegrator ()
 

Protected Member Functions

 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 

Protected Attributes

const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Detailed Description

This class is used to express the local action of a general nonlinear finite element operator. In addition it may provide the capability to assemble the local gradient operator and to compute the local energy.

Definition at line 27 of file nonlininteg.hpp.

Member Enumeration Documentation

◆ Mode

Enumerator
ELEMENTWISE 

Element-wise integration (default)

PATCHWISE 

Patch-wise integration (NURBS meshes)

PATCHWISE_REDUCED 

Patch-wise integration (NURBS meshes) with reduced integration rules.

Definition at line 30 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ NonlinearFormIntegrator()

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

Definition at line 51 of file nonlininteg.hpp.

◆ ~NonlinearFormIntegrator()

virtual mfem::NonlinearFormIntegrator::~NonlinearFormIntegrator ( )
inlinevirtual

Definition at line 169 of file nonlininteg.hpp.

Member Function Documentation

◆ AddMultGradPA()

void mfem::NonlinearFormIntegrator::AddMultGradPA ( const Vector x,
Vector y 
) const
virtual

Method for partially assembled gradient action.

All arguments are E-vectors. This method can be called only after the method AssembleGradPA() has been called.

Parameters
[in]xThe gradient Operator is applied to the Vector x.
[in,out]yThe result Vector: \( y += G x \).

Reimplemented in mfem::TMOPComboIntegrator, and mfem::TMOP_Integrator.

Definition at line 51 of file nonlininteg.cpp.

◆ AddMultMF()

void mfem::NonlinearFormIntegrator::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 in mfem::VectorDiffusionIntegrator, mfem::VectorMassIntegrator, mfem::ConvectionIntegrator, mfem::MassIntegrator, mfem::DiffusionIntegrator, mfem::SumIntegrator, mfem::VectorConvectionNLFIntegrator, and mfem::BilinearFormIntegrator.

Definition at line 69 of file nonlininteg.cpp.

◆ AddMultPA()

void mfem::NonlinearFormIntegrator::AddMultPA ( const Vector x,
Vector y 
) const
virtual

◆ AssembleElementGrad()

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

◆ AssembleElementVector()

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

◆ AssembleFaceGrad()

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

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

Reimplemented in mfem::BilinearFormIntegrator.

Definition at line 99 of file nonlininteg.cpp.

◆ AssembleFaceVector()

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

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

Reimplemented in mfem::BilinearFormIntegrator, and FaceIntegrator.

Definition at line 83 of file nonlininteg.cpp.

◆ AssembleGradDiagonalPA()

void mfem::NonlinearFormIntegrator::AssembleGradDiagonalPA ( Vector diag) const
virtual

Method for computing the diagonal of the gradient with partial assembly.

The result Vector diag is an E-Vector. This method can be called only after the method AssembleGradPA() has been called.

Parameters
[in,out]diagThe result Vector: \( diag += diag(G) \).

Reimplemented in mfem::TMOPComboIntegrator, and mfem::TMOP_Integrator.

Definition at line 57 of file nonlininteg.cpp.

◆ AssembleGradPA()

void mfem::NonlinearFormIntegrator::AssembleGradPA ( const Vector x,
const FiniteElementSpace fes 
)
virtual

Prepare the integrator for partial assembly (PA) gradient evaluations on the given FE space fes at the state x.

The result of the partial assembly is stored internally so that it can be used later in the methods AddMultGradPA() and AssembleGradDiagonalPA(). The state Vector x is an E-vector.

Reimplemented in mfem::TMOPComboIntegrator, and mfem::TMOP_Integrator.

Definition at line 38 of file nonlininteg.cpp.

◆ AssembleMF()

void mfem::NonlinearFormIntegrator::AssembleMF ( const FiniteElementSpace fes)
virtual

◆ AssemblePA() [1/2]

void mfem::NonlinearFormIntegrator::AssemblePA ( const FiniteElementSpace fes)
virtual

◆ AssemblePA() [2/2]

void mfem::NonlinearFormIntegrator::AssemblePA ( const FiniteElementSpace trial_fes,
const FiniteElementSpace test_fes 
)
virtual

◆ GetCeedOp()

ceed::Operator& mfem::NonlinearFormIntegrator::GetCeedOp ( )
inline

Definition at line 167 of file nonlininteg.hpp.

◆ GetElementEnergy()

double mfem::NonlinearFormIntegrator::GetElementEnergy ( const FiniteElement el,
ElementTransformation Tr,
const Vector elfun 
)
virtual

◆ GetIntegrationRule()

const IntegrationRule* mfem::NonlinearFormIntegrator::GetIntegrationRule ( ) const
inline

Get the integration rule of the integrator (possibly NULL).

Definition at line 75 of file nonlininteg.hpp.

◆ GetLocalStateEnergyPA()

double mfem::NonlinearFormIntegrator::GetLocalStateEnergyPA ( const Vector x) const
virtual

Compute the local (to the MPI rank) energy with partial assembly.

Here the state x is an E-vector. This method can be called only after the method AssemblePA() has been called.

Reimplemented in mfem::TMOPComboIntegrator, and mfem::TMOP_Integrator.

Definition at line 18 of file nonlininteg.cpp.

◆ HasNURBSPatchIntRule()

bool mfem::NonlinearFormIntegrator::HasNURBSPatchIntRule ( ) const
inline

Definition at line 63 of file nonlininteg.hpp.

◆ Patchwise()

bool mfem::NonlinearFormIntegrator::Patchwise ( ) const
inline

Definition at line 65 of file nonlininteg.hpp.

◆ SetIntegrationMode()

void mfem::NonlinearFormIntegrator::SetIntegrationMode ( Mode  m)
inline

Definition at line 59 of file nonlininteg.hpp.

◆ SetIntegrationRule()

void mfem::NonlinearFormIntegrator::SetIntegrationRule ( const IntegrationRule ir)
inline

Prescribe a fixed IntegrationRule to use.

Definition at line 68 of file nonlininteg.hpp.

◆ SetIntRule()

virtual void mfem::NonlinearFormIntegrator::SetIntRule ( const IntegrationRule ir)
inlinevirtual

Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == NULL).

Reimplemented in mfem::SumIntegrator, mfem::InverseIntegrator, mfem::LumpedIntegrator, and mfem::TransposeIntegrator.

Definition at line 57 of file nonlininteg.hpp.

◆ SetNURBSPatchIntRule()

void mfem::NonlinearFormIntegrator::SetNURBSPatchIntRule ( NURBSMeshRules pr)
inline

For patchwise integration, SetNURBSPatchIntRule must be called.

Definition at line 62 of file nonlininteg.hpp.

◆ SetPAMemoryType()

void mfem::NonlinearFormIntegrator::SetPAMemoryType ( MemoryType  mt)
inline

Set the memory type used for GeometricFactors and other large allocations in PA extensions.

Definition at line 72 of file nonlininteg.hpp.

◆ SupportsCeed()

virtual bool mfem::NonlinearFormIntegrator::SupportsCeed ( ) const
inlinevirtual

Indicates whether this integrator can use a Ceed backend.

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

Definition at line 154 of file nonlininteg.hpp.

Member Data Documentation

◆ ceedOp

ceed::Operator* mfem::NonlinearFormIntegrator::ceedOp
protected

Definition at line 47 of file nonlininteg.hpp.

◆ integrationMode

Mode mfem::NonlinearFormIntegrator::integrationMode = Mode::ELEMENTWISE
protected

Definition at line 41 of file nonlininteg.hpp.

◆ IntRule

const IntegrationRule* mfem::NonlinearFormIntegrator::IntRule
protected

Definition at line 39 of file nonlininteg.hpp.

◆ pa_mt

MemoryType mfem::NonlinearFormIntegrator::pa_mt = MemoryType::DEFAULT
protected

Definition at line 49 of file nonlininteg.hpp.

◆ patchRules

NURBSMeshRules* mfem::NonlinearFormIntegrator::patchRules = nullptr
protected

Definition at line 44 of file nonlininteg.hpp.


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