MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
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
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, Vector &elvect) override
 Perform the local action of the NonlinearFormIntegrator.
 
void AssembleElementGrad (const FiniteElement &el, ElementTransformation &trans, const Vector &elfun, DenseMatrix &elmat) override
 Assemble the local gradient matrix.
 
void AssemblePA (const FiniteElementSpace &fes) override
 Method defining partial assembly.
 
void AssembleMF (const FiniteElementSpace &fes) override
 Method defining fully unassembled operator.
 
void AddMultPA (const Vector &x, Vector &y) const override
 Method for partially assembled action.
 
void AddMultMF (const Vector &x, Vector &y) const override
 
virtual void AssemblePA (const FiniteElementSpace &trial_fes, const FiniteElementSpace &test_fes)
 
- Public Member Functions inherited from mfem::NonlinearFormIntegrator
void SetIntegrationMode (Mode m)
 
bool Patchwise () const
 
void SetPAMemoryType (MemoryType mt)
 
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.
 
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.
 
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.
 
virtual bool SupportsCeed () const
 Indicates whether this integrator can use a Ceed backend.
 
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 &fe, const ElementTransformation &T)
 

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::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)
 

Additional Inherited Members

- Public Types inherited from mfem::NonlinearFormIntegrator
enum  Mode { ELEMENTWISE = 0 , PATCHWISE = 1 , PATCHWISE_REDUCED = 2 }
 
- 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

Definition at line 380 of file nonlininteg.hpp.

Constructor & Destructor Documentation

◆ VectorConvectionNLFIntegrator() [1/2]

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

Definition at line 393 of file nonlininteg.hpp.

◆ VectorConvectionNLFIntegrator() [2/2]

mfem::VectorConvectionNLFIntegrator::VectorConvectionNLFIntegrator ( )
default

Member Function Documentation

◆ AddMultMF()

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

Definition at line 44 of file nonlininteg_vecconvection_mf.cpp.

◆ AddMultPA()

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

Definition at line 807 of file nonlininteg_vecconvection_pa.cpp.

◆ AssembleElementGrad()

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

Assemble the local gradient matrix.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 781 of file nonlininteg.cpp.

◆ AssembleElementVector()

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

Perform the local action of the NonlinearFormIntegrator.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 744 of file nonlininteg.cpp.

◆ AssembleMF()

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

Method defining fully unassembled operator.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 18 of file nonlininteg_vecconvection_mf.cpp.

◆ AssemblePA() [1/2]

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

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 19 of file nonlininteg_vecconvection_pa.cpp.

◆ AssemblePA() [2/2]

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

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.

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 99 of file nonlininteg.cpp.

◆ GetDefaultIntegrationRule()

const IntegrationRule * mfem::VectorConvectionNLFIntegrator::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 422 of file nonlininteg.hpp.

◆ GetRule()

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

Definition at line 737 of file nonlininteg.cpp.


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