DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Euler::Gas Namespace Reference

Classes

struct  RoePreamble
 Pre-computed Roe-averaged quantities shared by all Riemann solvers. More...
 

Typedefs

using tVec = Eigen::Vector3d
 Convenience alias for 3-D real vector.
 
using tVec2 = Eigen::Vector2d
 Convenience alias for 2-D real vector.
 

Enumerations

enum  RiemannSolverType {
  UnknownRS = 0 , Roe = 1 , HLLC = 2 , HLLEP = 3 ,
  HLLEP_V1 = 21 , Roe_M1 = 11 , Roe_M2 = 12 , Roe_M3 = 13 ,
  Roe_M4 = 14 , Roe_M5 = 15 , Roe_M6 = 16 , Roe_M7 = 17 ,
  Roe_M8 = 18 , Roe_M9 = 19
}
 Selects the approximate Riemann solver and its entropy-fix variant. More...
 

Functions

 DNDS_DEFINE_ENUM_JSON (RiemannSolverType, { {UnknownRS, "UnknownRS"}, {Roe, "Roe"}, {HLLC, "HLLC"}, {HLLEP, "HLLEP"}, {HLLEP_V1, "HLLEP_V1"}, {Roe_M1, "Roe_M1"}, {Roe_M2, "Roe_M2"}, {Roe_M3, "Roe_M3"}, {Roe_M4, "Roe_M4"}, {Roe_M5, "Roe_M5"}, {Roe_M6, "Roe_M6"}, {Roe_M7, "Roe_M7"}, {Roe_M8, "Roe_M8"}, {Roe_M9, "Roe_M9"}, }) template< int dim
 Fills the right eigenvector matrix for the 1-D Euler system in the x-direction.
 
class TeV inline void EulerGasRightEigenVector (const TVec &velo, real Vsqr, real H, real a, TeV &ReV)
 
template<int dim = 3, class TVec , class TeV >
void EulerGasLeftEigenVector (const TVec &velo, real Vsqr, real H, real a, real gamma, TeV &LeV)
 Fills the left eigenvector matrix (inverse of the right eigenvector matrix) for the 1-D Euler system in the x-direction.
 
void IdealGasThermal (real E, real rho, real vSqr, real gamma, real &p, real &asqr, real &H)
 Thin wrapper delegating to IdealGas::IdealGasThermal.
 
template<int dim = 3, typename TULm , typename TURm , typename TFdumpInfo >
RoePreamble< dim > ComputeRoePreamble (const TULm &ULm, const TURm &URm, real gamma, const TFdumpInfo &dumpInfo)
 Compute Roe-averaged quantities from mean-state L/R vectors.
 
template<int dim = 3, class TCons , class TPrim >
void IdealGasThermalConservative2Primitive (const TCons &U, TPrim &prim, real gamma)
 Converts conservative variables to primitive variables for an ideal gas.
 
template<int dim = 3, class TCons , class TPrim >
void IdealGasThermalPrimitive2Conservative (const TPrim &prim, TCons &U, real gamma)
 Converts primitive variables to conservative variables for an ideal gas.
 
template<int dim = 3, class TPrim >
std::tuple< real, realIdealGasThermalPrimitiveGetP0T0 (const TPrim &prim, real gamma, real rg)
 Computes total (stagnation) pressure p0 and temperature T0 from a primitive state using isentropic relations.
 
template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG >
void GasInviscidFlux (const TU &U, const TVec &velo, const TVecVG &vg, real p, TF &F)
 Computes the inviscid (Euler) flux in the x-direction for a moving grid.
 
template<int dim = 3, typename TU , typename TF , class TVec , class TVecN , class TVecVG >
void GasInviscidFlux_XY (const TU &U, const TVec &velo, const TVecVG &vg, const TVecN &n, real p, TF &F)
 Computes the inviscid flux projected onto an arbitrary face normal n, accounting for grid motion vg.
 
template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG , class TP >
void GasInviscidFlux_Batch (const TU &U, const TVec &velo, const TVecVG &vg, TP &&p, TF &F)
 Batched x-direction inviscid flux for column-major state matrices.
 
template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG , class TVecN , class TP >
void GasInviscidFlux_XY_Batch (const TU &U, const TVec &velo, const TVecVG &vg, const TVecN &n, TP &&p, TF &F)
 Batched face-normal inviscid flux for column-major state matrices.
 
template<int dim = 3, typename TU , class TVec >
void IdealGasUIncrement (const TU &U, const TU &dU, const TVec &velo, real gamma, TVec &dVelo, real &dp)
 Computes velocity and pressure increments from a conservative-state increment, used for the Lax-flux Jacobian computation.
 
template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG >
void GasInviscidFluxFacialIncrement (const TU &U, const TU &dU, const TVec &unitNorm, const TVecVG &velo, const TVec &dVelo, const TVec &vg, real dp, real p, TF &F)
 Computes the increment of the facial inviscid flux from state, velocity, and pressure increments.
 
template<int dim = 3, typename TU >
auto IdealGas_EulerGasRightEigenVector (const TU &U, real gamma)
 Convenience wrapper that computes the right eigenvector matrix directly from a conservative state vector U.
 
template<int dim = 3, typename TU >
auto IdealGas_EulerGasLeftEigenVector (const TU &U, real gamma)
 Convenience wrapper that computes the left eigenvector matrix directly from a conservative state vector U.
 
template<int dim = 3, int type = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void HLLEPFlux_IdealGas (const TUL &UL, const TUR &UR, const TULm &ULm, const TURm &URm, const TVecVG &vg, const TVecN &n, real gamma, TF &F, real dLambda, real fixScale, const TFdumpInfo &dumpInfo, real &lam0, real &lam123, real &lam4)
 HLLEP (HLL with Enhanced Pressure) approximate Riemann solver for an ideal gas on a moving grid.
 
template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void HLLCFlux_IdealGas_HartenYee (const TUL &UL, const TUR &UR, const TULm &ULm, const TURm &URm, const TVecVG &vg, const TVecN &n, real gamma, TF &F, real dLambda, real fixScale, const TFdumpInfo &dumpInfo, real &lam0, real &lam123, real &lam4)
 HLLC (Harten-Lax-van Leer-Contact) approximate Riemann solver for an ideal gas on a moving grid.
 
