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

RANS turbulence model functions (eddy viscosity, viscous flux, source terms). More...

Namespaces

namespace  SA
 Spalart-Allmaras model constants used across multiple files.
 

Functions

template<int dim, class TU , class TDiffU >
real GetMut_RealizableKe (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d)
 Compute turbulent viscosity mu_t from the Shih et al. realizable k-epsilon model.
 
template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void GetVisFlux_RealizableKe (TU &&UMeanXy, TDiffU &&DiffUxyPrim, TN &&uNorm, real mut, real d, real muPhy, TVFlux &vFlux)
 Compute the RANS viscous flux for the realizable k-epsilon transport equations.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSource_RealizableKe (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d, TSource &source, int mode)
 Compute source terms for the realizable k-epsilon transport equations.
 
template<int dim, class TU >
std::tuple< real, realSolveZeroGradEquilibrium (TU &u, real muPhy)
 Attempt to find equilibrium epsilon for zero-gradient (freestream) conditions via Newton iteration.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSourceJacobianDiag_RealizableKe (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, TSource &source)
 Compute the analytical diagonal of the source Jacobian for the realizable k-epsilon model.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSourceJacobianDiag_RealizableKe_ND (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, TSource &source)
 Compute the numerical (finite-difference) diagonal of the source Jacobian for the realizable k-epsilon model.
 
template<int dim, class TU , class TDiffU >
real GetMut_SST (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d)
 Compute turbulent viscosity mu_t from Menter's k-omega SST model.
 
template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void GetVisFlux_SST (TU &&UMeanXy, TDiffU &&DiffUxyPrim, TN &&uNorm, real mutIn, real d, real muf, TVFlux &vFlux)
 Compute the RANS viscous flux for the k-omega SST transport equations.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSource_SST (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d, real lLES, TSource &source, int mode)
 Compute source terms for the k-omega SST transport equations.
 
template<int dim, class TU , class TDiffU >
real GetMut_KOWilcox (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d)
 Compute turbulent viscosity mu_t from the Wilcox k-omega model.
 
template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void GetVisFlux_KOWilcox (TU &&UMeanXy, TDiffU &&DiffUxyPrim, TN &&uNorm, real mutIn, real d, real muf, TVFlux &vFlux)
 Compute the RANS viscous flux for the Wilcox k-omega transport equations.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSource_KOWilcox (TU &&UMeanXy, TDiffU &&DiffUxy, real muf, real d, TSource &source, int mode)
 Compute source terms for the Wilcox k-omega transport equations.
 
template<int dim, class TU , class TDiffU , class TSource >
void GetSource_SA (TU &&UMeanXy, TDiffU &&DiffUxy, real muRef, real mufPhy, real gamma, real d, real lLES, real hMax, int DESMode, TSource &source, int rotCor, int mode)
 Compute source terms for the Spalart-Allmaras one-equation turbulence model.
 

Detailed Description

RANS turbulence model functions (eddy viscosity, viscous flux, source terms).

Function Documentation

◆ GetMut_KOWilcox()

template<int dim, class TU , class TDiffU >
real DNDS::Euler::RANS::GetMut_KOWilcox ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d 
)

Compute turbulent viscosity mu_t from the Wilcox k-omega model.

Evaluates the eddy viscosity mu_t = rho * k / omega_tilde, where omega_tilde depends on the selected model version (KW_WILCOX_VER):

  • Version 0 (original 1988): omega_tilde = omega (no stress limiter).
  • Version 1 (2006): omega_tilde = max(omega, C_lim * sqrt(S_{ij}S_{ij} / beta*)), applying the stress limiter with C_lim = 7/8.
  • Version 2 (simplified): omega_tilde = omega (no stress limiter).

Optional upper bound of 1e5 * mu_laminar when KW_WILCOX_LIMIT_MUT is enabled.

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
Parameters
UMeanXyConservative state vector [rho, rhoU, ..., rhoE, rho*k, rho*omega].
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance (unused in this model but kept for interface consistency).
Returns
Turbulent dynamic viscosity mu_t (>= 0).

Definition at line 715 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetMut_RealizableKe()

template<int dim, class TU , class TDiffU >
real DNDS::Euler::RANS::GetMut_RealizableKe ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d 
)

