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

Component-wise upwinded flux. More...

#include <hyperbolic.hpp>

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

Public Member Functions

 ComponentwiseUpwindFlux (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

Component-wise upwinded flux.

Upwinded flux for scalar equations, a special case of Godunov or Engquist-Osher flux, is defined as follows: F̂ n = F(u⁺)n for dF(u)/du < 0 on [u⁻,u⁺] F̂ n = F(u⁻)n for dF(u)/du > 0 on [u⁻,u⁺]

Note
This construction assumes monotonous F(u,x) in u
Systems of equations are treated component-wise

Definition at line 539 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ ComponentwiseUpwindFlux()

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

Constructor for a flux function.

Parameters
fluxFunctionflux function F(u,x)

Definition at line 598 of file hyperbolic.cpp.

Member Function Documentation

◆ Average()

real_t mfem::ComponentwiseUpwindFlux::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()
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]fluxF̂ n = min(F(u⁻)n, F̄(u⁺,x)n) for u⁻ ≤ u⁺ or F̂ n = max(F(u⁻)n, F̄(u⁺,x)n) for u⁻ > u⁺
Returns
max(|dF(u⁺,x)/du⁺⋅n|, |dF(u⁻,x)/du⁻⋅n|)

Reimplemented from mfem::NumericalFlux.

Definition at line 670 of file hyperbolic.cpp.

◆ AverageGrad()

void mfem::ComponentwiseUpwindFlux::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()
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⁺) - F̄(u⁻,u⁺))n / (u⁺ - u⁻) when negative J(u⁻,x) n otherwise side = 2: min((F(u⁺) - F̄(u⁻,u⁺))n / (u⁺ - u⁻), 0)

Reimplemented from mfem::NumericalFlux.

Definition at line 697 of file hyperbolic.cpp.

◆ Eval()

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

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

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 = min(F(u⁻,x)n, F(u⁺,x)n) for u⁻ ≤ u⁺ or F̂ n = max(F(u⁻,x)n, F(u⁺,x)n) for u⁻ > u⁺
Returns
max(|dF(u⁺,x)/du⁺⋅n|, |dF(u⁻,x)/du⁻⋅n|)

Implements mfem::NumericalFlux.

Definition at line 610 of file hyperbolic.cpp.

◆ Grad()

void mfem::ComponentwiseUpwindFlux::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: max(J(u⁻,x)n, 0) side = 2: min(J(u⁺,x)n, 0)

Reimplemented from mfem::NumericalFlux.

Definition at line 635 of file hyperbolic.cpp.

Member Data Documentation

◆ fluxN1

Vector mfem::ComponentwiseUpwindFlux::fluxN1
mutableprotected

Definition at line 634 of file hyperbolic.hpp.

◆ fluxN2

Vector mfem::ComponentwiseUpwindFlux::fluxN2
protected

Definition at line 634 of file hyperbolic.hpp.

◆ JDotN

DenseMatrix mfem::ComponentwiseUpwindFlux::JDotN
mutableprotected

Definition at line 635 of file hyperbolic.hpp.


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