template<int eigScheme>
void Roe_EntropyFixer (const real aL, const real aR, const real aAve, const real uL, const real uR, const real uAve, const real VL, const real VR, const real VAve, real dLambda, real fixScale, real incFScale, real &lam0, real &lam123, real &lam4)
 Template-dispatched entropy-fix for Roe-type Riemann solvers.
 
template<int dim = 3, int eigScheme = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void RoeFlux_IdealGas_HartenYee (const TUL &UL, const TUR &UR, const TULm &ULm, const TURm &URm, const TVecVG &vg, const TVecN &n, real gamma, TF &F, real dLambda, real fixScale, real incFScale, const TFdumpInfo &dumpInfo, real &lam0, real &lam123, real &lam4)
 Core Roe approximate Riemann solver with a selectable entropy-fix scheme for an ideal gas on a moving grid.
 
template<int dim = 3, typename TUL , typename TUR , typename TVecV , typename TUOut >
void GetRoeAverage (const TUL &UL, const TUR &UR, real gamma, TVecV &veloRoe, real &vsqrRoe, real &aRoe, real &asqrRoe, real &HRoe, TUOut &UOut)
 Computes the Roe-averaged state vector including passive scalars (e.g. RANS turbulence variables).
 
template<int dim = 3, typename TDU , typename TDF , typename TVecV , typename TVecN >
void RoeFluxIncFDiff (const TDU &incU, const TVecN &n, const TVecV &veloRoe, real vsqrRoe, real aRoe, real asqrRoe, real HRoe, real lam0, real lam123, real lam4, real gamma, TDF &incF)
 Computes the dissipation part of the Roe flux |A|·dU, given pre-computed Roe averages and entropy-fixed eigenvalues.
 
template<int dim = 3, int eigScheme = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecVGm , typename TVecN , typename TVecNm , typename TF , typename TFdumpInfo >
void RoeFlux_IdealGas_HartenYee_Batch (const TUL &UL, const TUR &UR, const TULm &ULm, const TURm &URm, const TVecVG &vg, const TVecVGm &vgm, const TVecN &n, const TVecNm &nm, real gamma, TF &F, real dLambda, real fixScale, real incFScale, const TFdumpInfo &dumpInfo, real &lam0, real &lam123, real &lam4)
 Batched Roe flux with selectable entropy fix for column-major state matrices.
 
template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void InviscidFlux_IdealGas_Dispatcher (RiemannSolverType type, TUL &&UL, TUR &&UR, TULm &&ULm, TURm &&URm, TVecVG &&vg, TVecN &&n, real gamma, TF &&F, real dLambda, real fixScale, real incFScale, TFdumpInfo &&dumpInfo, real &lam0, real &lam123, real &lam4)
 Runtime dispatcher from RiemannSolverType to the correct Roe / HLLC / HLLEP template instantiation.
 
template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecVGm , typename TVecN , typename TVecNm , typename TF , typename TFdumpInfo >
void InviscidFlux_IdealGas_Batch_Dispatcher (RiemannSolverType type, TUL &&UL, TUR &&UR, TULm &&ULm, TURm &&URm, TVecVG &&vg, TVecVGm &&vgm, TVecN &&n, TVecNm &&nm, real gamma, TF &&F, real dLambda, real fixScale, real incFScale, TFdumpInfo &dumpInfo, real &lam0, real &lam123, real &lam4)
 Runtime dispatcher for batched Roe flux from RiemannSolverType to the correct RoeFlux_IdealGas_HartenYee_Batch template instantiation.
 
template<int dim = 3, typename TU , typename TGradU , typename TFlux , typename TNorm >
void ViscousFlux_IdealGas (const TU &U, const TGradU &GradUPrim, TNorm norm, bool adiabatic, real gamma, real mu, real mutRatio, bool mutQCRFix, real k, real Cp, TFlux &Flux)
 Computes the viscous (Navier-Stokes) flux projected onto a face normal for an ideal gas.
 
template<int dim = 3, typename TU , typename TGradU , typename TGradUPrim >
void GradientCons2Prim_IdealGas (const TU &U, const TGradU &GradU, TGradUPrim &GradUPrim, real gamma)
 Converts the gradient of conservative variables to the gradient of primitive variables for an ideal gas.
 
template<int dim, typename TU , typename TGradU >
auto GetGradVelo (const TU &U, const TGradU &GradU)
 Extracts the velocity gradient tensor from the conservative-variable gradient using the quotient rule.
 
template<int dim = 3, int scheme = 0, int nVarsFixed, typename TU , typename TUInc >
real IdealGasGetCompressionRatioPressure (const TU &u, const TUInc &uInc, real newrhoEinteralNew)
 Computes the maximum safe compression ratio (scaling factor α ∈ [0,1]) for a conservative-variable increment that keeps internal energy above a prescribed positive floor.
 

Typedef Documentation

◆ tVec

using DNDS::Euler::Gas::tVec = typedef Eigen::Vector3d

Convenience alias for 3-D real vector.

Definition at line 30 of file Gas.hpp.

◆ tVec2

using DNDS::Euler::Gas::tVec2 = typedef Eigen::Vector2d

Convenience alias for 2-D real vector.

Definition at line 31 of file Gas.hpp.

Enumeration Type Documentation

◆ RiemannSolverType

Selects the approximate Riemann solver and its entropy-fix variant.

