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

Rusanov flux, also known as local Lax-Friedrichs, F̂ n = ½(F(u⁺,x)n + F(u⁻,x)n) - ½λ(u⁺ - u⁻) where λ is the maximum characteristic speed. More...

#include <hyperbolic.hpp>

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

Public Member Functions

 RusanovFlux (const FluxFunction &fluxFunction)
 Constructor for a flux function.
 
real_t Eval (const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
 Normal numerical flux F̂(u⁻,u⁺,x) n.
 
void Grad (int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
 Jacobian of normal numerical flux F̂(u⁻,u⁺,x) n.
 
real_t Average (const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
 Average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,u,x) n.
 
void AverageGrad (int side, const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, DenseMatrix &grad) const override
 Jacobian of average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,u,x) n.
 
- Public Member Functions inherited from mfem::NumericalFlux
 NumericalFlux (const FluxFunction &fluxFunction)
 Constructor for a flux function.
 
virtual ~NumericalFlux ()=default
 
const FluxFunctionGetFluxFunction () const
 Get flux function F.
 

Protected Attributes

Vector fluxN1
 
Vector fluxN2
 
DenseMatrix JDotN
 
- Protected Attributes inherited from mfem::NumericalFlux
const FluxFunctionfluxFunction
 

Detailed Description

Rusanov flux, also known as local Lax-Friedrichs, F̂ n = ½(F(u⁺,x)n + F(u⁻,x)n) - ½λ(u⁺ - u⁻) where λ is the maximum characteristic speed.

Note
The implementation assumes monotonous |dF(u,x)/du⋅n| in u, so the maximum characteristic speed λ for any interval [u⁻, u⁺] is given by max(|dF(u⁺,x)/du⁺⋅n|, |dF(u⁻,x)/du⁻⋅n|).

Definition at line 428 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ RusanovFlux()

mfem::RusanovFlux::RusanovFlux ( const FluxFunction & fluxFunction)

Constructor for a flux function.

Parameters
fluxFunctionflux function F(u,x)

Definition at line 441 of file hyperbolic.cpp.

Member Function Documentation

◆ Average()

real_t mfem::RusanovFlux::Average ( const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
Vector & flux ) const
overridevirtual

Average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,u,x) n.

Note
The average normal flux F̄ n is required to be implemented in FluxFunction::ComputeAvgFluxDotN()
Systems of equations are treated component-wise
Parameters
[in]state1state value (u⁻) of the beginning of the interval (num_equations)
[in]state2state value (u⁺) of the end of the interval (num_equations)
[in]nornormal vector (not a unit vector) (dim)
[in]Trface element transformation
[out]flux½(F̄(u⁻,u⁺,x)n + F(u⁻,x)n) - ¼λ(u⁺ - u⁻)
Returns
max(|dF(u⁺,x)/du⁺⋅n|, |dF(u⁻,x)/du⁻⋅n|)

Reimplemented from mfem::NumericalFlux.

Definition at line 506 of file hyperbolic.cpp.

◆ AverageGrad()

void mfem::RusanovFlux::AverageGrad ( int side,
const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
DenseMatrix & grad ) const
overridevirtual

Jacobian of average normal numerical flux over the interval [u⁻, u⁺] in the second argument of the flux F̂(u⁻,u,x) n.

Note
The average normal flux F̄ n is required to be implemented in FluxFunction::ComputeAvgFluxDotN() and the Jacobian of flux J n in FluxFunction::ComputeFluxJacobianDotN()
Only the diagonal terms of the J n are considered, i.e., systems are treated as a set of independent equations
Parameters
[in]sidegradient w.r.t the first (u⁻) or second argument (u⁺)
[in]state1state value (u⁻) of the beginning of the interval (num_equations)
[in]state2state value (u⁺) of the end of the interval (num_equations)
[in]nornormal vector (not a unit vector) (dim)
[in]Trface element transformation
[out]gradJacobian of F̄(u⁻,u⁺,x) n side = 1: ½(F̄(u⁻,u⁺,x)n - F(u⁻,x)n) / (u⁺ - u⁻) - ½J(u⁻,x)n + ¼λ side = 2: ½(F(u⁺,x)n - F̄(u⁻,u⁺,x)n) / (u⁺ - u⁻) - ¼λ

Reimplemented from mfem::NumericalFlux.

Definition at line 527 of file hyperbolic.cpp.

◆ Eval()

real_t mfem::RusanovFlux::Eval ( const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
Vector & flux ) const
overridevirtual

Normal numerical flux F̂(u⁻,u⁺,x) n.

Note
Systems of equations are treated component-wise
Parameters
[in]state1state value (u⁻) at a point from the first element (num_equations)
[in]state2state value (u⁺) at a point from the second element (num_equations)
[in]nornormal vector (not a unit vector) (dim)
[in]Trface element transformation
[out]fluxF̂ n = ½(F(u⁺,x)n + F(u⁻,x)n) - ½λ(u⁺ - u⁻)
Returns
max(|dF(u⁺,x)/du⁺⋅n|, |dF(u⁻,x)/du⁻⋅n|)

Implements mfem::NumericalFlux.

Definition at line 450 of file hyperbolic.cpp.

◆ Grad()

void mfem::RusanovFlux::Grad ( int side,
const Vector & state1,
const Vector & state2,
const Vector & nor,
FaceElementTransformations & Tr,
DenseMatrix & grad ) const
overridevirtual

Jacobian of normal numerical flux F̂(u⁻,u⁺,x) n.

Note
The Jacobian of flux J n is required to be implemented in FluxFunction::ComputeFluxJacobianDotN()
Parameters
[in]sidegradient w.r.t the first (u⁻) or second argument (u⁺)
[in]state1state value (u⁻) of the beginning of the interval (num_equations)
[in]state2state value (u⁺) of the end of the interval (num_equations)
[in]nornormal vector (not a unit vector) (dim)
[in]Trface element transformation
[out]gradJacobian of F(u⁻,u⁺,x) n side = 1: ½J(u⁻,x)n + ½λ side = 2: ½J(u⁺,x)n - ½λ

Reimplemented from mfem::NumericalFlux.

Definition at line 470 of file hyperbolic.cpp.

Member Data Documentation

◆ fluxN1

Vector mfem::RusanovFlux::fluxN1
mutableprotected

Definition at line 524 of file hyperbolic.hpp.

◆ fluxN2

Vector mfem::RusanovFlux::fluxN2
protected

Definition at line 524 of file hyperbolic.hpp.

◆ JDotN

DenseMatrix mfem::RusanovFlux::JDotN
mutableprotected

Definition at line 525 of file hyperbolic.hpp.


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