Compute turbulent viscosity mu_t from the Shih et al. realizable k-epsilon model.

Evaluates the eddy viscosity using the realizable k-epsilon formulation with:

  • f_mu damping function based on the turbulent Reynolds number Re_t = rho * k^2 / (mu * epsilon).
  • Realizability constraint: mu_t <= phi * rho * k / S, where S is the strain-rate magnitude.
  • Optional upper bound of 1e5 * mu_laminar when KE_LIMIT_MUT is enabled.
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
Parameters
UMeanXyConservative state vector [rho, rhoU, ..., rhoE, rho*k, rho*epsilon].
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance (unused in this model but kept for interface consistency).
Returns
Turbulent dynamic viscosity mu_t (>= 0).

Definition at line 91 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetMut_SST()

template<int dim, class TU , class TDiffU >
real DNDS::Euler::RANS::GetMut_SST ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d 
)

Compute turbulent viscosity mu_t from Menter's k-omega SST model.

Evaluates the SST eddy viscosity with:

  • Vorticity magnitude Omega = ||(grad U)^T - grad U|| / sqrt(2).
  • Blending function F2 based on wall distance and turbulent Reynolds number.
  • Bradshaw's structural constraint: mu_t = a1 * k / max(Omega * F2, a1 * omega), where a1 = 0.31.
  • Optional upper bound of 1e5 * mu_laminar when KW_SST_LIMIT_MUT is enabled.
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
Parameters
UMeanXyConservative state vector [rho, rhoU, ..., rhoE, rho*k, rho*omega].
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance.
Returns
Turbulent dynamic viscosity mu_t (>= 0).

Definition at line 443 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSource_KOWilcox()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSource_KOWilcox ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d,
TSource &  source,
int  mode 
)

Compute source terms for the Wilcox k-omega transport equations.

Evaluates production, destruction, and cross-diffusion source contributions for the k and omega equations following the Wilcox k-omega model. The version-dependent behavior is controlled by KW_WILCOX_VER:

  • Version 0 (1988): alpha = 5/9, beta = 3/40, no stress limiter, no cross-diffusion.
  • Version 1 (2006): alpha = 13/25, beta = f_beta * beta_0, with the f_beta correction based on Chi_omega (mean rotation / strain invariant), stress limiter C_lim = 7/8, and cross-diffusion sigma_d term (sigma_d = 0.125 when grad(k) . grad(omega) > 0).
  • Version 2 (simplified): alpha = 13/25, beta = 0.075, vorticity-based production, no stress limiter, no cross-diffusion.

Production of k is optionally capped at 20 * beta* * rho * k * omega (CFL3D-style) when KW_WILCOX_PROD_LIMITS is enabled.

When mode == 0, the full source vector is returned:

  • source(I4+1) = P_k - beta* * rho * k * omega
  • source(I4+2) = P_omega - beta * rho * omega^2 + sigma_d / omega * grad(k) . grad(omega)

When mode == 1, the implicit diagonal contribution (positive) is returned:

  • source(I4+1) = beta* * omega
  • source(I4+2) = 2 * beta * omega
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output source vector.
Parameters
UMeanXyConservative state vector.
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance (unused here, kept for interface consistency).
[out]sourceOutput source vector; entries at I4+1 and I4+2 are set.
modeSource computation mode: 0 = full source, 1 = implicit diagonal (destruction only).

Definition at line 850 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSource_RealizableKe()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSource_RealizableKe ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d,
TSource &  source,
int  mode 
)

Compute source terms for the realizable k-epsilon transport equations.

Evaluates production, destruction, and extra source contributions for the k and epsilon equations following the Shih et al. realizable k-epsilon model.

Key terms computed:

  • Reynolds stress tensor tau_{ij} via Boussinesq hypothesis with realizability clipping.
  • Production of k: P_k = -tau_{ij} * dU_i/dx_j, capped by the Shih pphi limiter.
  • Turbulent time scale T_t with low-Re correction via zeta = sqrt(Re_t / 2).
  • Extra dissipation term E involving the gradient of the turbulent time scale.

When mode == 0, the full source vector is returned:

  • source(I4+1) = P_k - rho*epsilon
  • source(I4+2) = (c_{e1} * P_k - c_{e2} * rho*epsilon + E) / T_t