Enumerator Solver / fix scheme
Roe Standard Roe with Harten-Yee entropy fix (eigScheme 0).
HLLC Harten-Lax-van Leer-Contact (known accuracy issues, see IV).
HLLEP HLL with Enhanced Pressure (type 0).
HLLEP_V1 HLLEP variant 1 (type 1).
Roe_M1 Roe + cLLF entropy fix (eigScheme 1).
Roe_M2 Roe + Lax-Friedrichs / vanilla Lax (eigScheme 2, early exit).
Roe_M3 Roe + LD (low-dissipation) Roe fix (eigScheme 3).
Roe_M4 Roe + ID (intermediate-dissipation) Roe fix (eigScheme 4).
Roe_M5 Roe + LD cLLF fix (eigScheme 5).
Roe_M6 Roe + H-correction only (eigScheme 6).
Roe_M7 Roe + Harten-Yee fix only, no H-correction (eigScheme 7).
Roe_M8 Roe + H-correction + Harten-Yee fix (eigScheme 8).
Roe_M9 Reserved (eigScheme 9, currently asserts false).
Enumerator
UnknownRS 
Roe 
HLLC 
HLLEP 
HLLEP_V1 
Roe_M1 
Roe_M2 
Roe_M3 
Roe_M4 
Roe_M5 
Roe_M6 
Roe_M7 
Roe_M8 
Roe_M9 

Definition at line 61 of file Gas.hpp.

Function Documentation

◆ ComputeRoePreamble()

template<int dim = 3, typename TULm , typename TURm , typename TFdumpInfo >
RoePreamble< dim > DNDS::Euler::Gas::ComputeRoePreamble ( const TULm &  ULm,
const TURm &  URm,
real  gamma,
const TFdumpInfo &  dumpInfo 
)

Compute Roe-averaged quantities from mean-state L/R vectors.

Extracts velocities, calls IdealGasThermal for both sides, then computes the Roe-averaged velocity, enthalpy, and speed of sound.

Parameters
ULmLeft mean-state conservative vector
URmRight mean-state conservative vector
gammaRatio of specific heats
dumpInfoCallable invoked before assertion on invalid asqrRoe
Returns
RoePreamble<dim> with all fields populated

Definition at line 230 of file Gas.hpp.

Here is the call graph for this function:

◆ DNDS_DEFINE_ENUM_JSON()

DNDS::Euler::Gas::DNDS_DEFINE_ENUM_JSON ( RiemannSolverType  ,
{ {UnknownRS, "UnknownRS"}, {Roe, "Roe"}, {HLLC, "HLLC"}, {HLLEP, "HLLEP"}, {HLLEP_V1, "HLLEP_V1"}, {Roe_M1, "Roe_M1"}, {Roe_M2, "Roe_M2"}, {Roe_M3, "Roe_M3"}, {Roe_M4, "Roe_M4"}, {Roe_M5, "Roe_M5"}, {Roe_M6, "Roe_M6"}, {Roe_M7, "Roe_M7"}, {Roe_M8, "Roe_M8"}, {Roe_M9, "Roe_M9"}, }   
)

Fills the right eigenvector matrix for the 1-D Euler system in the x-direction.

The matrix has (dim+2) columns, one per characteristic wave:

  • column 0: left-running acoustic wave (λ = u − a)
  • column 1: entropy wave (λ = u)
  • columns 2..dim: shear waves (λ = u)
  • column dim+1: right-running acoustic wave (λ = u + a)
Warning
ReV must be pre-allocated to size (dim+2)×(dim+2).
Template Parameters
dimSpatial dimension (2 or 3).
TVecVelocity vector type (Eigen expression, length dim).
TeVEigenvector matrix type (Eigen expression, (dim+2)×(dim+2)).
Parameters
veloRoe-averaged (or cell-centred) velocity vector.
VsqrSquared velocity magnitude (= velo.squaredNorm()).
HTotal specific enthalpy.
aSpeed of sound.
[out]ReVRight eigenvector matrix, overwritten on output.

◆ EulerGasLeftEigenVector()

template<int dim = 3, class TVec , class TeV >
void DNDS::Euler::Gas::EulerGasLeftEigenVector ( const TVec &  velo,
real  Vsqr,
real  H,
real  a,
real  gamma,
TeV &  LeV 
)
inline

Fills the left eigenvector matrix (inverse of the right eigenvector matrix) for the 1-D Euler system in the x-direction.

Satisfies LeV * ReV = I for the (dim+2)×(dim+2) system.

Warning
LeV must be pre-allocated to size (dim+2)×(dim+2).
Template Parameters
dimSpatial dimension (2 or 3).
TVecVelocity vector type.
TeVEigenvector matrix type.
Parameters
veloRoe-averaged (or cell-centred) velocity vector.
VsqrSquared velocity magnitude.
HTotal specific enthalpy.
aSpeed of sound.
gammaRatio of specific heats (γ).
[out]LeVLeft eigenvector matrix, overwritten on output.

Definition at line 160 of file Gas.hpp.

◆ EulerGasRightEigenVector()

class TeV inline void DNDS::Euler::Gas::EulerGasRightEigenVector ( const TVec &  velo,
real  Vsqr,
real  H,
real  a,
TeV &  ReV 
)

Definition at line 120 of file Gas.hpp.

Here is the call graph for this function:

◆ GasInviscidFlux()

template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG >
void DNDS::Euler::Gas::GasInviscidFlux ( const TU &  U,
const TVec &  velo,
const TVecVG &  vg,
real  p,
TF &  F 
)
inline

Computes the inviscid (Euler) flux in the x-direction for a moving grid.

The flux vector is written to the first (dim+2) entries of F: F = [ρ(u−ug), ρu(u−ug)+p, ρv(u−ug), (ρw(u−ug),) (E+p)u − E·ug] where ug = vg(0) is the grid velocity in the x-direction.

Note
Passive-scalar (RANS) entries beyond index dim+1 are not filled.
Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UConservative state vector.
veloFluid velocity vector.
vgGrid velocity vector (for ALE / moving-mesh formulations).
pPressure.
[out]FFlux vector (first dim+2 entries overwritten).

Definition at line 365 of file Gas.hpp.

Here is the call graph for this function:

◆ GasInviscidFlux_Batch()

template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG , class TP >
void DNDS::Euler::Gas::GasInviscidFlux_Batch ( const TU &  U,
const TVec &  velo,
const TVecVG &  vg,
TP &&  p,
TF &  F 
)
inline

Batched x-direction inviscid flux for column-major state matrices.

Each column of U, velo, vg represents one face in the batch; p is a row-vector of pressures. Output columns of F are filled independently.

