13#ifndef MFEM_HYPERBOLIC
14#define MFEM_HYPERBOLIC
121 { MFEM_ABORT(
"Not Implemented."); }
162 { MFEM_ABORT(
"Not Implemented."); }
181#ifndef MFEM_THREAD_SAFE
241 { MFEM_ABORT(
"Not implemented."); }
261 { MFEM_ABORT(
"Not implemented."); }
283 { MFEM_ABORT(
"Not implemented."); }
313 const int IntOrderOffset;
315#ifndef MFEM_THREAD_SAFE
343 const int IntOrderOffset = 0,
352 max_char_speed = 0.0;
357 return max_char_speed;
452 Vector &flux)
const override;
494 Vector &flux)
const override;
523#ifndef MFEM_THREAD_SAFE
563 Vector &flux)
const override;
605 Vector &flux)
const override;
633#ifndef MFEM_THREAD_SAFE
644#ifndef MFEM_THREAD_SAFE
658#ifndef MFEM_THREAD_SAFE
685 Vector &fluxDotN)
const override;
711 Vector &fluxDotN)
const override;
773 Vector &fluxDotN)
const override;
803 Vector &fluxDotN)
const override;
869 Vector &fluxN)
const override;
876 const real_t specific_heat_ratio;
889 specific_heat_ratio(specific_heat_ratio) {}
913 Vector &fluxN)
const override;
AdvectionFlux(VectorCoefficient &b)
Construct AdvectionFlux FluxFunction with given velocity.
void ComputeFluxJacobian(const Vector &state, ElementTransformation &Tr, DenseTensor &J) const override
Compute J(u)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(u)
void ComputeFluxJacobianDotN(const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const override
Compute J(u) n.
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.
real_t ComputeAvgFluxDotN(const Vector &state1, const Vector &state2, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute average flux F̄(u) n.
BurgersFlux(const int dim)
Construct BurgersFlux FluxFunction with given spatial dimension.
void ComputeFluxJacobianDotN(const Vector &state, const Vector &normal, ElementTransformation &Tr, DenseMatrix &JDotN) const override
Compute J(u) n.
real_t ComputeAvgFlux(const Vector &state1, const Vector &state2, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute average flux F̄(u)
real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxDotN) const override
Compute F(u) n.
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(u)
void ComputeFluxJacobian(const Vector &state, ElementTransformation &Tr, DenseTensor &J) const override
Compute J(u)
Component-wise upwinded flux.
ComponentwiseUpwindFlux(const FluxFunction &fluxFunction)
Constructor for a flux function.
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 fl...
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⁻,...
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 Eval(const Vector &state1, const Vector &state2, const Vector &nor, FaceElementTransformations &Tr, Vector &flux) const override
Normal numerical flux F̂(u⁻,u⁺,x) n.
Data type dense matrix using column-major storage.
Rank 3 tensor (array of matrices)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(ρ, ρu, E)
EulerFlux(const int dim, const real_t specific_heat_ratio)
Construct a new EulerFlux FluxFunction with given spatial dimension and specific heat ratio.
real_t ComputeFluxDotN(const Vector &x, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxN) const override
Compute normal flux, F(ρ, ρu, E)n.
Abstract class for all finite elements.
Abstract class for hyperbolic flux for a system of hyperbolic conservation laws.
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 cla...
FluxFunction(const int num_equations, const int dim)
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 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...
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 dens...
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 ...
Abstract class for numerical flux for a system of hyperbolic conservation laws on a face with states,...
NumericalFlux(const FluxFunction &fluxFunction)
Constructor for a flux function.
virtual ~NumericalFlux()=default
const FluxFunction & fluxFunction
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 cla...
const FluxFunction & GetFluxFunction() const
Get flux function F.
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 ...
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 state...
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 overloade...
Rusanov flux, also known as local Lax-Friedrichs, F̂ n = ½(F(u⁺,x)n + F(u⁻,x)n) - ½λ(u⁺ - u⁻) where λ...
RusanovFlux(const FluxFunction &fluxFunction)
Constructor for a flux function.
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⁻,...
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.
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 fl...
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.
ShallowWaterFlux(const int dim, const real_t g=9.8)
Construct a new ShallowWaterFlux FluxFunction with given spatial dimension and gravity constant.
real_t ComputeFluxDotN(const Vector &state, const Vector &normal, FaceElementTransformations &Tr, Vector &fluxN) const override
Compute normal flux, F(h, hu)
real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const override
Compute F(h, hu)
Base class for vector Coefficients that optionally depend on time and space.
void SetSize(int s)
Resize the vector to size s.
MFEM_DEPRECATED typedef NumericalFlux RiemannSolver