When mode == 1, the implicit diagonal contribution (positive) is returned for use in implicit time stepping:

  • source(I4+1) = 0
  • source(I4+2) = c_{e2} / T_t
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output source vector.
Parameters
UMeanXyConservative state vector.
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance (unused here, kept for interface consistency).
[out]sourceOutput source vector; entries at I4+1 and I4+2 are set.
modeSource computation mode: 0 = full source, 1 = implicit diagonal (destruction only).

Definition at line 192 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSource_SA()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSource_SA ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muRef,
real  mufPhy,
real  gamma,
real  d,
real  lLES,
real  hMax,
int  DESMode,
TSource &  source,
int  rotCor,
int  mode 
)

Compute source terms for the Spalart-Allmaras one-equation turbulence model.

Evaluates the full SA source term including production, destruction, diffusion, and optional DES/DDES/IDDES/WMLES length-scale modification. The SA model solves a single transport equation for the modified kinematic viscosity nuTilde (stored as rho * nuTilde at index I4+1 in the conservative state vector).

Key terms computed:

  • Vorticity magnitude S = ||Omega|| * sqrt(2), with optional rotation correction (SRotCor = c_rot * min(0, S - SS)) when rotCor is enabled.
  • Modified vorticity Sh = S + Sbar, where Sbar = nuTilde * f_{v2} / (kappa^2 * d^2).
  • Production P = c_{b1} * (1 - f_{t2}) * Sh * nuTilde.
  • Destruction D = (c_{w1} * f_w - c_{b1}/kappa^2 * f_{t2}) * (nuTilde / d)^2.
  • Diffusion term: c_{b2}/sigma * |grad(nuTilde)|^2 (conservative form contribution).
  • f_{t2} laminar suppression term controlled by SA_USE_FT2_TERM.

DES/DDES/IDDES length-scale modification:

  • DESMode 0 (RANS): l_DES = d (pure RANS, wall distance).
  • DESMode 1 (IDDES): Blends l_RANS and l_LES using the IDDES shielding functions f_d, f_B, f_e, with psi correction for grid-induced separation avoidance.
  • DESMode 2 (WMLES): Uses the wall-modeled LES length scale l_WMLES = f_B * (1 + f_e) * l_RANS + (1 - f_B) * l_LES * psi.

When mode == 0, the full source is returned at source(I4+1). When mode == 1, the implicit diagonal contribution -dS/dU (positive, suitable for augmenting the implicit operator) is returned at source(I4+1).

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output source vector.
Parameters
UMeanXyConservative state vector [rho, rhoU, ..., rhoE, rho*nuTilde].
DiffUxyGradient tensor of conservative variables, size dim × nVars.
muRefReference viscosity scale used to non-dimensionalise nuTilde (nuTilde_physical = UMeanXy(I4+1) * muRef / rho).
mufPhyPhysical (dimensional) laminar dynamic viscosity.
gammaRatio of specific heats (used for pressure gradient in optional helicity correction).
dWall distance.
lLESLES length scale for DES shielding (set to a large value to disable DES).
hMaxMaximum cell size (used in IDDES alpha parameter: alpha = 0.25 - d/hMax).
DESModeDES operating mode: 0 = pure RANS, 1 = IDDES, 2 = WMLES.
[out]sourceOutput source vector; entry at I4+1 is set.
rotCorRotation correction flag: 0 = disabled, nonzero = enabled (c_rot = 2.0).
modeSource computation mode: 0 = full source, 1 = implicit diagonal (positive).

note that psi has lower bound of 1

super subs

Definition at line 986 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSource_SST()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSource_SST ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
real  d,
real  lLES,
TSource &  source,
int  mode 
)

Compute source terms for the k-omega SST transport equations.

Evaluates production, destruction, and cross-diffusion source contributions for the k and omega equations following Menter's SST model. Supports DES/DDES capability via the RANS-to-LES length scale ratio.

