My Project
programmer's documentation
Data Structures | Enumerations | Functions | Variables
cs_wall_functions.h File Reference
#include <bft_printf.h>
#include "cs_base.h"
#include "cs_turbulence_model.h"
Include dependency graph for cs_wall_functions.h:

Go to the source code of this file.

Data Structures

struct  cs_wall_functions_t
 wall functions descriptor. More...
 

Enumerations

enum  cs_wall_f_type_t {
  CS_WALL_F_DISABLED, CS_WALL_F_1SCALE_POWER, CS_WALL_F_1SCALE_LOG, CS_WALL_F_2SCALES_LOG,
  CS_WALL_F_SCALABLE_2SCALES_LOG, CS_WALL_F_2SCALES_VDRIEST, CS_WALL_F_2SCALES_SMOOTH_ROUGH, CS_WALL_F_2SCALES_CONTINUOUS
}
 
enum  cs_wall_f_s_type_t { CS_WALL_F_S_ARPACI_LARSEN, CS_WALL_F_S_VDRIEST }
 

Functions

static void cs_wall_functions_1scale_power (cs_real_t l_visc, cs_real_t vel, cs_real_t y, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp)
 Power law: Werner & Wengle. More...
 
static void cs_wall_functions_1scale_log (cs_lnum_t ifac, cs_real_t l_visc, cs_real_t vel, cs_real_t y, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp)
 Log law: piecewise linear and log, with one velocity scale based on the friction. More...
 
static cs_real_t _uplus (cs_real_t yp, cs_real_t ka, cs_real_t B, cs_real_t cuv, cs_real_t y0, cs_real_t n)
 
static cs_real_t _dupdyp (cs_real_t yp, cs_real_t ka, cs_real_t B, cs_real_t cuv, cs_real_t y0, cs_real_t n)
 
static void cs_wall_functions_2scales_continuous (cs_real_t rnnb, cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp)
 Continuous law of the wall between the linear and log law, with two velocity scales based on the friction and the turbulent kinetic energy. Can be used with RANS model either in high Reynolds or in low Reynolds (if the underlying RANS model is compatible). More...
 
static void cs_wall_functions_2scales_log (cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp)
 Log law: piecewise linear and log, with two velocity scales based on the friction and the turbulent kinetic energy. More...
 
static void cs_wall_functions_2scales_scalable (cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *dplus, cs_real_t *ypup, cs_real_t *cofimp)
 Scalable wall function: shift the wall if $ y^+ < y^+_{lim} $. More...
 
static cs_real_t _vdriest_dupdyp_integral (cs_real_t yplus)
 
static void cs_wall_functions_2scales_vdriest (cs_real_t rnnb, cs_real_t l_visc, cs_real_t vel, cs_real_t y, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp, cs_real_t *lmk, cs_real_t kr, bool wf)
 Two velocity scales wall function using Van Driest mixing length. More...
 
static void cs_wall_functions_2scales_smooth_rough (cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t roughness, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *dplus, cs_real_t *ypup, cs_real_t *cofimp)
 Two velocity scales wall function with automatic switch from rough to smooth. More...
 
static void cs_wall_functions_disabled (cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *dplus, cs_real_t *ypup, cs_real_t *cofimp)
 No wall function. More...
 
static void cs_wall_functions_s_arpaci_larsen (cs_real_t prl, cs_real_t prt, cs_real_t yplus, cs_real_t dplus, cs_real_t *htur, cs_real_t *yplim)
 The correction of the exchange coefficient is computed thanks to a similarity model between dynamic viscous sub-layer and themal sub-layer. More...
 
static void cs_wall_functions_s_vdriest (cs_real_t prl, cs_real_t prt, cs_real_t yplus, cs_real_t *htur)
 The correction of the exchange coefficient $ h_{tur} = \sigma \dfrac{y^+}{t^+} $ is computed thanks to a numerical integration of: More...
 
