[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real > Class Template Reference

#include <convective_numerical_flux.hpp>

Inheritance diagram for PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >:
Collaboration diagram for PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >:

Public Member Functions

 L2RoeRiemannSolverDissipation (std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input)
 Constructor.
 
void evaluate_entropy_fix (const std::array< real, 3 > &eig_L, const std::array< real, 3 > &eig_R, std::array< real, 3 > &eig_RoeAvg, const real vel2_ravg, const real sound_ravg) const
 
void evaluate_additional_modifications (const std::array< real, nstate > &soln_int, const std::array< real, nstate > &soln_ext, const std::array< real, 3 > &eig_L, const std::array< real, 3 > &eig_R, real &dV_normal, dealii::Tensor< 1, dim, real > &dV_tangent) const
 
- Public Member Functions inherited from PHiLiP::NumericalFlux::RoeBaseRiemannSolverDissipation< dim, nstate, real >
 RoeBaseRiemannSolverDissipation (std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input)
 Constructor.
 
std::array< real, nstate > evaluate_riemann_solver_dissipation (const std::array< real, nstate > &soln_int, const std::array< real, nstate > &soln_ext, const dealii::Tensor< 1, dim, real > &normal1) const
 
- Public Member Functions inherited from PHiLiP::NumericalFlux::RiemannSolverDissipation< dim, nstate, real >
virtual ~RiemannSolverDissipation ()=default
 < Base class destructor required for abstract classes.
 

Protected Member Functions

void evaluate_shock_indicator (const std::array< real, 3 > &eig_L, const std::array< real, 3 > &eig_R, int &ssw_LEFT, int &ssw_RIGHT) const
 

Additional Inherited Members

- Protected Attributes inherited from PHiLiP::NumericalFlux::RoeBaseRiemannSolverDissipation< dim, nstate, real >
const std::shared_ptr< Physics::Euler< dim, nstate, real > > euler_physics
 Numerical flux requires physics to evaluate convective eigenvalues.
 

Detailed Description

template<int dim, int nstate, typename real>
class PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >

L2Roe flux with entropy fix. Derived from RoeBase. — Reference: Osswald et al. (2016 L2Roe)

Definition at line 192 of file convective_numerical_flux.hpp.

Member Function Documentation

◆ evaluate_additional_modifications()

template<int dim, int nstate, typename real >
void PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >::evaluate_additional_modifications ( const std::array< real, nstate > &  soln_int,
const std::array< real, nstate > &  soln_ext,
const std::array< real, 3 > &  eig_L,
const std::array< real, 3 > &  eig_R,
real &  dV_normal,
dealii::Tensor< 1, dim, real > &  dV_tangent 
) const
virtual

Osswald's two modifications to Roe-Pike scheme –> L2Roe — Scale jump in (1) normal and (2) tangential velocities using a blending factor

Implements PHiLiP::NumericalFlux::RoeBaseRiemannSolverDissipation< dim, nstate, real >.

Definition at line 302 of file convective_numerical_flux.cpp.

◆ evaluate_entropy_fix()

template<int dim, int nstate, typename real >
void PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >::evaluate_entropy_fix ( const std::array< real, 3 > &  eig_L,
const std::array< real, 3 > &  eig_R,
std::array< real, 3 > &  eig_RoeAvg,
const real  vel2_ravg,
const real  sound_ravg 
) const
virtual

(1) Van Leer et al. (1989 Sonic) entropy fix for acoustic waves (i.e. i=1,5) (2) For waves (i=2,3,4) –> Entropy fix of Liou (2000 Mass) — See p.74 of Osswald et al. (2016 L2Roe)

Implements PHiLiP::NumericalFlux::RoeBaseRiemannSolverDissipation< dim, nstate, real >.

Definition at line 272 of file convective_numerical_flux.cpp.

◆ evaluate_shock_indicator()

template<int dim, int nstate, typename real >
void PHiLiP::NumericalFlux::L2RoeRiemannSolverDissipation< dim, nstate, real >::evaluate_shock_indicator ( const std::array< real, 3 > &  eig_L,
const std::array< real, 3 > &  eig_R,
int &  ssw_LEFT,
int &  ssw_RIGHT 
) const
protected

Shock indicator of Wada & Liou (1994 Flux) – Eq.(39) — See also p.74 of Osswald et al. (2016 L2Roe)

Definition at line 248 of file convective_numerical_flux.cpp.


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