8.1
general documentation
Module for turbulence constants
+ Collaboration diagram for Module for turbulence constants:

Variables

double precision, save xkappa = 0.42d0
 $ \kappa $ Karman constant. (= 0.42) Useful if and only if iturb >= 10. (mixing length, $k-\varepsilon$, $R_{ij}-\varepsilon$, LES, v2f or $k-\omega$) More...
 
real(c_double), pointer, save ypluli
 limit value of $y^+$ for the viscous sublayer. ypluli depends on the chosen wall function: it is initialized to 10.88 for the scalable wall function (iwallf=4), otherwise it is initialized to $1/\kappa\approx 2,38$. In LES, ypluli is taken by default to be 10.88. More...
 
real(c_double), pointer, save cmu
 constant $C_\mu$ for all the RANS turbulence models Warning, different values for the v2f model Useful if and only if iturb = 20, 21, 30, 31, 50, 51 or 60 ( $k-\varepsilon$, $R_{ij}-\varepsilon$ or $k-\omega$) More...
 
real(c_double), pointer, save csrij
 Coefficient of interfacial coefficient in k-eps, used in Lagrange treatment. More...
 
double precision, save xcl = 0.122d0
 constant of the Rij-epsilon EBRSM More...
 
double precision, save xa1 = 0.1d0
 constant in the expression of Ce1' for the Rij-epsilon EBRSM More...
 
double precision, save xct = 6.d0
 constant of the Rij-epsilon EBRSM More...
 
real(c_double), pointer, save almax
 is a characteristic macroscopic length of the domain, used for the initialization of the turbulence and the potential clipping (with iclkep=1) More...
 
real(c_double), pointer, save uref
 the characteristic flow velocity, used for the initialization of the turbulence. Negative value: not initialized. More...
 
real(c_double), pointer, save xlesfl
 constant used in the definition of LES filtering diameter: $ \delta = \text{xlesfl} . (\text{ales} . volume)^{\text{bles}}$ xlesfl is a constant used to define, for each cell $\omega_i$, the width of the (implicit) filter: $\overline{\Delta}=xlesfl(ales*|\Omega_i|)^{bles}$
Useful if and only if iturb = 40 or 41 More...
 
real(c_double), pointer, save ales
 constant used to define, for each cell $Omega_i$, the width of the (implicit) filter: More...
 
real(c_double), pointer, save bles
 constant used to define, for each cell $Omega_i$, More...
 
real(c_double), pointer, save csmago
 Smagorinsky constant used in the Smagorinsky model for LES. The sub-grid scale viscosity is calculated by $\displaystyle\mu_{sg}= \rho C_{smago}^2\bar{\Delta}^2\sqrt{2\bar{S}_{ij}\bar{S}_{ij}}$ where $\bar{\Delta}$ is the width of the filter and $\bar{S}_{ij}$ the filtered strain rate. More...
 
real(c_double), pointer, save xlesfd
 ratio between explicit and explicit filter width for a dynamic model constant used to define, for each cell $\Omega_i$, the width of the explicit filter used in the framework of the LES dynamic model: $\widetilde{\overline{\Delta}}=xlesfd\overline{\Delta}$. More...
 
real(c_double), pointer, save cdries
 van Driest constant appearing in the van Driest damping function applied to the Smagorinsky constant: More...
 
double precision, save volmin
 minimal control volume More...
 
double precision, save volmax
 maximal control volume More...
 
double precision, save voltot
 total domain volume More...
 
real(c_double), pointer, save xclt
 constant of EB-AFM and EB-DFM More...
 

Detailed Description

Variable Documentation

◆ ales

real(c_double), pointer, save ales

constant used to define, for each cell $Omega_i$, the width of the (implicit) filter:

  • $\overline{\Delta}=xlesfl(ales*|Omega_i|)^{bles}$

Useful if and only if iturb = 40 or 41.

◆ almax

real(c_double), pointer, save almax

is a characteristic macroscopic length of the domain, used for the initialization of the turbulence and the potential clipping (with iclkep=1)

  • Negative value: not initialized (the code then uses the cubic root of the domain volume).