void CS_PROCF (wallfunctions, WALLFUNCTIONS)(const cs_int_t *const iwallf
 
void CS_PROCF (hturbp, HTURBP)(const cs_int_t *const iwalfs
 
cs_wall_functions_tcs_get_glob_wall_functions (void)
 
void cs_wall_functions_velocity (cs_wall_f_type_t iwallf, cs_lnum_t ifac, cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t roughness, cs_real_t rnnb, cs_real_t kinetic_en, int *iuntur, cs_lnum_t *nsubla, cs_lnum_t *nlogla, cs_real_t *ustar, cs_real_t *uk, cs_real_t *yplus, cs_real_t *ypup, cs_real_t *cofimp, cs_real_t *dplus)
 Compute the friction velocity and $y^+$ / $u^+$. More...
 
void cs_wall_functions_scalar (cs_wall_f_s_type_t iwalfs, cs_real_t prl, cs_real_t prt, cs_real_t yplus, cs_real_t dplus, cs_real_t *htur, cs_real_t *yplim)
 Compute the correction of the exchange coefficient between the fluid and the wall for a turbulent flow. More...
 

Variables

const cs_wall_functions_tcs_glob_wall_functions
 
void const cs_lnum_t *const ifac
 
void const cs_lnum_t *const const cs_real_t *const viscosity
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const t_visc
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const vel
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const y
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const roughness
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const rnnb
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const kinetic_en
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_tiuntur
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_tnsubla
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_tnlogla
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_tustar
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_tuk
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_typlus
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_t cs_real_typup
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_t cs_real_t cs_real_tcofimp
 
void const cs_lnum_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_t cs_real_t cs_real_t cs_real_tdplus
 
void const cs_real_t *const prl
 
void const cs_real_t *const const cs_real_t *const prt
 
void const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_real_thtur
 
void const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const cs_real_t cs_real_typlim
 

Enumeration Type Documentation

◆ cs_wall_f_s_type_t

Enumerator
CS_WALL_F_S_ARPACI_LARSEN 
CS_WALL_F_S_VDRIEST 

◆ cs_wall_f_type_t

Enumerator
CS_WALL_F_DISABLED 
CS_WALL_F_1SCALE_POWER 
CS_WALL_F_1SCALE_LOG 
CS_WALL_F_2SCALES_LOG 
CS_WALL_F_SCALABLE_2SCALES_LOG 
CS_WALL_F_2SCALES_VDRIEST 
CS_WALL_F_2SCALES_SMOOTH_ROUGH 
CS_WALL_F_2SCALES_CONTINUOUS 

Function Documentation

◆ _dupdyp()

static cs_real_t _dupdyp ( cs_real_t  yp,
cs_real_t  ka,
cs_real_t  B,
cs_real_t  cuv,
cs_real_t  y0,
cs_real_t  n 
)
inlinestatic

◆ _uplus()

static cs_real_t _uplus ( cs_real_t  yp,
cs_real_t  ka,
cs_real_t  B,
cs_real_t  cuv,
cs_real_t  y0,
cs_real_t  n 
)
inlinestatic

◆ _vdriest_dupdyp_integral()

static cs_real_t _vdriest_dupdyp_integral ( cs_real_t  yplus)
inlinestatic

◆ cs_get_glob_wall_functions()

cs_wall_functions_t* cs_get_glob_wall_functions ( void  )

◆ CS_PROCF() [1/2]

void CS_PROCF ( hturbp  ,
HTURBP   
) const

◆ CS_PROCF() [2/2]

void CS_PROCF ( wallfunctions  ,
WALLFUNCTIONS   
) const

◆ cs_wall_functions_1scale_log()

static void cs_wall_functions_1scale_log ( cs_lnum_t  ifac,
cs_real_t  l_visc,
cs_real_t  vel,
cs_real_t  y,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Log law: piecewise linear and log, with one velocity scale based on the friction.

Parameters
[in]ifacface number
[in]l_visckinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_1scale_power()

static void cs_wall_functions_1scale_power ( cs_real_t  l_visc,
cs_real_t  vel,
cs_real_t  y,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Power law: Werner & Wengle.

Parameters
[in]l_visckinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_2scales_continuous()

static void cs_wall_functions_2scales_continuous ( cs_real_t  rnnb,
cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Continuous law of the wall between the linear and log law, with two velocity scales based on the friction and the turbulent kinetic energy. Can be used with RANS model either in high Reynolds or in low Reynolds (if the underlying RANS model is compatible).

Parameters
[in]rnnb$\vec{n}.(\tens{R}\vec{n})$
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]kinetic_enturbulent kinetic energy
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_2scales_log()

static void cs_wall_functions_2scales_log ( cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Log law: piecewise linear and log, with two velocity scales based on the friction and the turbulent kinetic energy.

Parameters
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]kinetic_enturbulent kinetic energy
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_2scales_scalable()

static void cs_wall_functions_2scales_scalable ( cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t dplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Scalable wall function: shift the wall if $ y^+ < y^+_{lim} $.

Parameters
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]kinetic_enturbulent kinetic energy
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]dplusdimensionless shift to the wall for scalable wall functions
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_2scales_smooth_rough()

static void cs_wall_functions_2scales_smooth_rough ( cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  roughness,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t dplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

Two velocity scales wall function with automatic switch from rough to smooth.

$ u^+ $ is computed as follows:

\[ u^+ = \dfrac{1}{\kappa} \ln \left(\dfrac{(y+\xi) u_k}{\nu + \alpha \xi u_k} \right) + Cst_{smooth} \]

with $ \alpha = \exp \left(- \kappa(Cst_{rough}-Cst_{smooth})\right) $.

Parameters
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]roughnessroughness
[in]kinetic_enturbulent kinetic energy
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]dplusdimensionless shift to the wall for scalable wall functions
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_2scales_vdriest()

static void cs_wall_functions_2scales_vdriest ( cs_real_t  rnnb,
cs_real_t  l_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp,
cs_real_t lmk,
cs_real_t  kr,
bool  wf 
)
inlinestatic

Two velocity scales wall function using Van Driest mixing length.

$ u^+ $ is computed as follows:

\[ u^+ = \int_0^{y_k^+} \dfrac{dy_k^+}{1+L_m^k} \]

with $ L_m^k $ standing for Van Driest mixing length:

\[ L_m^k = \kappa y_k^+ (1 - exp(\dfrac{-y_k^+}{A})) \]

.

A polynome fitting the integral is used for $ y_k^+ < 200 $, and a log law is used for $ y_k^+ >= 200 $.

A wall roughness can be taken into account in the mixing length as proposed by Rotta (1962) with Cebeci & Chang (1978) correlation.

Parameters
[in]rnnb$\vec{n}.(\tens{R}\vec{n})$
[in]l_visckinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]kinetic_enturbulent kinetic energy
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production
[in]lmkdimensionless mixing length
[in]krwall roughness
[in]wfenable full wall function computation, if false, uk is not recomputed and uplus is the only output

◆ cs_wall_functions_disabled()

static void cs_wall_functions_disabled ( cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t dplus,
cs_real_t ypup,
cs_real_t cofimp 
)
inlinestatic

No wall function.

Parameters
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[out]iunturindicator: 0 in the viscous sublayer
[out]nsublacounter of cell in the viscous sublayer
[out]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]dplusdimensionless shift to the wall for scalable wall functions
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production

