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

Abstract class for numerical flux for a system of hyperbolic conservation laws on a face with states, fluxes and characteristic speed. More...

#include <hyperbolic.hpp>

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

Public Member Functions

 NumericalFlux (const FluxFunction &fluxFunction)
 Constructor for a flux function.
 
virtual real_t Eval (const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const =0
 Evaluates normal numerical flux for the given states and normal. Must be implemented in a derived class.
 
virtual void Grad (int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const
 Evaluates Jacobian of the normal numerical flux for the given states and normal. Optionally overloaded in a derived class.
 
virtual real_t Average (const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const
 Evaluates average normal numerical flux over the interval between the given end states in the second argument and for the given normal. Optionally overloaded in a derived class.
 
virtual void AverageGrad (int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const
 Evaluates Jacobian of the average normal numerical flux over the interval between the given end states in the second argument and for the given normal. Optionally overloaded in a derived class.
 
virtual ~NumericalFlux ()=default
 
const FluxFunctionGetFluxFunction () const
 Get flux function F.
 

Protected Attributes

const FluxFunctionfluxFunction
 

Detailed Description

Abstract class for numerical flux for a system of hyperbolic conservation laws on a face with states, fluxes and characteristic speed.

Definition at line 193 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ NumericalFlux()

mfem::NumericalFlux::NumericalFlux ( const FluxFunction & fluxFunction)
inline

Constructor for a flux function.

Parameters
fluxFunctionflux function F(u,x)

Definition at line 200 of file hyperbolic.hpp.

◆ ~NumericalFlux()

virtual mfem::NumericalFlux::~NumericalFlux ( )
virtualdefault

Member Function Documentation

◆ Average()

virtual real_t mfem::NumericalFlux::Average ( const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
Vector & flux ) const
inlinevirtual

Evaluates average normal numerical flux over the interval between the given end states in the second argument and for the given normal. Optionally overloaded in a derived class.

Presently, not used. Reserved for future use.

Parameters
[in]state1state value of the beginning of the interval (num_equations)
[in]state2state value of the end of the interval (num_equations)
[in]norscaled normal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]fluxnumerical flux (num_equations)
Returns
real_t maximum characteristic speed |dF(u,x)/du⋅n|

Reimplemented in mfem::ComponentwiseUpwindFlux, and mfem::RusanovFlux.

Definition at line 258 of file hyperbolic.hpp.

◆ AverageGrad()

virtual void mfem::NumericalFlux::AverageGrad ( int side,
const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
DenseMatrix & grad ) const
inlinevirtual

Evaluates Jacobian of the average normal numerical flux over the interval between the given end states in the second argument and for the given normal. Optionally overloaded in a derived class.

Presently, not used. Reserved for future use.

Parameters
[in]sideindicates gradient w.r.t. the first (side = 1) or second (side = 2) state
[in]state1state value of the beginning of the interval (num_equations)
[in]state2state value of the end of the interval (num_equations)
[in]norscaled normal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]gradJacobian of the average normal numerical flux (num_equations, dim)

Reimplemented in mfem::ComponentwiseUpwindFlux, and mfem::RusanovFlux.

Definition at line 280 of file hyperbolic.hpp.

◆ Eval()

virtual real_t mfem::NumericalFlux::Eval ( const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
Vector & flux ) const
pure virtual

Evaluates normal numerical flux for the given states and normal. Must be implemented in a derived class.

Used in HyperbolicFormIntegrator::AssembleFaceVector() for evaluation of <F̂(u⁻,u⁺,x) n, [v]> term at the face.

Parameters
[in]state1state value at a point from the first element (num_equations)
[in]state2state value at a point from the second element (num_equations)
[in]norscaled normal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]fluxnumerical flux (num_equations)
Returns
real_t maximum characteristic speed |dF(u,x)/du⋅n|

Implemented in mfem::ComponentwiseUpwindFlux, and mfem::RusanovFlux.

◆ GetFluxFunction()

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

Get flux function F.

Returns
constant reference to the flux function.

Definition at line 289 of file hyperbolic.hpp.

◆ Grad()

virtual void mfem::NumericalFlux::Grad ( int side,
const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
DenseMatrix & grad ) const
inlinevirtual

Evaluates Jacobian of the normal numerical flux for the given states and normal. Optionally overloaded in a derived class.

Used in HyperbolicFormIntegrator::AssembleFaceGrad() for Jacobian of the term <F̂(u⁻,u⁺,x) n, [v]> at the face.

Parameters
[in]sideindicates gradient w.r.t. the first (side = 1) or second (side = 2) state
[in]state1state value of the beginning of the interval (num_equations)
[in]state2state value of the end of the interval (num_equations)
[in]norscaled normal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]gradJacobian of normal numerical flux (num_equations, dim)

Reimplemented in mfem::ComponentwiseUpwindFlux, and mfem::RusanovFlux.

Definition at line 238 of file hyperbolic.hpp.

Member Data Documentation

◆ fluxFunction

const FluxFunction& mfem::NumericalFlux::fluxFunction
protected

Definition at line 292 of file hyperbolic.hpp.


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