Useful if and only if iturb = 20, 21, 30, 31, 50 or 60 (RANS models)

◆ bles

real(c_double), pointer, save bles

constant used to define, for each cell $Omega_i$,

the width of the (implicit) filter:

  • $\overline{\Delta}=xlesfl(ales*|Omega_i|)^{bles}$

Useful if and only if iturb = 40 or 41

◆ cdries

real(c_double), pointer, save cdries

van Driest constant appearing in the van Driest damping function applied to the Smagorinsky constant:

  • $ (1-\exp^{(-y^+/cdries}) $.

Useful if and only if iturb = 40 or 41

◆ cmu

real(c_double), pointer, save cmu

constant $C_\mu$ for all the RANS turbulence models Warning, different values for the v2f model Useful if and only if iturb = 20, 21, 30, 31, 50, 51 or 60 ( $k-\varepsilon$, $R_{ij}-\varepsilon$ or $k-\omega$)

◆ csmago

real(c_double), pointer, save csmago

Smagorinsky constant used in the Smagorinsky model for LES. The sub-grid scale viscosity is calculated by $\displaystyle\mu_{sg}= \rho C_{smago}^2\bar{\Delta}^2\sqrt{2\bar{S}_{ij}\bar{S}_{ij}}$ where $\bar{\Delta}$ is the width of the filter and $\bar{S}_{ij}$ the filtered strain rate.

Useful if and only if iturb = 40

Note
In theory Smagorinsky constant is 0.18. For a planar canal plan, 0.065 value is rather taken.

◆ csrij

real(c_double), pointer, save csrij

Coefficient of interfacial coefficient in k-eps, used in Lagrange treatment.

constant $C_s$ for the $R_{ij}-\varepsilon$ models.

◆ uref

real(c_double), pointer, save uref

the characteristic flow velocity, used for the initialization of the turbulence. Negative value: not initialized.

Useful if and only if iturb = 20, 21, 30, 31, 50 or 60 (RANS model) and the turbulence is not initialized somewhere else (restart file or subroutine cs_user_initialization)

◆ volmax

double precision, save volmax

maximal control volume

◆ volmin

double precision, save volmin

minimal control volume

◆ voltot

double precision, save voltot

total domain volume

◆ xa1

double precision, save xa1 = 0.1d0

constant in the expression of Ce1' for the Rij-epsilon EBRSM

◆ xcl

double precision, save xcl = 0.122d0

constant of the Rij-epsilon EBRSM

◆ xclt

real(c_double), pointer, save xclt

constant of EB-AFM and EB-DFM

◆ xct

double precision, save xct = 6.d0

constant of the Rij-epsilon EBRSM

◆ xkappa

double precision, save xkappa = 0.42d0

$ \kappa $ Karman constant. (= 0.42) Useful if and only if iturb >= 10. (mixing length, $k-\varepsilon$, $R_{ij}-\varepsilon$, LES, v2f or $k-\omega$)

◆ xlesfd

real(c_double), pointer, save xlesfd

ratio between explicit and explicit filter width for a dynamic model constant used to define, for each cell $\Omega_i$, the width of the explicit filter used in the framework of the LES dynamic model: $\widetilde{\overline{\Delta}}=xlesfd\overline{\Delta}$.

Useful if and only if iturb = 41

◆ xlesfl

real(c_double), pointer, save xlesfl

constant used in the definition of LES filtering diameter: $ \delta = \text{xlesfl} . (\text{ales} . volume)^{\text{bles}}$ xlesfl is a constant used to define, for each cell $\omega_i$, the width of the (implicit) filter: $\overline{\Delta}=xlesfl(ales*|\Omega_i|)^{bles}$
Useful if and only if iturb = 40 or 41

◆ ypluli

real(c_double), pointer, save ypluli

limit value of $y^+$ for the viscous sublayer. ypluli depends on the chosen wall function: it is initialized to 10.88 for the scalable wall function (iwallf=4), otherwise it is initialized to $1/\kappa\approx 2,38$. In LES, ypluli is taken by default to be 10.88.

Always useful