Template Parameters
dimSpatial dimension.
Parameters
UConservative states [(dim+2) × nBatch].
veloVelocity vectors [dim × nBatch].
vgGrid velocities [dim × nBatch].
pPressures [1 × nBatch].
[out]FFlux matrix [(dim+2) × nBatch].

Definition at line 414 of file Gas.hpp.

Here is the call graph for this function:

◆ GasInviscidFlux_XY()

template<int dim = 3, typename TU , typename TF , class TVec , class TVecN , class TVecVG >
void DNDS::Euler::Gas::GasInviscidFlux_XY ( const TU &  U,
const TVec &  velo,
const TVecVG &  vg,
const TVecN &  n,
real  p,
TF &  F 
)
inline

Computes the inviscid flux projected onto an arbitrary face normal n, accounting for grid motion vg.

This is the main face-flux function used during spatial discretisation. The normal velocity is V_n = (velo − vg) · n.

Note
Passive-scalar entries beyond index dim+1 are not filled.
Template Parameters
dimSpatial dimension.
Parameters
UConservative state vector.
veloFluid velocity vector.
vgGrid velocity vector.
nOutward face-normal vector (unit or area-weighted).
pPressure.
[out]FFlux vector (first dim+2 entries overwritten).

Definition at line 391 of file Gas.hpp.

Here is the call graph for this function:

◆ GasInviscidFlux_XY_Batch()

template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG , class TVecN , class TP >
void DNDS::Euler::Gas::GasInviscidFlux_XY_Batch ( const TU &  U,
const TVec &  velo,
const TVecVG &  vg,
const TVecN &  n,
TP &&  p,
TF &  F 
)
inline

Batched face-normal inviscid flux for column-major state matrices.

Each column represents one face in the batch. The flux at each column is projected onto the corresponding column of the normal matrix n.

Template Parameters
dimSpatial dimension.
Parameters
UConservative states [(dim+2) × nBatch].
veloVelocity vectors [dim × nBatch].
vgGrid velocities [dim × nBatch].
nFace normals [dim × nBatch].
pPressures [1 × nBatch].
[out]FFlux matrix [(dim+2) × nBatch].

Definition at line 437 of file Gas.hpp.

Here is the call graph for this function:

◆ GasInviscidFluxFacialIncrement()

template<int dim = 3, typename TU , typename TF , class TVec , class TVecVG >
void DNDS::Euler::Gas::GasInviscidFluxFacialIncrement ( const TU &  U,
const TU &  dU,
const TVec &  unitNorm,
const TVecVG &  velo,
const TVec &  dVelo,
const TVec &  vg,
real  dp,
real  p,
TF &  F 
)
inline

Computes the increment of the facial inviscid flux from state, velocity, and pressure increments.

Used in the implicit time-stepping Jacobian (e.g. LU-SGS). Writes the full flux vector F including passive-scalar entries beyond index dim+1 (valid for RANS and extra-equation models).

Note
now this function writes to the whole F, assuming SeqI52Last part is passive scalar this is valid for RANS and EX
Template Parameters
dimSpatial dimension.
Parameters
UConservative state vector.
dUConservative state increment.
unitNormFace unit-normal vector.
veloVelocity vector.
dVeloVelocity increment (from IdealGasUIncrement()).
vgGrid velocity vector.
dpPressure increment (from IdealGasUIncrement()).
pPressure.
[out]FIncremental facial flux (all entries written).

Definition at line 496 of file Gas.hpp.

Here is the call graph for this function:

◆ GetGradVelo()

template<int dim, typename TU , typename TGradU >
auto DNDS::Euler::Gas::GetGradVelo ( const TU &  U,
const TGradU &  GradU 
)

Extracts the velocity gradient tensor from the conservative-variable gradient using the quotient rule.

Returns diffVelo(i,j) = ∂u_j/∂x_i, a [dim × dim] matrix, computed as (ρ · ∂(ρu_j)/∂x_i − (ρu_j) · ∂ρ/∂x_i) / ρ².

Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UConservative state vector.
GradUGradient of conservative variables [dim × nVars].
Returns
Eigen::Matrix<real, dim, dim> Velocity gradient tensor (dU_j/dx_i).

Definition at line 1713 of file Gas.hpp.

Here is the call graph for this function:

◆ GetRoeAverage()

template<int dim = 3, typename TUL , typename TUR , typename TVecV , typename TUOut >
void DNDS::Euler::Gas::GetRoeAverage ( const TUL &  UL,
const TUR &  UR,
real  gamma,
TVecV &  veloRoe,
real vsqrRoe,
real aRoe,
real asqrRoe,
real HRoe,
TUOut &  UOut 
)

Computes the Roe-averaged state vector including passive scalars (e.g. RANS turbulence variables).

The Roe average for the Euler components (density, momentum, energy) is the standard √ρ-weighted mean. Passive scalars beyond index dim+1 are also averaged using the same √ρ weighting.

Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UL,URLeft/right conservative state vectors (full length including any passive scalars).
gammaRatio of specific heats (γ).
[out]veloRoeRoe-averaged velocity vector.
[out]vsqrRoeRoe-averaged squared velocity magnitude.
[out]aRoeRoe-averaged speed of sound.
[out]asqrRoeRoe-averaged speed of sound squared.
[out]HRoeRoe-averaged total specific enthalpy.
[out]UOutRoe-averaged state vector (same layout as UL/UR, including passive scalars).

Definition at line 1149 of file Gas.hpp.

Here is the call graph for this function:

◆ GradientCons2Prim_IdealGas()

template<int dim = 3, typename TU , typename TGradU , typename TGradUPrim >
void DNDS::Euler::Gas::GradientCons2Prim_IdealGas ( const TU &  U,
const TGradU &  GradU,
TGradUPrim &  GradUPrim,
real  gamma 
)

Converts the gradient of conservative variables to the gradient of primitive variables for an ideal gas.

Given GradU = ∂U/∂x (the spatial gradient of the conservative vector), computes GradUPrim = ∂(ρ,u,v,[w,]p,...)/∂x in-place. The density gradient column (column 0) is unchanged, velocity gradient columns are transformed via the quotient rule, and the pressure gradient column uses the chain rule through the equation of state.

