MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::HyperbolicFormIntegrator Class Reference

Abstract hyperbolic form integrator, assembling (F(u, x), ∇v) and <F̂(u⁻,u⁺,x) n, [v]> terms for scalar finite elements. More...

#include <hyperbolic.hpp>

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

Public Member Functions

 HyperbolicFormIntegrator (const NumericalFlux &numFlux, const int IntOrderOffset=0, const real_t sign=1.)
 Construct a new Hyperbolic Form Integrator object.
 
void ResetMaxCharSpeed ()
 Reset the Max Char Speed 0.
 
real_t GetMaxCharSpeed ()
 
const FluxFunctionGetFluxFunction ()
 
void AssembleElementVector (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, Vector &elvect) override
 Implements (F(u), ∇v) with abstract F computed by FluxFunction::ComputeFlux()
 
void AssembleElementGrad (const FiniteElement &el, ElementTransformation &Tr, const Vector &elfun, DenseMatrix &grad) override
 Implements (J(u), ∇v) with abstract J computed by FluxFunction::ComputeFluxJacobian()
 
void AssembleFaceVector (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, Vector &elvect) override
 Implements <-F̂(u⁻,u⁺,x) n, [v]> with abstract F̂ computed by NumericalFlux::Eval() of the numerical flux object.
 
void AssembleFaceGrad (const FiniteElement &el1, const FiniteElement &el2, FaceElementTransformations &Tr, const Vector &elfun, DenseMatrix &elmat) override
 Implements <-Ĵ(u⁻,u⁺,x) n, [v]> with abstract Ĵ computed by NumericalFlux::Grad() of the numerical flux object.
 
- 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 AssemblePA (const FiniteElementSpace &fes)
 Method defining partial assembly.
 
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.
 
virtual real_t GetLocalStateEnergyPA (const Vector &x) const
 Compute the local (to the MPI rank) energy with partial assembly.
 
virtual void AddMultPA (const Vector &x, Vector &y) const
 Method for partially assembled action.
 
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.
 
virtual void AssembleMF (const FiniteElementSpace &fes)
 Method defining fully unassembled operator.
 
virtual void AddMultMF (const Vector &x, Vector &y) const
 
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.
 

Public Attributes

const int num_equations
 

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

Abstract hyperbolic form integrator, assembling (F(u, x), ∇v) and <F̂(u⁻,u⁺,x) n, [v]> terms for scalar finite elements.

This form integrator is coupled with a NumericalFlux that implements the numerical flux F̂ at the faces. The flux F is obtained from the FluxFunction assigned to the aforementioned NumericalFlux.

Definition at line 306 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ HyperbolicFormIntegrator()

mfem::HyperbolicFormIntegrator::HyperbolicFormIntegrator ( const NumericalFlux & numFlux,
const int IntOrderOffset = 0,
const real_t sign = 1. )

Construct a new Hyperbolic Form Integrator object.

Parameters
[in]numFluxnumerical flux
[in]IntOrderOffsetintegration order offset
[in]signsign of the convection term

Definition at line 21 of file hyperbolic.cpp.

Member Function Documentation

◆ AssembleElementGrad()

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

Implements (J(u), ∇v) with abstract J computed by FluxFunction::ComputeFluxJacobian()

Parameters
[in]ellocal finite element
[in]Trelement transformation
[in]elfunlocal coefficient of basis
[out]gradevaluated Jacobian

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 107 of file hyperbolic.cpp.

◆ AssembleElementVector()

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

Implements (F(u), ∇v) with abstract F computed by FluxFunction::ComputeFlux()

Parameters
[in]ellocal finite element
[in]Trelement transformation
[in]elfunlocal coefficient of basis
[out]elvectevaluated dual vector

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 44 of file hyperbolic.cpp.

◆ AssembleFaceGrad()

void mfem::HyperbolicFormIntegrator::AssembleFaceGrad ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Tr,
const Vector & elfun,
DenseMatrix & elmat )
overridevirtual

Implements <-Ĵ(u⁻,u⁺,x) n, [v]> with abstract Ĵ computed by NumericalFlux::Grad() of the numerical flux object.

Parameters
[in]el1finite element of the first element
[in]el2finite element of the second element
[in]Trface element transformations
[in]elfunlocal coefficient of basis from both elements
[out]elmatevaluated Jacobian matrix <-Ĵ(u⁻,u⁺,x) n, [v]>

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 264 of file hyperbolic.cpp.

◆ AssembleFaceVector()

void mfem::HyperbolicFormIntegrator::AssembleFaceVector ( const FiniteElement & el1,
const FiniteElement & el2,
FaceElementTransformations & Tr,
const Vector & elfun,
Vector & elvect )
overridevirtual

Implements <-F̂(u⁻,u⁺,x) n, [v]> with abstract F̂ computed by NumericalFlux::Eval() of the numerical flux object.

Parameters
[in]el1finite element of the first element
[in]el2finite element of the second element
[in]Trface element transformations
[in]elfunlocal coefficient of basis from both elements
[out]elvectevaluated dual vector <-F̂(u⁻,u⁺,x) n, [v]>

Reimplemented from mfem::NonlinearFormIntegrator.

Definition at line 177 of file hyperbolic.cpp.

◆ GetFluxFunction()

const FluxFunction & mfem::HyperbolicFormIntegrator::GetFluxFunction ( )
inline

Definition at line 360 of file hyperbolic.hpp.

◆ GetMaxCharSpeed()

real_t mfem::HyperbolicFormIntegrator::GetMaxCharSpeed ( )
inline

Definition at line 355 of file hyperbolic.hpp.

◆ ResetMaxCharSpeed()

void mfem::HyperbolicFormIntegrator::ResetMaxCharSpeed ( )
inline

Reset the Max Char Speed 0.

Definition at line 350 of file hyperbolic.hpp.

Member Data Documentation

◆ num_equations

const int mfem::HyperbolicFormIntegrator::num_equations

Definition at line 333 of file hyperbolic.hpp.


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