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

Advection flux. More...

#include <hyperbolic.hpp>

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

Public Member Functions

 AdvectionFlux (VectorCoefficient &b)
 Construct AdvectionFlux FluxFunction with given velocity.
 
real_t ComputeFlux (const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
 Compute F(u)
 
real_t ComputeFluxDotN (const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
 Compute F(u) n.
 
real_t ComputeAvgFlux (const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux) const override
 Compute average flux F̄(u)
 
real_t ComputeAvgFluxDotN (const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
 Compute average flux F̄(u) n.
 
void ComputeFluxJacobian (const Vector &state, ElementTransformation &Tr, DenseTensor &J) const override
 Compute J(u)
 
void ComputeFluxJacobianDotN (const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const override
 Compute J(u) n.
 
- Public Member Functions inherited from mfem::FluxFunction
 FluxFunction (const int num_equations, const int dim)
 
virtual ~FluxFunction ()
 

Additional Inherited Members

- Public Attributes inherited from mfem::FluxFunction
const int num_equations
 
const int dim
 

Detailed Description

Advection flux.

Definition at line 640 of file hyperbolic.hpp.

Constructor & Destructor Documentation

◆ AdvectionFlux()

mfem::AdvectionFlux::AdvectionFlux ( VectorCoefficient & b)
inline

Construct AdvectionFlux FluxFunction with given velocity.

Parameters
bvelocity coefficient, possibly depends on space

Definition at line 655 of file hyperbolic.hpp.

Member Function Documentation

◆ ComputeAvgFlux()

real_t mfem::AdvectionFlux::ComputeAvgFlux ( const Vector & state1,
const Vector & state2,
ElementTransformation & Tr,
DenseMatrix & flux ) const
overridevirtual

Compute average flux F̄(u)

Parameters
state1state value (u⁻) of the beginning of the interval
state2state value (u⁺) of the end of the interval
Trcurrent element transformation with the integration point
fluxF̄(u) = (u⁻+u⁺)/2*bᵀ
Returns
real_t maximum characteristic speed, |b|

Reimplemented from mfem::FluxFunction.

Definition at line 802 of file hyperbolic.cpp.

◆ ComputeAvgFluxDotN()

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

Compute average flux F̄(u) n.

Parameters
state1state value (u⁻) of the beginning of the interval
state2state value (u⁺) of the end of the interval
normalnormal vector, usually not a unit vector
Trcurrent element transformation with the integration point
fluxDotNF̄(u) n = (u⁻+u⁺)/2*(bᵀn)
Returns
real_t maximum characteristic speed, |b|

Reimplemented from mfem::FluxFunction.

Definition at line 816 of file hyperbolic.cpp.

◆ ComputeFlux()

real_t mfem::AdvectionFlux::ComputeFlux ( const Vector & state,
ElementTransformation & Tr,
DenseMatrix & flux ) const
overridevirtual

Compute F(u)

Parameters
statestate (u) at current integration point
Trcurrent element transformation with the integration point
fluxF(u) = ubᵀ
Returns
real_t maximum characteristic speed, |b|

Implements mfem::FluxFunction.

Definition at line 777 of file hyperbolic.cpp.

◆ ComputeFluxDotN()

real_t mfem::AdvectionFlux::ComputeFluxDotN ( const Vector & state,
const Vector & normal,
FaceElementTransformations & Tr,
Vector & fluxDotN ) const
overridevirtual

Compute F(u) n.

Parameters
statestate (u) at current integration point
normalnormal vector, usually not a unit vector
Trcurrent element transformation with the integration point
fluxDotNF(u) n = u (bᵀn)
Returns
real_t maximum characteristic speed, |b|

Reimplemented from mfem::FluxFunction.

Definition at line 789 of file hyperbolic.cpp.

◆ ComputeFluxJacobian()

void mfem::AdvectionFlux::ComputeFluxJacobian ( const Vector & state,
ElementTransformation & Tr,
DenseTensor & J ) const
overridevirtual

Compute J(u)

Parameters
statestate (u) at current integration point
Trcurrent element transformation with the integration point
JJ(u) = diag(b)

Reimplemented from mfem::FluxFunction.

Definition at line 829 of file hyperbolic.cpp.

◆ ComputeFluxJacobianDotN()

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

Compute J(u) n.

Parameters
statestate (u) at current integration point
normalnormal vector, usually not a unit vector
Trcurrent element transformation with the integration point
JDotNJ(u) n = bᵀn

Reimplemented from mfem::FluxFunction.

Definition at line 844 of file hyperbolic.cpp.


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