1 #ifndef __CS_WALL_FUNCTIONS_H__
2 #define __CS_WALL_FUNCTIONS_H__
139 const double ydvisc = y / l_visc;
145 *
yplus = *ustar * ydvisc;
150 *ustar = sqrt(
vel / ydvisc);
151 *
yplus = *ustar * ydvisc;
210 double ustarwer, ustarmin, ustaro, ydvisc;
219 reynolds =
vel * ydvisc;
228 *ustar = sqrt(
vel / ydvisc);
229 *
yplus = *ustar * ydvisc;
245 ustaro =
CS_MAX(ustarwer, ustarmin);
250 for (iter = 0; iter < niter_max
251 && fabs(*ustar - ustaro) >=
eps * ustaro; iter++) {
257 if (iter >= niter_max) {
258 bft_printf(
_(
"WARNING: non-convergence in the computation\n"
259 "******** of the friction velocity\n\n"
260 "friction vel: %f \n" ), *ustar);
264 *
yplus = *ustar * ydvisc;
295 f_blend = exp(-0.25*cuv*pow(yp,3));
296 uplus = f_blend*yp + (log(yp)/
ka +B)*(1.-exp(-pow(yp/y0,n)))*(1-f_blend);
319 dupdyp = exp(-0.25*cuv*pow(yp,3))
320 - 0.75*cuv*pow(yp,3.)*exp(-0.25*cuv*pow(yp,3.))
321 + n*(1.-exp(-0.25*cuv*pow(yp,3.)))*(pow(yp,n-1.)/pow(y0,n))
322 *exp(-pow(yp/y0,n))*((1./
ka)*log(yp)+B)
323 + 0.75*cuv*pow(yp,2.)*exp(-0.25*cuv*pow(yp,3.))
324 *(1.-exp(-pow(yp/y0,n)))*((1./
ka)*log(yp)+B)
325 + (1./
ka/yp)*(1.-exp(-pow(yp/y0,n)))*(1-exp(-0.25*cuv*pow(yp,3.)));
372 double Re, g, t_visc_durb;
382 Re = sqrt(kinetic_en) * y / l_visc;
387 + g * l_visc *
vel/y);
390 *
yplus = *uk * y / l_visc;
393 *ustar =
vel / uplus;
395 if (*
yplus < 1.e-1) {
408 cstcuv, csty0, cstN);
413 *ypup = *
yplus / uplus;
427 assert(turb_model != NULL);
428 if (turb_model->
itytur == 3)
429 t_visc_durb = t_visc / (kinetic_en *
cs_turb_cmu ) * rnnb * 0.22;
431 t_visc_durb = t_visc;
434 = 1. - *ypup * (2. * sqrt(l_visc / t_visc_durb * dup1 * (1. - dup1))
487 double rcprod, ml_visc, Re, g;
492 Re = sqrt(kinetic_en) * y / l_visc;
496 + g * l_visc *
vel / y);
498 *
yplus = *uk * y / l_visc;
515 if (*
yplus > 1.e-12) {
573 double rcprod, ml_visc, Re, g;
577 Re = sqrt(kinetic_en) * y / l_visc;
581 + g * l_visc *
vel / y);
583 *
yplus = *uk * y / l_visc;
587 *
yplus = *uk * y / l_visc;
632 static double aa[11] = {-0.0091921, 3.9577, 0.031578,
633 -0.51013, -2.3254, -0.72665,
634 2.969, 0.48506, -1.5944,
637 cs_real_t y1,y2,y3,y4,y5,y6,y7,y8,y9,y10, uplus;
639 y1 = 0.25 * log(
yplus);
720 double urplus, d_up, lmk15;
733 cs_real_t dyrp = 0.9 * (sqrt(krp) - krp * exp(-krp / 6.));
738 else if (dyrp <= 200.)
743 if (yrplus <= 1.e-1) {
758 }
else if (yrplus <= 200.) {
765 *ypup = *
yplus / (urplus-d_up);
774 *cofimp = 1. - (2. / (1. + *lmk) - 1. / (1. + lmk15)) * *ypup;
784 *ypup = *
yplus / (urplus-d_up);
793 *cofimp = 1. - (2. / *lmk - 1. / lmk15) * *ypup;
798 *ustar =
vel / (urplus-d_up);
857 double rcprod, ml_visc, Re, g;
880 Re = sqrt(kinetic_en) * (y + y0) / l_visc;
884 + g * l_visc *
vel / (y + y0));
887 *
yplus = *uk * y / l_visc;
890 *dplus = *uk * y0 / l_visc;
913 *ustar =
vel / uplus;
915 bft_printf(
"uet=%f, u=%f, uplus=%f, yk=%f, duplus=%f\n", *ustar,
vel, uplus, *
yplus, 1./uplus);
917 *ypup = *
yplus / uplus;
923 * ( 2. * rcprod - 1. / (2. * *
yplus + *dplus));
972 *ustar = sqrt(
vel * l_visc / y);
973 *
yplus = *ustar * y / l_visc;
1071 if (
yplus > (*yplim)) {
1073 (*htur) = prl *
yplus / tplus;
1082 (*yplim) = pow(1000./prl, 1./3.);
1084 a2 = 15.*pow(prl, 2./3.);
1088 (*htur) = prl *
yplus / tplus;
1094 (*htur) = prl *
yplus / tplus;
1135 const int ninter_max = 100;
1147 int npeff =
CS_MAX((
int)(ypint / ypmax * (
double)(ninter_max)), 1);
1149 double dy = ypint / (double)(npeff);
1154 for (
int ip = 1; ip <= npeff; ip++) {
1155 double yp = ypint * (double)(ip) / (double)(npeff);
1157 stplus += dy / (1. + prlrat * 0.5 * (nut1 + nut2));
1161 if (
yplus > ypint) {
1163 stplus += log( (1. + r*
yplus) / (1. + r*ypint)) / r;
1166 if (stplus >= 1.e-6)
1167 *htur =
yplus / stplus;
1223 (*htur) = prl *
yplus / tplus;
int bft_printf(const char *const format,...)
Replacement for printf() with modifiable behavior.
Definition: bft_printf.c:140
#define BEGIN_C_DECLS
Definition: cs_defs.h:509
double cs_real_t
Floating-point value.
Definition: cs_defs.h:319
#define CS_MIN(a, b)
Definition: cs_defs.h:472
#define _(String)
Definition: cs_defs.h:63
#define CS_MAX(a, b)
Definition: cs_defs.h:473
#define CS_PROCF(x, y)
Definition: cs_defs.h:523
unsigned long cs_gnum_t
global mesh entity number
Definition: cs_defs.h:298
#define CS_UNUSED(x)
Definition: cs_defs.h:495
#define END_C_DECLS
Definition: cs_defs.h:510
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:313
@ eps
Definition: cs_field_pointer.h:71
@ vel
Definition: cs_field_pointer.h:68
@ yplus
Definition: cs_field_pointer.h:236
const cs_real_t cs_math_epzero
double cs_turb_vdriest
Definition: cs_turbulence_model.c:446
double cs_turb_crij2
Definition: cs_turbulence_model.c:536
double cs_turb_crij1
Definition: cs_turbulence_model.c:530
double cs_turb_cmu025
Definition: cs_turbulence_model.c:499
double cs_turb_cstlog
Definition: cs_turbulence_model.c:457
cs_turb_model_t * cs_get_glob_turb_model(void)
Provide write access to turbulence model structure.
Definition: cs_turbulence_model.c:1494
double cs_turb_dpow
Definition: cs_turbulence_model.c:488
double cs_turb_apow
Definition: cs_turbulence_model.c:482
double cs_turb_cmu
Definition: cs_turbulence_model.c:496
double cs_turb_bpow
Definition: cs_turbulence_model.c:485
double cs_turb_xkappa
Definition: cs_turbulence_model.c:437
cs_wall_f_s_type_t
Definition: cs_wall_functions.h:68
@ CS_WALL_F_S_SMOOTH_ROUGH
Definition: cs_wall_functions.h:75
@ CS_WALL_F_S_MONIN_OBUKHOV
Definition: cs_wall_functions.h:74
@ CS_WALL_F_S_VDRIEST
Definition: cs_wall_functions.h:72
@ CS_WALL_F_S_ARPACI_LARSEN
Definition: cs_wall_functions.h:71
@ CS_WALL_F_S_UNSET
Definition: cs_wall_functions.h:70
@ CS_WALL_F_S_LOUIS
Definition: cs_wall_functions.h:73
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)
Definition: cs_wall_functions.h:310
void cs_wall_functions_scalar(cs_wall_f_s_type_t iwalfs, cs_real_t l_visc, cs_real_t prl, cs_real_t prt, cs_real_t rough_t, cs_real_t uk, 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 flo...
Definition: cs_wall_functions.c:537
void cs_immersed_boundary_wall_functions(int f_id, cs_real_t *st_exp, cs_real_t *st_imp)
Compute boundary contributions for all immersed boundaries.
Definition: cs_wall_functions.c:606
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_gnum_t *nsubla, cs_gnum_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 fric...
Definition: cs_wall_functions.h:356
static void cs_wall_functions_1scale_log(cs_real_t l_visc, cs_real_t vel, cs_real_t y, int *iuntur, cs_gnum_t *nsubla, cs_gnum_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.
Definition: cs_wall_functions.h:196
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 rough_d, cs_real_t kinetic_en, int *iuntur, cs_gnum_t *nsubla, cs_gnum_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.
Definition: cs_wall_functions.h:837
static void cs_wall_functions_s_smooth_rough(cs_real_t l_visc, cs_real_t prl, cs_real_t prt, cs_real_t rough_t, cs_real_t uk, cs_real_t yplus, cs_real_t dplus, cs_real_t *htur)
Rough Smooth Thermal Wall Function - Prototype.
Definition: cs_wall_functions.h:1192
static void cs_wall_functions_1scale_power(cs_real_t l_visc, cs_real_t vel, cs_real_t y, int *iuntur, cs_gnum_t *nsubla, cs_gnum_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.
Definition: cs_wall_functions.h:125
const cs_wall_functions_t * cs_glob_wall_functions
static cs_real_t _vdriest_dupdyp_integral(cs_real_t yplus)
Definition: cs_wall_functions.h:629
static void cs_wall_functions_s_arpaci_larsen(cs_real_t l_visc, cs_real_t prl, cs_real_t prt, cs_real_t rough_t, cs_real_t uk, 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 v...
Definition: cs_wall_functions.h:1025
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_gnum_t *nsubla, cs_gnum_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.
Definition: cs_wall_functions.h:950
void cs_wall_functions_velocity(cs_wall_f_type_t iwallf, cs_real_t l_visc, cs_real_t t_visc, cs_real_t vel, cs_real_t y, cs_real_t rough_d, cs_real_t rnnb, cs_real_t kinetic_en, int *iuntur, cs_gnum_t *nsubla, cs_gnum_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 / .
Definition: cs_wall_functions.c:350
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_gnum_t *nsubla, cs_gnum_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.
Definition: cs_wall_functions.h:703
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 is computed thanks to a numerical integration of:
Definition: cs_wall_functions.h:1127
void hturbp(const int *const iwalfs, const cs_real_t *const l_visc, const cs_real_t *const prl, const cs_real_t *const prt, const cs_real_t *const rough_t, const cs_real_t *const uk, const cs_real_t *const yplus, const cs_real_t *const dplus, cs_real_t *htur, cs_real_t *yplim)
Definition: cs_wall_functions.c:281
cs_wall_f_type_t
Definition: cs_wall_functions.h:54
@ CS_WALL_F_1SCALE_LOG
Definition: cs_wall_functions.h:59
@ CS_WALL_F_1SCALE_POWER
Definition: cs_wall_functions.h:58
@ CS_WALL_F_2SCALES_SMOOTH_ROUGH
Definition: cs_wall_functions.h:63
@ CS_WALL_F_2SCALES_LOG
Definition: cs_wall_functions.h:60
@ CS_WALL_F_2SCALES_VDRIEST
Definition: cs_wall_functions.h:62
@ CS_WALL_F_DISABLED
Definition: cs_wall_functions.h:57
@ CS_WALL_F_UNSET
Definition: cs_wall_functions.h:56
@ CS_WALL_F_2SCALES_CONTINUOUS
Definition: cs_wall_functions.h:64
@ CS_WALL_F_SCALABLE_2SCALES_LOG
Definition: cs_wall_functions.h:61
void wallfunctions(const int *const iwallf, const cs_lnum_t *const ifac, const cs_real_t *const viscosity, const cs_real_t *const t_visc, const cs_real_t *const vel, const cs_real_t *const y, const cs_real_t *const rough_d, const cs_real_t *const rnnb, const cs_real_t *const 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)
Definition: cs_wall_functions.c:216
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_gnum_t *nsubla, cs_gnum_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 .
Definition: cs_wall_functions.h:554
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_gnum_t *nsubla, cs_gnum_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 k...
Definition: cs_wall_functions.h:471
cs_wall_functions_t * cs_get_glob_wall_functions(void)
Definition: cs_wall_functions.c:317
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)
Definition: cs_wall_functions.h:286
double precision epzero
epsilon
Definition: cstnum.f90:40
real(c_double), pointer, save ypluli
limit value of for the viscous sublayer. ypluli depends on the chosen wall function: it is initializ...
Definition: cstphy.f90:303
integer(c_int), pointer, save iwallf
Wall functions Indicates the type of wall function used for the velocity boundary conditions on a fri...
Definition: optcal.f90:434
integer(c_int), pointer, save iwalfs
Wall functions for scalar.
Definition: optcal.f90:439
double precision, dimension(ncharm), save a2
Definition: cpincl.f90:233
double precision, dimension(4, npot), save ka
Definition: ppcpfu.f90:159
Turbulence model general options descriptor.
Definition: cs_turbulence_model.h:115
int itytur
Definition: cs_turbulence_model.h:139
wall functions descriptor.
Definition: cs_wall_functions.h:82
double ypluli
Definition: cs_wall_functions.h:88
cs_wall_f_s_type_t iwalfs
Definition: cs_wall_functions.h:86
cs_wall_f_type_t iwallf
Definition: cs_wall_functions.h:84