Note
The input GradU has size [dim × nVars], where nVars ≥ dim+2. Passive-scalar gradient columns beyond index dim+1 are also converted (divided by ρ after removing the density contribution).
Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UConservative state vector.
GradUGradient of conservative variables [dim × nVars].
[out]GradUPrimGradient of primitive variables [dim × nVars].
gammaRatio of specific heats (γ).

Definition at line 1678 of file Gas.hpp.

Here is the call graph for this function:

◆ HLLCFlux_IdealGas_HartenYee()

template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::HLLCFlux_IdealGas_HartenYee ( const TUL &  UL,
const TUR &  UR,
const TULm &  ULm,
const TURm &  URm,
const TVecVG &  vg,
const TVecN &  n,
real  gamma,
TF &  F,
real  dLambda,
real  fixScale,
const TFdumpInfo &  dumpInfo,
real lam0,
real lam123,
real lam4 
)

HLLC (Harten-Lax-van Leer-Contact) approximate Riemann solver for an ideal gas on a moving grid.

Warning
This solver has known accuracy issues (observed in the IV test problem). Prefer Roe or HLLEP for production runs.
Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UL,URLeft/right conservative face states.
ULm,URmLeft/right mean states for Roe averaging.
vgGrid velocity vector.
nFace-normal vector.
gammaRatio of specific heats (γ).
[out]FNumerical flux (first dim+2 entries written).
dLambdaH-correction transverse eigenvalue estimate.
fixScaleEntropy-fix scaling factor.
dumpInfoDebug dump callable for assertion failures.
[out]lam0|V_n − a| (Roe-averaged, before HLLC wave selection).
[out]lam123|V_n| (Roe-averaged).
[out]lam4|V_n + a| (Roe-averaged).

warning: has accuracy issue (see IV test)

is this right for moving mesh?

Definition at line 740 of file Gas.hpp.

Here is the call graph for this function:

◆ HLLEPFlux_IdealGas()

template<int dim = 3, int type = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::HLLEPFlux_IdealGas ( const TUL &  UL,
const TUR &  UR,
const TULm &  ULm,
const TURm &  URm,
const TVecVG &  vg,
const TVecN &  n,
real  gamma,
TF &  F,
real  dLambda,
real  fixScale,
const TFdumpInfo &  dumpInfo,
real lam0,
real lam123,
real lam4 
)

HLLEP (HLL with Enhanced Pressure) approximate Riemann solver for an ideal gas on a moving grid.

Uses Roe-decomposed wave strengths combined with HLL wave-speed estimates. When type = 0, computes the standard HLLEP flux; when type = 1, computes the HLLEP_V1 variant which uses modified eigenvalue scaling factors δ1, δ2, δ3 instead of the HLL blending formula.

Left/right face states (UL, UR) are used for the actual flux and jump, while the mean states (ULm, URm) are used for Roe averaging (these may differ under MUSCL or WENO reconstruction).

Template Parameters
dimSpatial dimension (2 or 3).
type0 = HLLEP, 1 = HLLEP_V1.
Parameters
UL,URLeft/right conservative face states.
ULm,URmLeft/right mean-state conservative vectors for Roe averaging.
vgGrid velocity vector.
nFace-normal vector (unit or area-weighted).
gammaRatio of specific heats (γ).
[out]FNumerical flux (first dim+2 entries written).
dLambdaH-correction transverse eigenvalue estimate from neighbouring faces.
fixScaleUser-configurable scaling factor for entropy fix (from settings.rsFixScale).
dumpInfoCallable invoked before assertion on invalid state (negative density, NaN, etc.).
[out]lam0Absolute eigenvalue |V_n − a| after entropy fixing.
[out]lam123Absolute eigenvalue |V_n| after entropy fixing.
[out]lam4Absolute eigenvalue |V_n + a| after entropy fixing.

not using m, for this is accuracy-limited!

Definition at line 605 of file Gas.hpp.

Here is the call graph for this function:

◆ IdealGas_EulerGasLeftEigenVector()

template<int dim = 3, typename TU >
auto DNDS::Euler::Gas::IdealGas_EulerGasLeftEigenVector ( const TU &  U,
real  gamma 
)
inline

Convenience wrapper that computes the left eigenvector matrix directly from a conservative state vector U.

Extracts velocity and thermodynamic quantities from U, then delegates to EulerGasLeftEigenVector().

Template Parameters
dimSpatial dimension.
Parameters
UConservative state vector.
gammaRatio of specific heats (γ).
Returns
Eigen::Matrix<real, dim+2, dim+2> Left eigenvector matrix.

Definition at line 555 of file Gas.hpp.

Here is the call graph for this function:

◆ IdealGas_EulerGasRightEigenVector()

template<int dim = 3, typename TU >
auto DNDS::Euler::Gas::IdealGas_EulerGasRightEigenVector ( const TU &  U,
real  gamma 
)
inline

Convenience wrapper that computes the right eigenvector matrix directly from a conservative state vector U.

Extracts velocity and thermodynamic quantities from U, then delegates to EulerGasRightEigenVector().

Template Parameters
dimSpatial dimension.
Parameters
UConservative state vector.
gammaRatio of specific heats (γ).
Returns
Eigen::Matrix<real, dim+2, dim+2> Right eigenvector matrix.

Definition at line 528 of file Gas.hpp.

Here is the call graph for this function:

◆ IdealGasGetCompressionRatioPressure()

template<int dim = 3, int scheme = 0, int nVarsFixed, typename TU , typename TUInc >
real DNDS::Euler::Gas::IdealGasGetCompressionRatioPressure ( const TU &  u,
const TUInc &  uInc,
real  newrhoEinteralNew 
)

Computes the maximum safe compression ratio (scaling factor α ∈ [0,1]) for a conservative-variable increment that keeps internal energy above a prescribed positive floor.

Given a state u and an increment uInc, finds the largest α such that the internal energy e = E − ½ρ|v|² of (u + α·uInc) remains above newrhoEinteralNew. Two schemes are available:

  • scheme 0: analytic quadratic solve with iterative safety fallback.
  • scheme 1: convex (linear) estimation only.