Key terms:

  • Production P_k from the Boussinesq Reynolds stress, optionally capped at 20 * beta* * rho * k * omega (when KW_SST_PROD_LIMITS is enabled).
  • Omega production P_omega with version-dependent formulation controlled by KW_SST_PROD_OMEGA_VERSION (0 = classical, 1 = strain-rate, 2 = vorticity).
  • Destruction terms: beta* * rho * k * omega for k, beta * rho * omega^2 for omega.
  • Cross-diffusion term: 2 * (1 - F1) * rho * sigma_omega2 / omega * grad(k) . grad(omega).
  • DES shielding: F_DES = max(1, l_RANS / l_LES * (1 - F2)), reducing the destruction of k outside the RANS region.

When mode == 0, the full source vector is returned. When mode == 1, the implicit diagonal contribution (positive) is returned:

  • source(I4+1) = beta* * omega * F_DES
  • source(I4+2) = 2 * beta * omega
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output source vector.
Parameters
UMeanXyConservative state vector.
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
dWall distance.
lLESLES length scale for DES/DDES shielding (set to a large value to disable DES).
[out]sourceOutput source vector; entries at I4+1 and I4+2 are set.
modeSource computation mode: 0 = full source, 1 = implicit diagonal (destruction only).

Definition at line 609 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSourceJacobianDiag_RealizableKe()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSourceJacobianDiag_RealizableKe ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
TSource &  source 
)

Compute the analytical diagonal of the source Jacobian for the realizable k-epsilon model.

Evaluates -dS/dU diagonal entries for implicit time stepping of the k and epsilon transport equations. The diagonal entries approximate the linearised destruction terms:

  • source(I4+1) = -P_k / (rho*k) (destruction rate of k)
  • source(I4+2) = c_{e2} / T_t (destruction rate of epsilon)

These positive values are used to augment the implicit operator diagonal, improving stability of the coupled RANS system.

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output Jacobian diagonal vector.
Parameters
UMeanXyConservative state vector.
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
[out]sourceOutput Jacobian diagonal vector; entries at I4+1 and I4+2 are set.

Definition at line 335 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetSourceJacobianDiag_RealizableKe_ND()

template<int dim, class TU , class TDiffU , class TSource >
void DNDS::Euler::RANS::GetSourceJacobianDiag_RealizableKe_ND ( TU &&  UMeanXy,
TDiffU &&  DiffUxy,
real  muf,
TSource &  source 
)

Compute the numerical (finite-difference) diagonal of the source Jacobian for the realizable k-epsilon model.

Evaluates -dS_i/dU_i by forward finite-difference perturbation of rho*k and rho*epsilon independently. The perturbation step is proportional to sqrt(smallReal) times the variable magnitude. Results are clamped to non-negative values and amplified by a factor of 10 for enhanced implicit stability.

This is a fallback when the analytical Jacobian diagonal (GetSourceJacobianDiag_RealizableKe) is suspected of inaccuracy.

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, conservative variables).
TSourceEigen-compatible type for the output Jacobian diagonal vector.
Parameters
UMeanXyConservative state vector.
DiffUxyGradient tensor of conservative variables, size dim × nVars.
mufLaminar (molecular) dynamic viscosity.
[out]sourceOutput numerical Jacobian diagonal vector; entries at I4+1 and I4+2 are set.

Definition at line 400 of file RANS_ke.hpp.

◆ GetVisFlux_KOWilcox()

template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void DNDS::Euler::RANS::GetVisFlux_KOWilcox ( TU &&  UMeanXy,
TDiffU &&  DiffUxyPrim,
TN &&  uNorm,
real  mutIn,
real  d,
real  muf,
TVFlux &  vFlux 
)

Compute the RANS viscous flux for the Wilcox k-omega transport equations.

Evaluates the diffusion (viscous) flux contributions for the k and omega equations projected onto the face normal direction. The Wilcox model uses constant diffusion coefficients:

  • sigma_k = 0.5
  • sigma_omega = 0.5

Effective diffusivity: mu_laminar + mu_t * sigma_{k,omega}.

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TNEigen-compatible type for the face unit normal vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, primitive variables).
TVFluxEigen-compatible type for the output viscous flux vector.
Parameters
UMeanXyConservative state vector.
DiffUxyPrimGradient tensor of primitive variables, size dim × nVars.
uNormOutward-pointing face unit normal vector (length dim).
mutInTurbulent dynamic viscosity (precomputed).
dWall distance (unused here, kept for interface consistency).
mufLaminar (molecular) dynamic viscosity.
[out]vFluxOutput viscous flux vector; entries at I4+1 and I4+2 are set.

