[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Base class of Roe (Roe-Pike) flux with entropy fix. Derived from RiemannSolverDissipation. More...
#include <convective_numerical_flux.hpp>
Public Member Functions | |
RoeBaseRiemannSolverDissipation (std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input) | |
Constructor. | |
virtual 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 =0 |
Virtual member function for evaluating the entropy fix for a Roe-Pike flux. | |
virtual 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 =0 |
Virtual member function for evaluating additional modifications/corrections for a Roe-Pike flux. | |
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 |
![]() | |
virtual | ~RiemannSolverDissipation ()=default |
< Base class destructor required for abstract classes. | |
Protected Attributes | |
const std::shared_ptr< Physics::Euler< dim, nstate, real > > | euler_physics |
Numerical flux requires physics to evaluate convective eigenvalues. | |
Base class of Roe (Roe-Pike) flux with entropy fix. Derived from RiemannSolverDissipation.
Definition at line 122 of file convective_numerical_flux.hpp.
|
virtual |
Returns the convective flux at an interface — See Blazek 2015, p.103-105 — Note: Modified calculation of alpha_{3,4} to use dVt (jump in tangential velocities); expressions are equivalent.
Implements PHiLiP::NumericalFlux::RiemannSolverDissipation< dim, nstate, real >.
Definition at line 331 of file convective_numerical_flux.cpp.