Used during implicit time stepping and limiting to prevent negative pressures.

TODO: vectorize

Template Parameters
dimSpatial dimension (2 or 3).
scheme0 = quadratic solve, 1 = convex estimate.
nVarsFixedCompile-time size of the state vector.
Parameters
uCurrent conservative state.
uIncProposed conservative state increment.
newrhoEinteralNewDesired minimum internal energy floor ρ·e_min (= p_min / (γ−1)).
Returns
The safe scaling factor α in [0, 1].

using convex estimation

Definition at line 1750 of file Gas.hpp.

Here is the call graph for this function:

◆ IdealGasThermal()

void DNDS::Euler::Gas::IdealGasThermal ( real  E,
real  rho,
real  vSqr,
real  gamma,
real p,
real asqr,
real H 
)
inline

Thin wrapper delegating to IdealGas::IdealGasThermal.

Definition at line 190 of file Gas.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ IdealGasThermalConservative2Primitive()

template<int dim = 3, class TCons , class TPrim >
void DNDS::Euler::Gas::IdealGasThermalConservative2Primitive ( const TCons &  U,
TPrim &  prim,
real  gamma 
)
inline

Converts conservative variables to primitive variables for an ideal gas.

Conservative layout: U = [ρ, ρu, ρv, (ρw,) E, ...] Primitive layout: prim = [ρ, u, v, (w,) p, ...]

The primitive variable at index I4 = dim+1 stores pressure (Euler module convention), not temperature.

Template Parameters
dimSpatial dimension (2 or 3).
TConsConservative vector type (Eigen expression, length ≥ dim+2).
TPrimPrimitive vector type (Eigen expression, same size as TCons).
Parameters
UConservative state vector (input).
[out]primPrimitive state vector (output).
gammaRatio of specific heats (γ).

Definition at line 279 of file Gas.hpp.

◆ IdealGasThermalPrimitive2Conservative()

template<int dim = 3, class TCons , class TPrim >
void DNDS::Euler::Gas::IdealGasThermalPrimitive2Conservative ( const TPrim &  prim,
TCons &  U,
real  gamma 
)
inline

Converts primitive variables to conservative variables for an ideal gas.

Inverse of IdealGasThermalConservative2Primitive().

Template Parameters
dimSpatial dimension (2 or 3).
TConsConservative vector type.
TPrimPrimitive vector type.
Parameters
primPrimitive state vector (input).
[out]UConservative state vector (output).
gammaRatio of specific heats (γ).

Definition at line 306 of file Gas.hpp.

◆ IdealGasThermalPrimitiveGetP0T0()

template<int dim = 3, class TPrim >
std::tuple< real, real > DNDS::Euler::Gas::IdealGasThermalPrimitiveGetP0T0 ( const TPrim &  prim,
real  gamma,
real  rg 
)

Computes total (stagnation) pressure p0 and temperature T0 from a primitive state using isentropic relations.

Uses p0 = p * (1 + (γ−1)/2 · M²)^(γ/(γ−1)) and T0 = T * (1 + (γ−1)/2 · M²).

Template Parameters
dimSpatial dimension (2 or 3).
TPrimPrimitive vector type.
Parameters
primPrimitive state [ρ, u, v, (w,) p, ...].
gammaRatio of specific heats (γ).
rgSpecific gas constant R_gas = Cp − Cv.
Returns
std::tuple<real,real> (p0, T0).

Definition at line 333 of file Gas.hpp.

◆ IdealGasUIncrement()

template<int dim = 3, typename TU , class TVec >
void DNDS::Euler::Gas::IdealGasUIncrement ( const TU &  U,
const TU &  dU,
const TVec &  velo,
real  gamma,
TVec &  dVelo,
real dp 
)
inline

Computes velocity and pressure increments from a conservative-state increment, used for the Lax-flux Jacobian computation.

Given a state U and its increment dU, computes dVelo = d(momentum/ρ) and dp = (γ−1)(dE − ½(dMomentum·v + momentum·dv)).

Template Parameters
dimSpatial dimension.
Parameters
UConservative state vector.
dUConservative state increment.
veloVelocity vector (= U[1..dim] / U[0]).
gammaRatio of specific heats (γ).
[out]dVeloVelocity increment (output).
[out]dpPressure increment (output).

Definition at line 463 of file Gas.hpp.

◆ InviscidFlux_IdealGas_Batch_Dispatcher()

template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecVGm , typename TVecN , typename TVecNm , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::InviscidFlux_IdealGas_Batch_Dispatcher ( RiemannSolverType  type,
TUL &&  UL,
TUR &&  UR,
TULm &&  ULm,
TURm &&  URm,
TVecVG &&  vg,
TVecVGm &&  vgm,
TVecN &&  n,
TVecNm &&  nm,
real  gamma,
TF &&  F,
real  dLambda,
real  fixScale,
real  incFScale,
TFdumpInfo &  dumpInfo,
real lam0,
real lam123,
real lam4 
)

Runtime dispatcher for batched Roe flux from RiemannSolverType to the correct RoeFlux_IdealGas_HartenYee_Batch template instantiation.

Only Roe variants (Roe, Roe_M1..M9) are supported in batch mode; HLLC and HLLEP are not available and will trigger an assertion failure.

Template Parameters
dimSpatial dimension (2 or 3).
Parameters
typeRiemann solver selector (runtime).
UL,URLeft/right face states [(dim+2) × nBatch].
ULm,URmMean states for Roe averaging (single column).
vg,vgmPer-face and mean-state grid velocities.
n,nmPer-face and mean-state face normals.
gammaRatio of specific heats (γ).
[out]FNumerical flux matrix [(dim+2) × nBatch].
dLambdaH-correction transverse eigenvalue estimate.
fixScaleEntropy-fix scaling factor.
incFScaleIncremental flux dissipation scale.
dumpInfoDebug dump callable.
[out]lam0|V_n − a| after entropy fixing.
[out]lam123|V_n| after entropy fixing.
[out]lam4|V_n + a| after entropy fixing.

Definition at line 1520 of file Gas.hpp.

Here is the call graph for this function:

◆ InviscidFlux_IdealGas_Dispatcher()

