My Project
programmer's documentation
Loading...
Searching...
No Matches
cs_wall_functions.h
Go to the documentation of this file.
1#ifndef __CS_WALL_FUNCTIONS_H__
2#define __CS_WALL_FUNCTIONS_H__
3
4/*============================================================================
5 * Wall functions
6 *============================================================================*/
7
8/*
9 This file is part of Code_Saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2019 EDF S.A.
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU General Public License as published by the Free Software
15 Foundation; either version 2 of the License, or (at your option) any later
16 version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21 details.
22
23 You should have received a copy of the GNU General Public License along with
24 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25 Street, Fifth Floor, Boston, MA 02110-1301, USA.
26*/
27
28/*----------------------------------------------------------------------------*/
29
30/*----------------------------------------------------------------------------
31 * BFT library headers
32 *----------------------------------------------------------------------------*/
33
34#include <bft_printf.h>
35
36/*----------------------------------------------------------------------------
37 * Local headers
38 *----------------------------------------------------------------------------*/
39
40#include "cs_base.h"
41#include "cs_turbulence_model.h"
42
43/*----------------------------------------------------------------------------*/
44
46
47/*=============================================================================
48 * Local Macro definitions
49 *============================================================================*/
50
51/*============================================================================
52 * Type definition
53 *============================================================================*/
54
55/* Wall function type */
56/*--------------------*/
57
70
77
78/* Wall functions descriptor */
79/*---------------------------*/
80
81typedef struct {
82
83 cs_wall_f_type_t iwallf; /* wall function type */
84
85 cs_wall_f_s_type_t iwalfs; /* wall function type for scalars */
86
87 int iwallt; /* exchange coefficient correlation
88 - 0: not used by default
89 - 1: exchange coefficient computed with a
90 correlation */
91
92 double ypluli; /* limit value of y+ for the viscous
93 sublayer */
94
96
97/*============================================================================
98 * Global variables
99 *============================================================================*/
100
101/* Pointer to wall functions descriptor structure */
102
104
105/*============================================================================
106 * Private function definitions
107 *============================================================================*/
108
109/*----------------------------------------------------------------------------*/
126/*----------------------------------------------------------------------------*/
127
128inline static void
131 cs_real_t y,
132 int *iuntur,
133 cs_lnum_t *nsubla,
134 cs_lnum_t *nlogla,
135 cs_real_t *ustar,
136 cs_real_t *uk,
138 cs_real_t *ypup,
139 cs_real_t *cofimp)
140{
141 const double ypluli = cs_glob_wall_functions->ypluli;
142
143 const double ydvisc = y / l_visc;
144
145 /* Compute the friction velocity ustar */
146
147 *ustar = pow((vel/(cs_turb_apow * pow(ydvisc, cs_turb_bpow))), cs_turb_dpow);
148 *uk = *ustar;
149 *yplus = *ustar * ydvisc;
150
151 /* In the viscous sub-layer: U+ = y+ */
152 if (*yplus <= ypluli) {
153
154 *ustar = sqrt(vel / ydvisc);
155 *yplus = *ustar * ydvisc;
156 *uk = *ustar;
157 *ypup = 1.;
158 *cofimp = 0.;
159
160 /* Disable the wall funcion count the cell in the viscous sub-layer */
161 *iuntur = 0;
162 *nsubla += 1;
163
164 /* In the log layer */
165 } else {
166 *ypup = pow(vel, 2. * cs_turb_dpow-1.)
167 / pow(cs_turb_apow, 2. * cs_turb_dpow);
168 *cofimp = 1. + cs_turb_bpow
169 * pow(*ustar, cs_turb_bpow + 1. - 1./cs_turb_dpow)
170 * (pow(2., cs_turb_bpow - 1.) - 2.);
171
172 /* Count the cell in the log layer */
173 *nlogla += 1;
174
175 }
176}
177
178/*----------------------------------------------------------------------------*/
197/*----------------------------------------------------------------------------*/
198
199inline static void
201 cs_real_t l_visc,
203 cs_real_t y,
204 int *iuntur,
205 cs_lnum_t *nsubla,
206 cs_lnum_t *nlogla,
207 cs_real_t *ustar,
208 cs_real_t *uk,
210 cs_real_t *ypup,
211 cs_real_t *cofimp)
212{
213 const double ypluli = cs_glob_wall_functions->ypluli;
214
215 double ustarwer, ustarmin, ustaro, ydvisc;
216 double eps = 0.001;
217 int niter_max = 100;
218 int iter = 0;
219 double reynolds;
220
221 /* Compute the local Reynolds number */
222
223 ydvisc = y / l_visc;
224 reynolds = vel * ydvisc;
225
226 /*
227 * Compute the friction velocity ustar
228 */
229
230 /* In the viscous sub-layer: U+ = y+ */
231 if (reynolds <= ypluli * ypluli) {
232
233 *ustar = sqrt(vel / ydvisc);
234 *yplus = *ustar * ydvisc;
235 *uk = *ustar;
236 *ypup = 1.;
237 *cofimp = 0.;
238
239 /* Disable the wall funcion count the cell in the viscous sub-layer */
240 *iuntur = 0;
241 *nsubla += 1;
242
243 /* In the log layer */
244 } else {
245
246 /* The initial value is Wener or the minimun ustar to ensure convergence */
247 ustarwer = pow(fabs(vel) / cs_turb_apow / pow(ydvisc, cs_turb_bpow),
249 ustarmin = exp(-cs_turb_cstlog * cs_turb_xkappa)/ydvisc;
250 ustaro = CS_MAX(ustarwer, ustarmin);
251 *ustar = (cs_turb_xkappa * vel + ustaro)
252 / (log(ydvisc * ustaro) + cs_turb_xkappa * cs_turb_cstlog + 1.);
253
254 /* Iterative solving */
255 for (iter = 0; iter < niter_max
256 && fabs(*ustar - ustaro) >= eps * ustaro; iter++) {
257 ustaro = *ustar;
258 *ustar = (cs_turb_xkappa * vel + ustaro)
259 / (log(ydvisc * ustaro) + cs_turb_xkappa * cs_turb_cstlog + 1.);
260 }
261
262 if (iter >= niter_max) {
263 bft_printf(_("WARNING: non-convergence in the computation\n"
264 "******** of the friction velocity\n\n"
265 "face number: %d \n"
266 "friction vel: %f \n" ), ifac, *ustar);
267 }
268
269 *uk = *ustar;
270 *yplus = *ustar * ydvisc;
271 *ypup = *yplus / (log(*yplus) / cs_turb_xkappa + cs_turb_cstlog);
272 *cofimp = 1. - *ypup / cs_turb_xkappa * 1.5 / *yplus;
273
274 /* Count the cell in the log layer */
275 *nlogla += 1;
276
277 }
278
279}
280
281/*----------------------------------------------------------------------------
282 * Compute du+/dy+ for a given yk+.
283 * parameters:
284 * yplus <-- dimensionless distance
285 * returns:
286 * the resulting dimensionless velocity.
287 *----------------------------------------------------------------------------*/
288
289inline static cs_real_t
291 cs_real_t ka ,
292 cs_real_t B ,
293 cs_real_t cuv ,
294 cs_real_t y0 ,
295 cs_real_t n )
296{
297 cs_real_t uplus, f_blend ;
298
299 f_blend = exp(-0.25*cuv*pow(yp,3)) ;
300 uplus = f_blend*yp + (log(yp)/ka +B)*(1.-exp(-pow(yp/y0,n)))*(1-f_blend);
301
302 return uplus ;
303}
304
305/*----------------------------------------------------------------------------
306 * Compute du+/dy+ for a given yk+.
307 * parameters:
308 * yplus <-- dimensionless distance
309 * returns:
310 * the resulting dimensionless velocity gradient.
311 *----------------------------------------------------------------------------*/
312
313inline static cs_real_t
315 cs_real_t ka ,
316 cs_real_t B ,
317 cs_real_t cuv ,
318 cs_real_t y0 ,
319 cs_real_t n )
320{
321 cs_real_t dupdyp ;
322
323 dupdyp = exp(-0.25*cuv*pow(yp,3))
324 - 0.75*cuv*pow(yp,3.)*exp(-0.25*cuv*pow(yp,3.))
325 + n*(1.-exp(-0.25*cuv*pow(yp,3.)))*(pow(yp,n-1.)/pow(y0,n))*exp(-pow(yp/y0,n))*((1./ka)*log(yp)+B)
326 + 0.75*cuv*pow(yp,2.)*exp(-0.25*cuv*pow(yp,3.))*(1.-exp(-pow(yp/y0,n)))*((1./ka)*log(yp)+B)
327 + (1./ka/yp)*(1.-exp(-pow(yp/y0,n)))*(1-exp(-0.25*cuv*pow(yp,3.))) ;
328
329 return dupdyp ;
330}
331
332/*----------------------------------------------------------------------------*/
355/*----------------------------------------------------------------------------*/
356
357inline static void
359 cs_real_t l_visc,
360 cs_real_t t_visc,
362 cs_real_t y,
363 cs_real_t kinetic_en,
364 int *iuntur,
365 cs_lnum_t *nsubla,
366 cs_lnum_t *nlogla,
367 cs_real_t *ustar,
368 cs_real_t *uk,
370 cs_real_t *ypup,
371 cs_real_t *cofimp)
372{
373 const double ypluli = cs_glob_wall_functions->ypluli;
374 double Re, g, t_visc_durb;
375 cs_real_t cstcuv, csty0, cstN ;
376 cs_real_t dup1, dup2, uplus ;
377
378 /* Local constants */
379 cstcuv = 1.0674e-3;
380 csty0 = 14.5e0;
381 cstN = 2.25e0;
382
383 /* Iterative process to determine uk through TKE law */
384 Re = sqrt(kinetic_en) * y / l_visc;
385 g = exp(-Re/11.);
386
387 /* Comutation of uk*/
388 *uk = sqrt( (1.-g) * cs_turb_cmu025 * cs_turb_cmu025 * kinetic_en
389 + g * l_visc * vel / y);
390
391 /* Local value of y+, estimated U+ */
392 *yplus = *uk * y / l_visc;
393 uplus = _uplus( *yplus, cs_turb_xkappa, cs_turb_cstlog, cstcuv, csty0, cstN);
394 /* Deduced velocity sclale uet*/
395 *ustar = vel / uplus ;
396
397 if( *yplus < 1.e-1 ) {
398
399 *ypup = 1.0 ;
400 *cofimp = 0.0 ;
401
402 *iuntur = 0;
403 *nsubla += 1;
404
405 } else {
406
407 /* Dimensionless velocity gradient in y+ */
408 dup1 = _dupdyp( *yplus, cs_turb_xkappa, cs_turb_cstlog, cstcuv, csty0, cstN);
409 /* Dimensionless velocity gradient in 2 x y+ */
410 dup2 = _dupdyp(2.0 * *yplus , cs_turb_xkappa, cs_turb_cstlog, cstcuv, csty0, cstN);
411
412 *ypup = *yplus / uplus ;
413
414 /* ------------------------------------------------------------
415 * Cofimp = U,F/U,I is built so that the theoretical expression
416 * of the production P_theo = dup1 * (1.0 - dup1) is equal to
417 * P_calc = mu_t,I * ((U,I - U,F + IF*dup2)/(2IF) )^2
418 * This is a generalization of the process implemented in the 2
419 * scales wall function (iwallf = 3).
420 * ------------------------------------------------------------*/
421
422 /* Turbulent viscocity is modified for RSM so that its expression
423 * remain valid down to the wall, according to Durbin :
424 * nu_t = 0.22 * v'2 * k / eps */
425 if (cs_glob_turb_model->itytur == 3) {
426 t_visc_durb = t_visc / (kinetic_en * cs_turb_cmu ) * rnnb * 0.22 ;
427 } else {
428 t_visc_durb = t_visc ;
429 }
430
431 *cofimp = 1. - *ypup * (2.0 * sqrt( l_visc / t_visc_durb * dup1 * (1.0 - dup1) ) - dup2) ;
432
433 /* log layer */
434 if (*yplus > ypluli) {
435 *nlogla += 1;
436 /* viscous sub-layer or buffer layer*/
437 } else {
438 *iuntur = 0;
439 *nsubla += 1;
440 }
441 }
442}
443
444/*----------------------------------------------------------------------------*/
464/*----------------------------------------------------------------------------*/
465
466inline static void
468 cs_real_t t_visc,
470 cs_real_t y,
471 cs_real_t kinetic_en,
472 int *iuntur,
473 cs_lnum_t *nsubla,
474 cs_lnum_t *nlogla,
475 cs_real_t *ustar,
476 cs_real_t *uk,
478 cs_real_t *ypup,
479 cs_real_t *cofimp)
480{
481 const double ypluli = cs_glob_wall_functions->ypluli;
482
483 double rcprod, ml_visc, Re, g;
484
485 /* Compute the friction velocity ustar */
486
487 /* Blending for very low values of k */
488 Re = sqrt(kinetic_en) * y / l_visc;
489 g = exp(-Re/11.);
490
491 *uk = sqrt( (1.-g) * cs_turb_cmu025 * cs_turb_cmu025 * kinetic_en
492 + g * l_visc * vel / y);
493
494 *yplus = *uk * y / l_visc;
495
496 /* log layer */
497 if (*yplus > ypluli) {
498
499 *ustar = vel / (log(*yplus) / cs_turb_xkappa + cs_turb_cstlog);
500 *ypup = *yplus / (log(*yplus) / cs_turb_xkappa + cs_turb_cstlog);
501 /* Mixing length viscosity */
502 ml_visc = cs_turb_xkappa * l_visc * *yplus;
503 rcprod = CS_MIN(cs_turb_xkappa, CS_MAX(1., sqrt(ml_visc / t_visc)) / *yplus);
504 *cofimp = 1. - *ypup / cs_turb_xkappa * ( 2. * rcprod - 1. / (2. * *yplus));
505
506 *nlogla += 1;
507
508 /* viscous sub-layer */
509 } else {
510
511 if (*yplus > 1.e-12) {
512 *ustar = fabs(vel / *yplus); /* FIXME remove that: its is here only to
513 be fully equivalent to the former code. */
514 } else {
515 *ustar = 0.;
516 }
517 *ypup = 1.;
518 *cofimp = 0.;
519
520 *iuntur = 0;
521 *nsubla += 1;
522
523 }
524}
525
526/*----------------------------------------------------------------------------*/
547/*----------------------------------------------------------------------------*/
548
549inline static void
551 cs_real_t t_visc,
553 cs_real_t y,
554 cs_real_t kinetic_en,
555 int *iuntur,
556 cs_lnum_t *nsubla,
557 cs_lnum_t *nlogla,
558 cs_real_t *ustar,
559 cs_real_t *uk,
561 cs_real_t *dplus,
562 cs_real_t *ypup,
563 cs_real_t *cofimp)
564{
565 CS_UNUSED(iuntur);
566
567 const double ypluli = cs_glob_wall_functions->ypluli;
568
569 double rcprod, ml_visc, Re, g;
570 /* Compute the friction velocity ustar */
571
572 /* Blending for very low values of k */
573 Re = sqrt(kinetic_en) * y / l_visc;
574 g = exp(-Re/11.);
575
576 *uk = sqrt( (1.-g) * cs_turb_cmu025 * cs_turb_cmu025 * kinetic_en
577 + g * l_visc * vel / y);
578
579 *yplus = *uk * y / l_visc;
580
581 /* Compute the friction velocity ustar */
582 *uk = cs_turb_cmu025 * sqrt(kinetic_en);
583 *yplus = *uk * y / l_visc;
584
585 /* Log layer */
586 if (*yplus > ypluli) {
587
588 *dplus = 0.;
589
590 *nlogla += 1;
591
592 /* Viscous sub-layer and therefore shift */
593 } else {
594
595 *dplus = ypluli - *yplus;
596 *yplus = ypluli;
597
598 /* Count the cell as if it was in the viscous sub-layer */
599 *nsubla += 1;
600
601 }
602
603 /* Mixing length viscosity */
604 ml_visc = cs_turb_xkappa * l_visc * *yplus;
605 rcprod = CS_MIN(cs_turb_xkappa, CS_MAX(1., sqrt(ml_visc / t_visc)) / *yplus);
606
607 *ustar = vel / (log(*yplus) / cs_turb_xkappa + cs_turb_cstlog);
608 *ypup = (*yplus - *dplus) / (log(*yplus) / cs_turb_xkappa + cs_turb_cstlog);
609 *cofimp = 1. - *ypup
610 / cs_turb_xkappa * (2. * rcprod - 1. / (2. * *yplus - *dplus));
611}
612
613/*----------------------------------------------------------------------------
614 * Compute u+ for a given yk+ between 0.1 and 200 according to the two
615 * scales wall functions using Van Driest mixing length.
616 * This function holds the coefficients of the polynome fitting log(u+).
617 *
618 * parameters:
619 * yplus <-- dimensionless distance
620 *
621 * returns:
622 * the resulting dimensionless velocity.
623 *----------------------------------------------------------------------------*/
624
625inline static cs_real_t
627{
628 /* Coefficients of the polynome fitting log(u+) for yk < 200 */
629 static double aa[11] = {-0.0091921, 3.9577, 0.031578,
630 -0.51013, -2.3254, -0.72665,
631 2.969, 0.48506, -1.5944,
632 0.087309, 0.1987 };
633
634 cs_real_t y1,y2,y3,y4,y5,y6,y7,y8,y9,y10, uplus;
635
636 y1 = 0.25 * log(yplus);
637 y2 = y1 * y1;
638 y3 = y2 * y1;
639 y4 = y3 * y1;
640 y5 = y4 * y1;
641 y6 = y5 * y1;
642 y7 = y6 * y1;
643 y8 = y7 * y1;
644 y9 = y8 * y1;
645 y10 = y9 * y1;
646
647 uplus = aa[0]
648 + aa[1] * y1
649 + aa[2] * y2
650 + aa[3] * y3
651 + aa[4] * y4
652 + aa[5] * y5
653 + aa[6] * y6
654 + aa[7] * y7
655 + aa[8] * y8
656 + aa[9] * y9
657 + aa[10] * y10;
658
659 return exp(uplus);
660}
661
662/*----------------------------------------------------------------------------*/
697/*----------------------------------------------------------------------------*/
698
699inline static void
701 cs_real_t l_visc,
703 cs_real_t y,
704 cs_real_t kinetic_en,
705 int *iuntur,
706 cs_lnum_t *nsubla,
707 cs_lnum_t *nlogla,
708 cs_real_t *ustar,
709 cs_real_t *uk,
711 cs_real_t *ypup,
712 cs_real_t *cofimp,
713 cs_real_t *lmk,
714 cs_real_t kr,
715 bool wf)
716{
717 double urplus, d_up, lmk15;
718
719 if (wf)
720 *uk = sqrt(sqrt((1.-cs_turb_crij2)/cs_turb_crij1 * rnnb * kinetic_en));
721
722 /* Set a low threshold value in case tangential velocity is zero */
723 *yplus = CS_MAX(*uk * y / l_visc, 1.e-4);
724
725 /* Dimensionless roughness */
726 cs_real_t krp = *uk * kr / l_visc;
727
728 /* Extension of Van Driest mixing length according to Rotta (1962) with
729 Cebeci & Chang (1978) correlation */
730 cs_real_t dyrp = 0.9 * (sqrt(krp) - krp * exp(-krp / 6.));
731 cs_real_t yrplus = *yplus + dyrp;
732
733 if (dyrp <= 1.e-1)
734 d_up = dyrp;
735 else if (dyrp <= 200.)
736 d_up = _vdriest_dupdyp_integral(dyrp);
737 else
738 d_up = 16.088739022054590 + log(dyrp/200.) / cs_turb_xkappa;
739
740 if (yrplus <= 1.e-1) {
741
742 urplus = yrplus;
743
744 if (wf) {
745 *iuntur = 0;
746 *nsubla += 1;
747
748 *lmk = 0.;
749
750 *ypup = 1.;
751
752 *cofimp = 0.;
753 }
754
755 } else if (yrplus <= 200.) {
756
757 urplus = _vdriest_dupdyp_integral(yrplus);
758
759 if (wf) {
760 *nlogla += 1;
761
762 *ypup = *yplus / (urplus-d_up);
763
764 /* Mixing length in y+ */
765 *lmk = cs_turb_xkappa * (*yplus) *(1-exp(- (*yplus) / cs_turb_vdriest));
766
767 /* Mixing length in 3/2*y+ */
768 lmk15 = cs_turb_xkappa * 1.5 * (*yplus) *(1-exp(- 1.5 * (*yplus)
769 / cs_turb_vdriest));
770
771 *cofimp = 1. - (2. / (1. + *lmk) - 1. / (1. + lmk15)) * *ypup;
772 }
773
774 } else {
775
776 urplus = 16.088739022054590 + log(yrplus/200) / cs_turb_xkappa;
777
778 if (wf) {
779 *nlogla += 1;
780
781 *ypup = *yplus / (urplus-d_up);
782
783 /* Mixing length in y+ */
784 *lmk = cs_turb_xkappa * (*yplus) *(1-exp(- (*yplus) / cs_turb_vdriest));
785
786 /* Mixing length in 3/2*y+ */
787 lmk15 = cs_turb_xkappa * 1.5 * (*yplus) *(1-exp(- 1.5 * (*yplus)
788 / cs_turb_vdriest));
789
790 *cofimp = 1. - (2. / *lmk - 1. / lmk15) * *ypup;
791 }
792
793 }
794
795 *ustar = vel / (urplus-d_up);
796}
797
798/*----------------------------------------------------------------------------*/
827/*----------------------------------------------------------------------------*/
828
829inline static void
831 cs_real_t t_visc,
833 cs_real_t y,
834 cs_real_t roughness,
835 cs_real_t kinetic_en,
836 int *iuntur,
837 cs_lnum_t *nsubla,
838 cs_lnum_t *nlogla,
839 cs_real_t *ustar,
840 cs_real_t *uk,
842 cs_real_t *dplus,
843 cs_real_t *ypup,
844 cs_real_t *cofimp)
845{
846 CS_UNUSED(iuntur);
847
848 const double ypluli = cs_glob_wall_functions->ypluli;
849
850 double rcprod, ml_visc, Re, g;
851
852 /* Compute the friction velocity ustar */
853
854 /* Shifting of the wall distance to be consistant with
855 * the fully rough wall function
856 *
857 * ln((y+y0)/y0) = ln((y+y0)/alpha xi) + kappa * 5.2
858 *
859 * y0 = roughness * exp(-kappa * 8.5)
860 */
861 double y0 = roughness*exp(-cs_turb_xkappa*cs_turb_cstlog_rough);
862
863 /* Blending for very low values of k */
864 Re = sqrt(kinetic_en) * (y + y0) / l_visc;
865 g = exp(-Re/11.);
866
867 *uk = sqrt( (1.-g) * cs_turb_cmu025 * cs_turb_cmu025 * kinetic_en
868 + g * l_visc * vel / (y + y0));
869
870 double effective_visc = (l_visc + cs_turb_cstlog_alpha * roughness * *uk);
871
872 /* NB: tends to "y/xi" in rough regime, to "y.uk/nu" in smooth regime */
873 *yplus = *uk * (y + y0) / effective_visc;
874
875 double yk = *uk * y / l_visc;
876
877 /* As for scalable wall functions, "y*uk/effective_visc+ dplus = yplus" */
878 *dplus = *uk * y0 / effective_visc;
879
880 /* Log layer and shifted with the roughness */
881 if (yk > ypluli) {
882
883 *nlogla += 1;
884
885 /* Viscous sub-layer and therefore shift again */
886 } else {
887
888 *dplus = ypluli - *yplus;
889 *yplus = ypluli;
890 /* Count the cell as if it was in the viscous sub-layer */
891 *nsubla += 1;
892
893 }
894
895 double uplus = log(*yplus) / cs_turb_xkappa + cs_turb_cstlog;
896 *ustar = vel / uplus;
897 *ypup = yk / uplus;
898
899 /* Mixing length viscosity, compatible with both regimes */
900 ml_visc = cs_turb_xkappa * *uk * (y + y0);
901 rcprod = CS_MIN(cs_turb_xkappa, CS_MAX(1., sqrt(ml_visc / t_visc)) / *yplus);
902 *cofimp = 1. - (*yplus - *dplus) / uplus
903 / cs_turb_xkappa * ( 2. * rcprod - 1. / (2. * *yplus - *dplus));
904
905}
906
907/*----------------------------------------------------------------------------*/
927/*----------------------------------------------------------------------------*/
928
929inline static void
931 cs_real_t t_visc,
933 cs_real_t y,
934 int *iuntur,
935 cs_lnum_t *nsubla,
936 cs_lnum_t *nlogla,
937 cs_real_t *ustar,
938 cs_real_t *uk,
940 cs_real_t *dplus,
941 cs_real_t *ypup,
942 cs_real_t *cofimp)
943{
944 CS_UNUSED(t_visc);
945 CS_UNUSED(nlogla);
946 CS_UNUSED(dplus);
947
948 const double ypluli = cs_glob_wall_functions->ypluli;
949
950 /* Compute the friction velocity ustar */
951
952 *ustar = sqrt(vel * l_visc / y);
953 *yplus = *ustar * y / l_visc;
954 *uk = *ustar;
955 *ypup = 1.;
956 *cofimp = 0.;
957 *iuntur = 0;
958
959 if (*yplus <= ypluli) {
960
961 /* Disable the wall funcion count the cell in the viscous sub-layer */
962 *nsubla += 1;
963
964 } else {
965
966 /* Count the cell as if it was in the viscous sub-layer */
967 *nsubla += 1;
968
969 }
970}
971
972/*----------------------------------------------------------------------------*/
998/*----------------------------------------------------------------------------*/
999
1000inline static void
1002 cs_real_t prt,
1004 cs_real_t dplus,
1005 cs_real_t *htur,
1006 cs_real_t *yplim)
1007{
1008 /* Local variables */
1009 double tplus;
1010 double beta2,a2;
1011 double yp2;
1012 double prlm1;
1013
1014 const double epzero = 1.e-12;
1015
1016 /*==========================================================================*/
1017
1018 /*==========================================================================
1019 1. Initializations
1020 ==========================================================================*/
1021
1022 /*==========================================================================*/
1023
1024 (*htur) = CS_MAX(yplus-dplus,epzero)/CS_MAX(yplus,epzero);
1025
1026 prlm1 = 0.1;
1027
1028 /*==========================================================================
1029 2. Compute htur for small Prandtl numbers
1030 ==========================================================================*/
1031
1032 if (prl <= prlm1) {
1033 (*yplim) = prt/(prl*cs_turb_xkappa);
1034 if (yplus > (*yplim)) {
1035 tplus = prl*(*yplim) + prt/cs_turb_xkappa * log(yplus/(*yplim));
1036 (*htur) = prl*(yplus-dplus)/tplus;
1037 }
1038
1039 /*========================================================================
1040 3. Compute htur for the model with three sub-layers
1041 ========================================================================*/
1042
1043 } else {
1044 yp2 = cs_turb_xkappa*1000./prt;
1045 yp2 = sqrt(yp2);
1046 (*yplim) = pow(1000./prl,1./3.);
1047
1048 a2 = 15.*pow(prl,2./3.);
1049 beta2 = a2 - 500./ pow(yp2,2);
1050
1051 if (yplus >= (*yplim) && yplus < yp2) {
1052 tplus = a2 - 500./(yplus*yplus);
1053 (*htur) = prl*(yplus-dplus)/tplus;
1054 }
1055
1056 if (yplus >= yp2) {
1057 tplus = beta2 + prt/cs_turb_xkappa*log(yplus/yp2);
1058 (*htur) = prl*(yplus-dplus)/tplus;
1059 }
1060
1061 }
1062}
1063
1064/*----------------------------------------------------------------------------*/
1087/*----------------------------------------------------------------------------*/
1088
1089inline static void
1091 cs_real_t prt,
1093 cs_real_t *htur)
1094{
1095 cs_real_t prlrat = prl / prt;
1096
1097 /* Parameters of the numerical quadrature */
1098 const int ninter_max = 100;
1099 const cs_real_t ypmax = 1.e2;
1100
1101 /* No correction for very small yplus */
1102 if (yplus <= 0.1)
1103 *htur = 1.;
1104 else {
1105 cs_real_t ypint = CS_MIN(yplus, ypmax);
1106
1107 /* The number of sub-intervals is taken proportional to yplus and equal to
1108 * ninter_max if yplus=ypmax */
1109
1110 int npeff = CS_MAX((int)(ypint / ypmax * (double)(ninter_max)), 1);
1111
1112 double dy = ypint / (double)(npeff);
1113 cs_real_t stplus = 0.;
1114 cs_real_t nut1 = 0.;
1115 cs_real_t nut2 = 0.;
1116
1117 for (int ip = 1; ip <= npeff; ip++) {
1118 double yp = ypint * (double)(ip) / (double)(npeff);
1119 nut2 = cs_turb_xkappa * yp * (1. - exp(-yp / cs_turb_vdriest));
1120 stplus += dy / (1. + prlrat * 0.5 * (nut1 + nut2));
1121 nut1 = nut2;
1122 }
1123
1124 if (yplus > ypint) {
1125 cs_real_t r = prlrat * cs_turb_xkappa;
1126 stplus += log( (1. + r*yplus) / (1. + r*ypint)) / r;
1127 }
1128
1129 if (stplus >= 1.e-6)
1130 *htur = yplus / stplus;
1131 else
1132 *htur = 1.;
1133 }
1134}
1135
1136/*============================================================================
1137 * Public function definitions for Fortran API
1138 *============================================================================*/
1139
1140/*----------------------------------------------------------------------------
1141 * Wrapper to cs_wall_functions_velocity.
1142 *----------------------------------------------------------------------------*/
1143
1144void CS_PROCF (wallfunctions, WALLFUNCTIONS)
1145(
1146 const cs_int_t *const iwallf,
1147 const cs_lnum_t *const ifac,
1148 const cs_real_t *const viscosity,
1149 const cs_real_t *const t_visc,
1150 const cs_real_t *const vel,
1151 const cs_real_t *const y,
1152 const cs_real_t *const roughness,
1153 const cs_real_t *const rnnb,
1154 const cs_real_t *const kinetic_en,
1155 cs_int_t *iuntur,
1156 cs_lnum_t *nsubla,
1157 cs_lnum_t *nlogla,
1158 cs_real_t *ustar,
1159 cs_real_t *uk,
1161 cs_real_t *ypup,
1162 cs_real_t *cofimp,
1163 cs_real_t *dplus
1164);
1165
1166/*----------------------------------------------------------------------------
1167 * Wrapper to cs_wall_functions_scalar.
1168 *----------------------------------------------------------------------------*/
1169
1170void CS_PROCF (hturbp, HTURBP)
1171(
1172 const cs_int_t *const iwalfs,
1173 const cs_real_t *const prl,
1174 const cs_real_t *const prt,
1175 const cs_real_t *const yplus,
1176 const cs_real_t *const dplus,
1177 cs_real_t *htur,
1178 cs_real_t *yplim
1179);
1180
1181/*=============================================================================
1182 * Public function prototypes
1183 *============================================================================*/
1184
1185/*----------------------------------------------------------------------------
1186 *! \brief Provide acces to cs_glob_wall_functions
1187 *
1188 * needed to initialize structure with GUI
1189 *----------------------------------------------------------------------------*/
1190
1193
1197/*-------------------------------------------------------------------------------
1198 Arguments
1199 ______________________________________________________________________________.
1200 mode name role !
1201 ______________________________________________________________________________*/
1227/*-------------------------------------------------------------------------------*/
1228
1229void
1231 cs_lnum_t ifac,
1232 cs_real_t l_visc,
1233 cs_real_t t_visc,
1234 cs_real_t vel,
1235 cs_real_t y,
1236 cs_real_t roughness,
1237 cs_real_t rnnb,
1238 cs_real_t kinetic_en,
1239 int *iuntur,
1240 cs_lnum_t *nsubla,
1241 cs_lnum_t *nlogla,
1242 cs_real_t *ustar,
1243 cs_real_t *uk,
1245 cs_real_t *ypup,
1246 cs_real_t *cofimp,
1247 cs_real_t *dplus);
1248
1249/*-------------------------------------------------------------------------------*/
1250
1264/*-------------------------------------------------------------------------------
1265 Arguments
1266 ______________________________________________________________________________.
1267 mode name role !
1268 ______________________________________________________________________________*/
1279/*-------------------------------------------------------------------------------*/
1280
1281void
1283 cs_real_t prl,
1284 cs_real_t prt,
1286 cs_real_t dplus,
1287 cs_real_t *htur,
1288 cs_real_t *yplim);
1289
1290/*----------------------------------------------------------------------------*/
1291
1293
1294#endif /* __CS_WALL_FUNCTIONS_H__ */
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:467
double cs_real_t
Floating-point value.
Definition cs_defs.h:302
#define CS_MIN(a, b)
Definition cs_defs.h:430
#define _(String)
Definition cs_defs.h:52
int cs_int_t
Fortran-compatible integer.
Definition cs_defs.h:301
#define CS_MAX(a, b)
Definition cs_defs.h:431
#define CS_PROCF(x, y)
Definition cs_defs.h:481
#define CS_UNUSED(x)
Definition cs_defs.h:453
#define END_C_DECLS
Definition cs_defs.h:468
int cs_lnum_t
local mesh entity id
Definition cs_defs.h:298
@ eps
Definition cs_field_pointer.h:71
@ vel
Definition cs_field_pointer.h:68
@ yplus
Definition cs_field_pointer.h:237
const double cs_turb_bpow
Definition cs_turbulence_model.c:420
const double cs_turb_cstlog_rough
Definition cs_turbulence_model.c:404
double cs_turb_cstlog_alpha
Definition cs_turbulence_model.c:414
double cs_turb_crij2
Definition cs_turbulence_model.c:484
const double cs_turb_vdriest
Definition cs_turbulence_model.c:382
const double cs_turb_cstlog
Definition cs_turbulence_model.c:393
double cs_turb_crij1
Definition cs_turbulence_model.c:478
double cs_turb_cmu025
Definition cs_turbulence_model.c:434
const double cs_turb_xkappa
Definition cs_turbulence_model.c:373
double cs_turb_dpow
Definition cs_turbulence_model.c:423
const double cs_turb_apow
Definition cs_turbulence_model.c:417
double cs_turb_cmu
Definition cs_turbulence_model.c:431
const cs_turb_model_t * cs_glob_turb_model
cs_wall_f_s_type_t
Definition cs_wall_functions.h:71
@ CS_WALL_F_S_VDRIEST
Definition cs_wall_functions.h:74
@ CS_WALL_F_S_ARPACI_LARSEN
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:314
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 / .
Definition cs_wall_functions.c:320
const cs_wall_functions_t * cs_glob_wall_functions
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.
Definition cs_wall_functions.h:930
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 fric...
Definition cs_wall_functions.h:358
static cs_real_t _vdriest_dupdyp_integral(cs_real_t yplus)
Definition cs_wall_functions.h:626
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 .
Definition cs_wall_functions.h:550
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.
Definition cs_wall_functions.h:700
cs_wall_functions_t * cs_get_glob_wall_functions(void)
Definition cs_wall_functions.c:280
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.
Definition cs_wall_functions.h:200
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:1090
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 v...
Definition cs_wall_functions.h:1001
cs_wall_f_type_t
Definition cs_wall_functions.h:58
@ CS_WALL_F_1SCALE_LOG
Definition cs_wall_functions.h:62
@ CS_WALL_F_1SCALE_POWER
Definition cs_wall_functions.h:61
@ CS_WALL_F_2SCALES_SMOOTH_ROUGH
Definition cs_wall_functions.h:66
@ CS_WALL_F_2SCALES_LOG
Definition cs_wall_functions.h:63
@ CS_WALL_F_2SCALES_VDRIEST
Definition cs_wall_functions.h:65
@ CS_WALL_F_DISABLED
Definition cs_wall_functions.h:60
@ CS_WALL_F_2SCALES_CONTINUOUS
Definition cs_wall_functions.h:67
@ CS_WALL_F_SCALABLE_2SCALES_LOG
Definition cs_wall_functions.h:64
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 k...
Definition cs_wall_functions.h:467
void hturbp(const cs_int_t *const iwalfs, const cs_real_t *const prl, const cs_real_t *const prt, 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:250
void wallfunctions(const cs_int_t *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 roughness, const cs_real_t *const rnnb, const cs_real_t *const kinetic_en, cs_int_t *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:202
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.
Definition cs_wall_functions.h:830
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 flo...
Definition cs_wall_functions.c:511
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.
Definition cs_wall_functions.h:129
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:290
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:614
integer(c_int), pointer, save iwalfs
Wall functions for scalar.
Definition optcal.f90:619
int itytur
Definition cs_turbulence_model.h:99
wall functions descriptor.
Definition cs_wall_functions.h:81
double ypluli
Definition cs_wall_functions.h:92
cs_wall_f_s_type_t iwalfs
Definition cs_wall_functions.h:85
int iwallt
Definition cs_wall_functions.h:87
cs_wall_f_type_t iwallf
Definition cs_wall_functions.h:83