Definition at line 776 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ GetVisFlux_RealizableKe()

template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void DNDS::Euler::RANS::GetVisFlux_RealizableKe ( TU &&  UMeanXy,
TDiffU &&  DiffUxyPrim,
TN &&  uNorm,
real  mut,
real  d,
real  muPhy,
TVFlux &  vFlux 
)

Compute the RANS viscous flux for the realizable k-epsilon transport equations.

Evaluates the diffusion (viscous) flux contributions for the k and epsilon equations projected onto the face normal direction. The effective diffusivities are:

  • k-equation: (mu_laminar + mu_t / sigma_k), sigma_k = 1.0
  • epsilon-equation: (mu_laminar + mu_t / sigma_e), sigma_e = 1.3
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TNEigen-compatible type for the face unit normal vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, primitive variables).
TVFluxEigen-compatible type for the output viscous flux vector.
Parameters
UMeanXyConservative state vector.
DiffUxyPrimGradient tensor of primitive variables, size dim × nVars.
uNormOutward-pointing face unit normal vector (length dim).
mutTurbulent dynamic viscosity (precomputed).
dWall distance (unused here, kept for interface consistency).
muPhyLaminar (molecular) dynamic viscosity.
[out]vFluxOutput viscous flux vector; entries at I4+1 and I4+2 are set.

Definition at line 148 of file RANS_ke.hpp.

◆ GetVisFlux_SST()

template<int dim, class TU , class TN , class TDiffU , class TVFlux >
void DNDS::Euler::RANS::GetVisFlux_SST ( TU &&  UMeanXy,
TDiffU &&  DiffUxyPrim,
TN &&  uNorm,
real  mutIn,
real  d,
real  muf,
TVFlux &  vFlux 
)

Compute the RANS viscous flux for the k-omega SST transport equations.

Evaluates the diffusion (viscous) flux contributions for the k and omega equations projected onto the face normal direction. The effective diffusion coefficients are blended between the inner-layer (set 1) and outer-layer (set 2) values using the SST blending function F1:

  • sigma_k = sigma_k1 * F1 + sigma_k2 * (1 - F1)
  • sigma_omega = sigma_omega1 * F1 + sigma_omega2 * (1 - F1)

Effective diffusivity for each equation: mu_laminar + mu_t * sigma_{k,omega}.

Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
TNEigen-compatible type for the face unit normal vector.
TDiffUEigen-compatible type for the gradient tensor (dim × nVars, primitive variables).
TVFluxEigen-compatible type for the output viscous flux vector.
Parameters
UMeanXyConservative state vector.
DiffUxyPrimGradient tensor of primitive variables, size dim × nVars.
uNormOutward-pointing face unit normal vector (length dim).
mutInTurbulent dynamic viscosity (precomputed, used for diffusion coefficients).
dWall distance.
mufLaminar (molecular) dynamic viscosity.
[out]vFluxOutput viscous flux vector; entries at I4+1 and I4+2 are set.

Definition at line 511 of file RANS_ke.hpp.

Here is the call graph for this function:

◆ SolveZeroGradEquilibrium()

template<int dim, class TU >
std::tuple< real, real > DNDS::Euler::RANS::SolveZeroGradEquilibrium ( TU &  u,
real  muPhy 
)

Attempt to find equilibrium epsilon for zero-gradient (freestream) conditions via Newton iteration.

Iteratively solves for the epsilon value that zeroes the k-equation source term when all spatial gradients are zero. Uses a simple Newton-Raphson loop with finite-difference Jacobian evaluation.

Warning
This function is noted as not applicable in practice — the zero-gradient equilibrium assumption does not hold for the realizable k-epsilon model. Retained for reference and experimentation only.
Template Parameters
dimSpatial dimension (2 or 3).
TUEigen-compatible type for the conservative state vector.
Parameters
[in,out]uConservative state vector; u(I4+2) (rho*epsilon) is updated in place.
muPhyLaminar (molecular) dynamic viscosity.
Returns
A tuple (rhs_initial, rhs_final) giving the k-equation residual before and after the Newton iteration, for convergence diagnostics.

Definition at line 275 of file RANS_ke.hpp.