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

#include <nonlininteg.hpp>

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

Public Member Functions

 VectorConvectionNLFIntegrator (Coefficient &q)
 
 VectorConvectionNLFIntegrator ()=default
 
virtual void AssembleElementVector (const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect)
 Perform the local action of the NonlinearFormIntegrator. More...
 
virtual void AssembleElementGrad (const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat)
 Assemble the local gradient matrix. More...
 
virtual void AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly. More...
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator. More...
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action. More...
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
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 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 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 &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 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 ()
 

Static Public Member Functions

static const IntegrationRuleGetRule (const FiniteElement &fe, ElementTransformation &T)
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0, PATCHWISE = 1, PATCHWISE_REDUCED = 2 }
 
- Protected Member Functions inherited from mfem::NonlinearFormIntegrator
 NonlinearFormIntegrator (const IntegrationRule *ir=NULL)
 
- Protected Attributes inherited from mfem::NonlinearFormIntegrator
const IntegrationRuleIntRule
 
Mode integrationMode = Mode::ELEMENTWISE
 
NURBSMeshRulespatchRules = nullptr
 
ceed::OperatorceedOp
 
MemoryType pa_mt = MemoryType::DEFAULT
 

Detailed Description

Definition at line 390 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ VectorConvectionNLFIntegrator() [1/2]

mfem::VectorConvectionNLFIntegrator::VectorConvectionNLFIntegrator ( Coefficient q)
inline

Definition at line 403 of file nonlininteg.hpp.

◆ VectorConvectionNLFIntegrator() [2/2]

mfem::VectorConvectionNLFIntegrator::VectorConvectionNLFIntegrator ( )
default

Member Function Documentation

◆ AddMultMF()

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

◆ AddMultPA()

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

◆ AssembleElementGrad()

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

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Reimplemented in mfem::SkewSymmetricVectorConvectionNLFIntegrator, and mfem::ConvectiveVectorConvectionNLFIntegrator.

Definition at line 779 of file nonlininteg.cpp.

◆ AssembleElementVector()

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

Perform the local action of the NonlinearFormIntegrator.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 742 of file nonlininteg.cpp.

◆ AssembleMF()

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

Method defining fully unassembled operator.

Reimplemented from mfem::NonlinearFormIntegrator.

◆ AssemblePA() [1/3]

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

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/3]

virtual void mfem::VectorConvectionNLFIntegrator::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().

Reimplemented from mfem::NonlinearFormIntegrator.

◆ GetRule()

const IntegrationRule & mfem::VectorConvectionNLFIntegrator::GetRule ( const FiniteElement fe,
ElementTransformation T 
)
static

Definition at line 735 of file nonlininteg.cpp.


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