◆ cs_wall_functions_s_arpaci_larsen()

static void cs_wall_functions_s_arpaci_larsen ( cs_real_t  prl,
cs_real_t  prt,
cs_real_t  yplus,
cs_real_t  dplus,
cs_real_t htur,
cs_real_t yplim 
)
inlinestatic

The correction of the exchange coefficient is computed thanks to a similarity model between dynamic viscous sub-layer and themal sub-layer.

$ T^+ $ is computed as follows:

  • For a laminar Prandtl number smaller than 0.1 (such as liquid metals), the standard model with two sub-layers (Prandtl-Taylor) is used.
  • For a laminar Prandtl number larger than 0.1 (such as liquids and gaz), a model with three sub-layers (Arpaci-Larsen) is used.

The final exchange coefficient is:

\[ h = \dfrac{K}{\centip \centf} h_{tur} \]

Parameters
[in]prllaminar Prandtl number
[in]prtturbulent Prandtl number
[in]yplusdimensionless distance to the wall
[out]dplusdimensionless shift to the wall for scalable wall functions
[out]hturcorrected exchange coefficient
[out]yplimvalue of the limit for $ y^+ $

◆ cs_wall_functions_s_vdriest()

static void cs_wall_functions_s_vdriest ( cs_real_t  prl,
cs_real_t  prt,
cs_real_t  yplus,
cs_real_t htur 
)
inlinestatic

The correction of the exchange coefficient $ h_{tur} = \sigma \dfrac{y^+}{t^+} $ is computed thanks to a numerical integration of:

\[ \dfrac{T^+}{\sigma} = \int_0^{y_k^+} \dfrac{dy_k^+}{1+\dfrac{\sigma}{\sigma_t}\nu_t^+} \]

with $ \nu_t^+ = L_m^k $ as assumed in the derivation of the two scales wall function using Van Driest mixing length. Therefore $ L_m^k = \kappa y_k^+(1 - exp(\dfrac{-y_k^+}{A})) $ is taken.

Notice that we integrate up to $ y^+=100 $ (YP100), beyond that value the profile is prolonged by a logarithm relying on the fact that $ \nu_t^+=\kappa y^+ $ beyond $ y^+=100 $.

Parameters
[in]prlmolecular Prandtl number ( $ Pr=\sigma=\frac{\mu C_p}{\lambda_f} $ )
[in]prtturbulent Prandtl number ( $ \sigma_t $ )
[in]yplusdimensionless distance to the wall
[out]hturcorrected exchange coefficient

◆ cs_wall_functions_scalar()

void cs_wall_functions_scalar ( cs_wall_f_s_type_t  iwalfs,
cs_real_t  prl,
cs_real_t  prt,
cs_real_t  yplus,
cs_real_t  dplus,
cs_real_t htur,
cs_real_t yplim 
)

Compute the correction of the exchange coefficient between the fluid and the wall for a turbulent flow.

This is function of the dimensionless distance to the wall $ y^+ = \dfrac{\centip \centf u_\star}{\nu}$.

Then the return coefficient reads:

\[ h_{tur} = Pr \dfrac{y^+}{T^+} \]

Parameters
[in]iwalfstype of wall functions for scalar
[in]prllaminar Prandtl number
[in]prtturbulent Prandtl number
[in]yplusdimensionless distance to the wall
[in]dplusdimensionless distance for scalable wall functions
[out]hturcorrected exchange coefficient
[out]yplimvalue of the limit for $ y^+ $

◆ cs_wall_functions_velocity()

void cs_wall_functions_velocity ( cs_wall_f_type_t  iwallf,
cs_lnum_t  ifac,
cs_real_t  l_visc,
cs_real_t  t_visc,
cs_real_t  vel,
cs_real_t  y,
cs_real_t  roughness,
cs_real_t  rnnb,
cs_real_t  kinetic_en,
int *  iuntur,
cs_lnum_t nsubla,
cs_lnum_t nlogla,
cs_real_t ustar,
cs_real_t uk,
cs_real_t yplus,
cs_real_t ypup,
cs_real_t cofimp,
cs_real_t dplus 
)

Compute the friction velocity and $y^+$ / $u^+$.

Parameters
[in]iwallfwall function type
[in]ifacface number
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]roughnessroughness
[in]rnnb$\vec{n}.(\tens{R}\vec{n})$
[in]kinetic_enturbulent kinetic energy
[in]iunturindicator: 0 in the viscous sublayer
[in]nsublacounter of cell in the viscous sublayer
[in]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusnon-dimension wall distance
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production
[out]dplusdimensionless shift to the wall for scalable wall functions
[in]iwallfwall function type
[in]ifacface number
[in]l_visckinematic viscosity
[in]t_viscturbulent kinematic viscosity
[in]velwall projected cell center velocity
[in]ywall distance
[in]roughnessroughness
[in]rnnb$\vec{n}.(\tens{R}\vec{n})$
[in]kinetic_enturbulente kinetic energy
[in]iunturindicator: 0 in the viscous sublayer
[in]nsublacounter of cell in the viscous sublayer
[in]nloglacounter of cell in the log-layer
[out]ustarfriction velocity
[out]ukfriction velocity
[out]yplusdimensionless distance to the wall
[out]ypupyplus projected vel ratio
[out]cofimp$\frac{|U_F|}{|U_I^p|}$ to ensure a good turbulence production
[out]dplusdimensionless shift to the wall for scalable wall functions

Variable Documentation

◆ cofimp

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_t cs_real_t cs_real_t* cofimp

◆ cs_glob_wall_functions

const cs_wall_functions_t* cs_glob_wall_functions

◆ dplus

void const cs_real_t *const const cs_real_t *const const cs_real_t *const const cs_real_t *const dplus

◆ htur

void const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_real_t* htur

◆ ifac

void const cs_lnum_t* const ifac

◆ iuntur

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t* iuntur

◆ kinetic_en

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const kinetic_en

◆ nlogla

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t cs_lnum_t* nlogla

◆ nsubla

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t* nsubla

◆ prl

void const cs_real_t* const prl

◆ prt

void const cs_real_t* const const cs_real_t* const prt

◆ rnnb

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const rnnb

◆ roughness

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const roughness

◆ t_visc

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const t_visc

◆ uk

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t* uk

◆ ustar

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t cs_lnum_t cs_real_t* ustar

◆ vel

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const vel

◆ viscosity

void const cs_lnum_t* const const cs_real_t* const viscosity

◆ y

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const y

◆ yplim

void const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_real_t cs_real_t* yplim

◆ yplus

void const cs_real_t* const const cs_real_t* const const cs_real_t* const yplus

◆ ypup

void const cs_lnum_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const const cs_real_t* const cs_int_t cs_lnum_t cs_lnum_t cs_real_t cs_real_t cs_real_t cs_real_t* ypup