template<int dim = 3, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::InviscidFlux_IdealGas_Dispatcher ( RiemannSolverType  type,
TUL &&  UL,
TUR &&  UR,
TULm &&  ULm,
TURm &&  URm,
TVecVG &&  vg,
TVecN &&  n,
real  gamma,
TF &&  F,
real  dLambda,
real  fixScale,
real  incFScale,
TFdumpInfo &&  dumpInfo,
real lam0,
real lam123,
real lam4 
)

Runtime dispatcher from RiemannSolverType to the correct Roe / HLLC / HLLEP template instantiation.

Maps RiemannSolverType enumerators to the corresponding compile-time eigScheme / type template arguments. For Roe variants, calls RoeFlux_IdealGas_HartenYee; for HLLEP, calls HLLEPFlux_IdealGas; for HLLC, calls HLLCFlux_IdealGas_HartenYee.

Template Parameters
dimSpatial dimension (2 or 3).
Parameters
typeRiemann solver selector (runtime).
UL,URLeft/right conservative face states.
ULm,URmMean states for Roe averaging.
vgGrid velocity.
nFace-normal.
gammaRatio of specific heats (γ).
[out]FNumerical flux.
dLambdaH-correction transverse eigenvalue estimate.
fixScaleEntropy-fix scaling factor.
incFScaleIncremental flux dissipation scale.
dumpInfoDebug dump callable.
[out]lam0|V_n − a| after entropy fixing.
[out]lam123|V_n| after entropy fixing.
[out]lam4|V_n + a| after entropy fixing.

Definition at line 1439 of file Gas.hpp.

Here is the call graph for this function:

◆ Roe_EntropyFixer()

template<int eigScheme>
void DNDS::Euler::Gas::Roe_EntropyFixer ( const real  aL,
const real  aR,
const real  aAve,
const real  uL,
const real  uR,
const real  uAve,
const real  VL,
const real  VR,
const real  VAve,
real  dLambda,
real  fixScale,
real  incFScale,
real lam0,
real lam123,
real lam4 
)

Template-dispatched entropy-fix for Roe-type Riemann solvers.

Modifies the absolute Roe eigenvalues (lam0, lam123, lam4) in-place to prevent the sonic glitch and control numerical dissipation. The fix strategy is selected at compile time via eigScheme:

eigScheme Strategy
0 Harten-Yee (H2) — threshold = max(0, λ−λ_L, λ_R−λ) · s.
1 cLLF — componentwise local Lax-Friedrichs.
3 LD Roe — low-dissipation (Fleischmann et al. 2020).
4 ID Roe — intermediate-dissipation (Fleischmann et al. 2020).
5 LD cLLF — low-dissipation componentwise LLF.
6 H-correction only (floor by dLambda · scaleHFix).
7 Harten-Yee only, no H-correction.
8 H-correction + Harten-Yee combined.
Template Parameters
eigSchemeCompile-time entropy-fix scheme selector (0–8).
Parameters
aL,aRLeft/right mean-state speed of sound.
aAveRoe-averaged speed of sound.
uL,uRLeft/right normal velocities relative to grid.
uAveRoe-averaged normal velocity relative to grid.
VL,VRLeft/right total velocity magnitudes relative to grid.
VAveRoe-averaged total velocity magnitude relative to grid.
dLambdaH-correction transverse eigenvalue estimate.
fixScaleUser scaling factor applied to the base entropy-fix constants (kScaleHartenYee, kScaleLD, kScaleHFix).
incFScaleIncremental flux dissipation scale (settings.rsIncFScale), controls overall dissipation level in the contact/shear wave.
[in,out]lam0|V_n − a| — modified in place.
[in,out]lam123|V_n| — modified in place.
[in,out]lam4|V_n + a| — modified in place.

Nico Fleischmann, Stefan Adami, Xiangyu Y. Hu, Nikolaus A. Adams, A low dissipation method to cure the grid-aligned shock instability, 2020

Nico Fleischmann, Stefan Adami, Xiangyu Y. Hu, Nikolaus A. Adams, A low dissipation method to cure the grid-aligned shock instability, 2020

why signM here?

Nico Fleischmann, Stefan Adami, Xiangyu Y. Hu, Nikolaus A. Adams, A low dissipation method to cure the grid-aligned shock instability, 2020

Definition at line 879 of file Gas.hpp.

Here is the call graph for this function:

◆ RoeFlux_IdealGas_HartenYee()

template<int dim = 3, int eigScheme = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecN , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::RoeFlux_IdealGas_HartenYee ( const TUL &  UL,
const TUR &  UR,
const TULm &  ULm,
const TURm &  URm,
const TVecVG &  vg,
const TVecN &  n,
real  gamma,
TF &  F,
real  dLambda,
real  fixScale,
real  incFScale,
const TFdumpInfo &  dumpInfo,
real lam0,
real lam123,
real lam4 
)

Core Roe approximate Riemann solver with a selectable entropy-fix scheme for an ideal gas on a moving grid.

Computes the numerical inter-cell flux as F = ½(F_L + F_R) − ½ |A_Roe| (U_R − U_L) where |A_Roe| is the Roe matrix with eigenvalues modified by the chosen entropy fix (see Roe_EntropyFixer()).

When eigScheme = 2 (Lax-Friedrichs), the function takes an early-exit path returning the vanilla Lax-Friedrichs flux without Roe decomposition.

Template Parameters
dimSpatial dimension (2 or 3).
eigSchemeCompile-time entropy-fix selector (0–8); default 0 (Harten-Yee).
Parameters
UL,URLeft/right conservative face states.
ULm,URmLeft/right mean states for Roe averaging.
vgGrid velocity vector.
nFace-normal vector.
gammaRatio of specific heats (γ).
[out]FNumerical flux (first dim+2 entries written).
dLambdaH-correction transverse eigenvalue estimate.
fixScaleEntropy-fix scaling factor (settings.rsFixScale).
incFScaleIncremental flux dissipation scale (settings.rsIncFScale).
dumpInfoDebug dump callable for assertion failures.
[out]lam0|V_n − a| after entropy fixing.
[out]lam123|V_n| after entropy fixing.
[out]lam4|V_n + a| after entropy fixing.

