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

Enumerations

enum class  PrimVariable { Pressure , InternalEnergy }
 

Functions

DNDS_DEVICE_CALLABLE void IdealGasThermal (real E, real rho, real vSqr, real gamma, real &p, real &asqr, real &H)
 Compute pressure, speed-of-sound squared, and specific enthalpy from total energy, density, and velocity squared.
 
DNDS_DEVICE_CALLABLE real Pressure_From_InternalEnergy (real e, real gamma)
 Pressure from internal energy: p = (gamma - 1) * e.
 
DNDS_DEVICE_CALLABLE real InternalEnergy_From_Pressure (real p, real gamma)
 Internal energy from pressure: e = p / (gamma - 1)
 
DNDS_DEVICE_CALLABLE real Enthalpy (real E, real rho, real p)
 Specific enthalpy from conservative state: H = (E + p) / rho.
 
DNDS_DEVICE_CALLABLE real SpeedOfSoundSqr (real gamma, real p, real rho)
 Speed of sound squared: a^2 = gamma * p / rho.
 
template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real Cons2PrimEnergy (real E, real rho, real vSqr, real gamma)
 Convert conservative energy to primitive energy-index value.
 
template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real Prim2ConsEnergy (real primE, real rho, real vSqr, real gamma)
 Convert primitive energy-index value to conservative total energy.
 
template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real PrimE2Pressure (real primE, real gamma)
 Get pressure from the primitive energy-index value.
 
DNDS_DEVICE_CALLABLE real RoeSpeedOfSoundSqr (real gamma, real HRoe, real vsqrRoe)
 Roe-averaged speed of sound squared: a^2 = (gamma-1)(H - 0.5*v^2).
 
DNDS_DEVICE_CALLABLE void RoeAlphaDecomposition (real incU0, real incU123N, real incU4b, real veloRoeN, real HRoe, real asqrRoe, real aRoe, real gamma, real &alpha0, real &alpha1, real &alpha4)
 Roe alpha-decomposition coefficients for the 1D wave structure.
 
DNDS_DEVICE_CALLABLE void EntropyFix_HCorrHY (real aL, real aR, real vnL, real vnR, real dLambda, real fixScale, real &lam0, real &lam123, real &lam4)
 H-correction + Harten-Yee entropy fix (scheme 8 in Euler module).
 

Enumeration Type Documentation

◆ PrimVariable

enum class DNDS::IdealGas::PrimVariable
strong

Which thermodynamic variable is stored at the energy index of the primitive state vector.

Enumerator
Pressure 

prim[I4] = p (Euler convention)

InternalEnergy 

prim[I4] = rho * e_internal (EulerP convention)

Definition at line 23 of file IdealGasPhysics.hpp.

Function Documentation

◆ Cons2PrimEnergy()

template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real DNDS::IdealGas::Cons2PrimEnergy ( real  E,
real  rho,
real  vSqr,
real  gamma 
)
inline

Convert conservative energy to primitive energy-index value.

Template Parameters
primWhether to store pressure or internal energy.
Parameters
ETotal energy (conservative).
rhoDensity.
vSqrVelocity squared.
gammaRatio of specific heats.
Returns
The value to store at prim[I4].

Definition at line 90 of file IdealGasPhysics.hpp.

◆ Enthalpy()

DNDS_DEVICE_CALLABLE real DNDS::IdealGas::Enthalpy ( real  E,
real  rho,
real  p 
)
inline

Specific enthalpy from conservative state: H = (E + p) / rho.

Definition at line 62 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ EntropyFix_HCorrHY()

DNDS_DEVICE_CALLABLE void DNDS::IdealGas::EntropyFix_HCorrHY ( real  aL,
real  aR,
real  vnL,
real  vnR,
real  dLambda,
real  fixScale,
real lam0,
real lam123,
real lam4 
)
inline

H-correction + Harten-Yee entropy fix (scheme 8 in Euler module).

This is the entropy fix scheme used by both Euler and EulerP.

Definition at line 194 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ IdealGasThermal()

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

Compute pressure, speed-of-sound squared, and specific enthalpy from total energy, density, and velocity squared.

Definition at line 38 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ InternalEnergy_From_Pressure()

DNDS_DEVICE_CALLABLE real DNDS::IdealGas::InternalEnergy_From_Pressure ( real  p,
real  gamma 
)
inline

Internal energy from pressure: e = p / (gamma - 1)

Definition at line 55 of file IdealGasPhysics.hpp.

◆ Pressure_From_InternalEnergy()

DNDS_DEVICE_CALLABLE real DNDS::IdealGas::Pressure_From_InternalEnergy ( real  e,
real  gamma 
)
inline

Pressure from internal energy: p = (gamma - 1) * e.

Definition at line 48 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ Prim2ConsEnergy()

template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real DNDS::IdealGas::Prim2ConsEnergy ( real  primE,
real  rho,
real  vSqr,
real  gamma 
)
inline

Convert primitive energy-index value to conservative total energy.

Template Parameters
primWhether prim[I4] stores pressure or internal energy.
Parameters
primEThe primitive energy-index value (p or e).
rhoDensity.
vSqrVelocity squared.
gammaRatio of specific heats.
Returns
Total energy E.

Definition at line 111 of file IdealGasPhysics.hpp.

◆ PrimE2Pressure()

template<PrimVariable prim>
DNDS_DEVICE_CALLABLE real DNDS::IdealGas::PrimE2Pressure ( real  primE,
real  gamma 
)
inline

Get pressure from the primitive energy-index value.

Template Parameters
primWhether prim[I4] stores pressure or internal energy.
Parameters
primEThe primitive energy-index value.
gammaRatio of specific heats.

Definition at line 128 of file IdealGasPhysics.hpp.

◆ RoeAlphaDecomposition()

DNDS_DEVICE_CALLABLE void DNDS::IdealGas::RoeAlphaDecomposition ( real  incU0,
real  incU123N,
real  incU4b,
real  veloRoeN,
real  HRoe,
real  asqrRoe,
real  aRoe,
real  gamma,
real alpha0,
real alpha1,
real alpha4 
)
inline

Roe alpha-decomposition coefficients for the 1D wave structure.

Computes the wave strengths alpha0, alpha1, alpha4 from the jump across the interface.

Parameters
incU0Jump in density (UR(0) - UL(0)).
incU123NJump in normal momentum.
incU4bJump in (E - alpha23VT . veloRoe).
veloRoeNRoe-averaged normal velocity.
HRoeRoe-averaged enthalpy.
asqrRoeRoe-averaged speed of sound squared.
aRoeRoe-averaged speed of sound.
gammaRatio of specific heats.
[out]alpha0Left-going acoustic wave strength.
[out]alpha1Entropy wave strength.
[out]alpha4Right-going acoustic wave strength.

Definition at line 168 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ RoeSpeedOfSoundSqr()

DNDS_DEVICE_CALLABLE real DNDS::IdealGas::RoeSpeedOfSoundSqr ( real  gamma,
real  HRoe,
real  vsqrRoe 
)
inline

Roe-averaged speed of sound squared: a^2 = (gamma-1)(H - 0.5*v^2).

Definition at line 144 of file IdealGasPhysics.hpp.

Here is the caller graph for this function:

◆ SpeedOfSoundSqr()

DNDS_DEVICE_CALLABLE real DNDS::IdealGas::SpeedOfSoundSqr ( real  gamma,
real  p,
real  rho 
)
inline

Speed of sound squared: a^2 = gamma * p / rho.

Definition at line 69 of file IdealGasPhysics.hpp.

Here is the caller graph for this function: