MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::FluxFunction Class Referenceabstract

Abstract class for hyperbolic flux for a system of hyperbolic conservation laws. More...

#include <hyperbolic.hpp>

Inheritance diagram for mfem::FluxFunction:
[legend]

Public Member Functions

 FluxFunction (const int num_equations, const int dim)
 
virtual ~FluxFunction ()
 
virtual real_t ComputeFlux (const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const =0
 Compute flux F(u, x). Must be implemented in a derived class.
 
virtual real_t ComputeFluxDotN (const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const
 Compute normal flux F(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a full dense matrix for flux.
 
virtual real_t ComputeAvgFlux (const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux_) const
 Compute average flux over the given interval of states. Optionally overloaded in a derived class.
 
virtual real_t ComputeAvgFluxDotN (const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const
 Compute average normal flux over the given interval of states. Optionally overloaded in a derived class.
 
virtual void ComputeFluxJacobian (const Vector &state, ElementTransformation &Tr, DenseTensor &J_) const
 Compute flux Jacobian J(u, x). Optionally overloaded in a derived class when Jacobian is necessary (e.g. Newton iteration, flux limiter)
 
virtual void ComputeFluxJacobianDotN (const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const
 Compute normal flux Jacobian J(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a full dense tensor for Jacobian.
 

Public Attributes

const int num_equations
 
const int dim
 

Detailed Description

Abstract class for hyperbolic flux for a system of hyperbolic conservation laws.

Definition at line 57 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ FluxFunction()

mfem::FluxFunction::FluxFunction ( const int num_equations,
const int dim )
inline

Definition at line 63 of file hyperbolic.hpp.

◆ ~FluxFunction()

virtual mfem::FluxFunction::~FluxFunction ( )
inlinevirtual

Definition at line 66 of file hyperbolic.hpp.

Member Function Documentation

◆ ComputeAvgFlux()

virtual real_t mfem::FluxFunction::ComputeAvgFlux ( const Vector & state1,
const Vector & state2,
ElementTransformation & Tr,
DenseMatrix & flux_ ) const
inlinevirtual

Compute average flux over the given interval of states. Optionally overloaded in a derived class.

The average flux is defined as F̄(u1,u2) = ∫ F(u) du / (u2 - u1) for u ∈ [u1,u2], where u1 is the first state (state1) and the u2 the second state (state2), while F(u) is the flux as defined in ComputeFlux().

Used in the default implementation of ComputeAvgFluxDotN().

Parameters
[in]state1state of the beginning of the interval (num_equations)
[in]state2state of the end of the interval (num_equations)
[in]Trelement transformation
[out]flux_average flux from the given element at the current integration point (num_equations, dim)
Returns
real_t maximum characteristic speed |dF(u,x)/du| over the interval [u1,u2]

Reimplemented in mfem::AdvectionFlux, and mfem::BurgersFlux.

Definition at line 118 of file hyperbolic.hpp.

◆ ComputeAvgFluxDotN()

real_t mfem::FluxFunction::ComputeAvgFluxDotN ( const Vector & state1,
const Vector & state2,
const Vector & normal,
FaceElementTransformations & Tr,
Vector & fluxDotN ) const
virtual

Compute average normal flux over the given interval of states. Optionally overloaded in a derived class.

The average normal flux is defined as F̄(u1,u2)n = ∫ F(u)n du / (u2 - u1) for u ∈ [u1,u2], where u1 is the first state (state1) and the u2 the second state (state2), while n is the normal and F(u) is the flux as defined in ComputeFlux().

Used in NumericalFlux::Average() and NumericalFlux::AverageGrad() for evaluation of the average normal flux on a face.

Parameters
[in]state1state of the beginning of the interval (num_equations)
[in]state2state of the end of the interval (num_equations)
[in]normalnormal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]fluxDotNaverage normal flux from the given element at the current integration point (num_equations)
Returns
real_t maximum (normal) characteristic speed |dF(u,x)/du⋅n| over the interval [u1,u2]

Reimplemented in mfem::AdvectionFlux, and mfem::BurgersFlux.

Definition at line 408 of file hyperbolic.cpp.

◆ ComputeFlux()

virtual real_t mfem::FluxFunction::ComputeFlux ( const Vector & state,
ElementTransformation & Tr,
DenseMatrix & flux ) const
pure virtual

Compute flux F(u, x). Must be implemented in a derived class.

Used in HyperbolicFormIntegrator::AssembleElementVector() for evaluation of (F(u), ∇v) and in the default implementation of ComputeFluxDotN() for evaluation of F(u)⋅n.

Parameters
[in]statestate at the current integration point (num_equations)
[in]Trelement transformation
[out]fluxflux from the given element at the current integration point (num_equations, dim)
Returns
real_t maximum characteristic speed |dF(u,x)/du|
Note
One can put assertion in here to detect non-physical solution

Implemented in mfem::AdvectionFlux, mfem::BurgersFlux, mfem::EulerFlux, and mfem::ShallowWaterFlux.

◆ ComputeFluxDotN()

real_t mfem::FluxFunction::ComputeFluxDotN ( const Vector & state,
const Vector & normal,
FaceElementTransformations & Tr,
Vector & fluxDotN ) const
virtual

Compute normal flux F(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a full dense matrix for flux.

Used in NumericalFlux for evaluation of the normal flux on a face.

Parameters
[in]statestate at the current integration point (num_equations)
[in]normalnormal vector, see mfem::CalcOrtho() (dim)
[in]Trface transformation
[out]fluxDotNnormal flux from the given element at the current integration point (num_equations)
Returns
real_t maximum (normal) characteristic speed |dF(u,x)/du⋅n|

Reimplemented in mfem::AdvectionFlux, mfem::BurgersFlux, mfem::EulerFlux, and mfem::ShallowWaterFlux.

Definition at line 393 of file hyperbolic.cpp.

◆ ComputeFluxJacobian()

virtual void mfem::FluxFunction::ComputeFluxJacobian ( const Vector & state,
ElementTransformation & Tr,
DenseTensor & J_ ) const
inlinevirtual

Compute flux Jacobian J(u, x). Optionally overloaded in a derived class when Jacobian is necessary (e.g. Newton iteration, flux limiter)

Used in HyperbolicFormIntegrator::AssembleElementGrad() for evaluation of Jacobian of the flux in an element and in the default implementation of ComputeFluxJacobianDotN().

Parameters
[in]statestate at the current integration point (num_equations)
[in]Trelement transformation
[out]J_flux Jacobian, \( J(i,j,d) = dF_{id} / du_j \)

Reimplemented in mfem::AdvectionFlux, and mfem::BurgersFlux.

Definition at line 159 of file hyperbolic.hpp.

◆ ComputeFluxJacobianDotN()

void mfem::FluxFunction::ComputeFluxJacobianDotN ( const Vector & state,
const Vector & normal,
ElementTransformation & Tr,
DenseMatrix & JDotN ) const
virtual

Compute normal flux Jacobian J(u, x)⋅n. Optionally overloaded in a derived class to avoid creating a full dense tensor for Jacobian.

Used in NumericalFlux for evaluation of Jacobian of the normal flux on a face.

Parameters
[in]statestate at the current integration point (num_equations)
[in]normalnormal vector, see mfem::CalcOrtho() (dim)
[in]Trelement transformation
[out]JDotNnormal flux Jacobian, \( JDotN(i,j) = d(F_{id} n_d) / du_j \)

Reimplemented in mfem::AdvectionFlux, and mfem::BurgersFlux.

Definition at line 423 of file hyperbolic.cpp.

Member Data Documentation

◆ dim

const int mfem::FluxFunction::dim

Definition at line 61 of file hyperbolic.hpp.

◆ num_equations

const int mfem::FluxFunction::num_equations

Definition at line 60 of file hyperbolic.hpp.


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