not using m, for this is accuracy-limited!

Definition at line 1043 of file Gas.hpp.

Here is the call graph for this function:

◆ RoeFlux_IdealGas_HartenYee_Batch()

template<int dim = 3, int eigScheme = 0, typename TUL , typename TUR , typename TULm , typename TURm , typename TVecVG , typename TVecVGm , typename TVecN , typename TVecNm , typename TF , typename TFdumpInfo >
void DNDS::Euler::Gas::RoeFlux_IdealGas_HartenYee_Batch ( const TUL &  UL,
const TUR &  UR,
const TULm &  ULm,
const TURm &  URm,
const TVecVG &  vg,
const TVecVGm &  vgm,
const TVecN &  n,
const TVecNm &  nm,
real  gamma,
TF &  F,
real  dLambda,
real  fixScale,
real  incFScale,
const TFdumpInfo &  dumpInfo,
real lam0,
real lam123,
real lam4 
)

Batched Roe flux with selectable entropy fix for column-major state matrices.

Evaluates the Roe flux for nB = UL.cols() faces simultaneously. The Roe averaging is performed from the single-column mean states ULm, URm, sharing one set of Roe eigenvalues across the batch. Per-face normal vectors n and grid velocities vg allow each column to have a distinct orientation and motion; nm and vgm are the mean-state normal/grid velocity used for the shared Roe average.

Template Parameters
dimSpatial dimension (2 or 3).
eigSchemeCompile-time entropy-fix selector (0–8); default 0.
Parameters
UL,URLeft/right face states [(dim+2) × nBatch].
ULm,URmMean-state vectors for Roe averaging (single column).
vg,vgmPer-face and mean-state grid velocities.
n,nmPer-face and mean-state face normals.
gammaRatio of specific heats (γ).
[out]FNumerical flux matrix [(dim+2) × nBatch].
dLambdaH-correction transverse eigenvalue estimate.
fixScaleEntropy-fix scaling factor.
incFScaleIncremental flux dissipation scale.
dumpInfoDebug dump callable.
[out]lam0|V_n − a| after entropy fixing.
[out]lam123|V_n| after entropy fixing.
[out]lam4|V_n + a| after entropy fixing.

not using m, for this is accuracy-limited!

Definition at line 1278 of file Gas.hpp.

Here is the call graph for this function:

◆ RoeFluxIncFDiff()

template<int dim = 3, typename TDU , typename TDF , typename TVecV , typename TVecN >
void DNDS::Euler::Gas::RoeFluxIncFDiff ( const TDU &  incU,
const TVecN &  n,
const TVecV &  veloRoe,
real  vsqrRoe,
real  aRoe,
real  asqrRoe,
real  HRoe,
real  lam0,
real  lam123,
real  lam4,
real  gamma,
TDF &  incF 
)

Computes the dissipation part of the Roe flux |A|·dU, given pre-computed Roe averages and entropy-fixed eigenvalues.

This function accumulates into incF (does not zero it first). It handles the Euler components (indices 0..dim+1) via the standard Roe eigenvector decomposition, and passive-scalar components (indices dim+2..end) with the contact/shear eigenvalue lam123.

Used by the LUSGS implicit solver where the dissipation term is needed separately from the central flux.

Template Parameters
dimSpatial dimension (2 or 3).
Parameters
incUState jump U_R − U_L (full length, including scalars).
nFace-normal vector.
veloRoeRoe-averaged velocity.
vsqrRoeRoe-averaged |v|².
aRoeRoe-averaged speed of sound.
asqrRoeRoe-averaged a².
HRoeRoe-averaged total enthalpy.
lam0Entropy-fixed eigenvalue |V_n − a|.
lam123Entropy-fixed eigenvalue |V_n|.
lam4Entropy-fixed eigenvalue |V_n + a|.
gammaRatio of specific heats (γ).
[in,out]incFDissipation flux; accumulated (not zeroed).

Definition at line 1210 of file Gas.hpp.

Here is the call graph for this function:

◆ ViscousFlux_IdealGas()

template<int dim = 3, typename TU , typename TGradU , typename TFlux , typename TNorm >
void DNDS::Euler::Gas::ViscousFlux_IdealGas ( const TU &  U,
const TGradU &  GradUPrim,
TNorm  norm,
bool  adiabatic,
real  gamma,
real  mu,
real  mutRatio,
bool  mutQCRFix,
real  k,
real  Cp,
TFlux &  Flux 
)

Computes the viscous (Navier-Stokes) flux projected onto a face normal for an ideal gas.

Evaluates the viscous stress tensor τ (including Stokes' hypothesis λ = −2/3 μ) and the heat flux vector q = −k ∇T, then projects onto norm. Supports:

  • Sutherland-law viscosity (via the externally provided mu).
  • Turbulent eddy viscosity through mutRatio (μ_t / μ).
  • QCR (Quadratic Constitutive Relation) correction when mutQCRFix is true, using ccr1 = 0.3.
  • Adiabatic wall treatment: when adiabatic is true, the wall-normal component of ∇T is removed.
Note
GradUPrim has size [dim × nVars]; it is the gradient of the primitive variables (ρ, u, v, [w,] p, ...) in physical space. "3×5 TGradU" in the original comment refers to dim=3 with 5 Euler variables.
Template Parameters
dimSpatial dimension (2 or 3).
Parameters
UConservative state vector.
GradUPrimGradient of primitive variables [dim × nVars].
normFace-normal vector (unit or area-weighted).
adiabaticIf true, zero the wall-normal heat flux component.
gammaRatio of specific heats (γ).
muDynamic (molecular + turbulent) viscosity μ + μ_t.
mutRatioTurbulent-to-molecular viscosity ratio μ_t / μ.
mutQCRFixEnable QCR stress correction.
kThermal conductivity (molecular + turbulent).
CpSpecific heat at constant pressure.
[out]FluxViscous flux vector (dim+2 entries, Flux[0] = 0).

is this fix reasonable?

Definition at line 1595 of file Gas.hpp.

Here is the call graph for this function: