1!-------------------------------------------------------------------------------
2
3! This file is part of Code_Saturne, a general-purpose CFD tool.
4!
5! Copyright (C) 1998-2021 EDF S.A.
6!
7! This program is free software; you can redistribute it and/or modify it under
8! the terms of the GNU General Public License as published by the Free Software
9! Foundation; either version 2 of the License, or (at your option) any later
10! version.
11!
12! This program is distributed in the hope that it will be useful, but WITHOUT
13! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14! FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
15! details.
16!
17! You should have received a copy of the GNU General Public License along with
18! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
19! Street, Fifth Floor, Boston, MA 02110-1301, USA.
20
21!-------------------------------------------------------------------------------
22
23!> \file cstphy.f90
24!> \brief Module for physical constants
25
26module cstphy
27
28  !=============================================================================
29
30  use, intrinsic :: iso_c_binding
31
32  use paramx
33
34  implicit none
35
36  !=============================================================================
37
38  !> \defgroup cstphy Module for physical constants
39
40  !> \addtogroup cstphy
41  !> \{
42
43  !> Temperature in Kelvin correponding to 0 degrees Celsius (= +273,15)
44  double precision :: tkelvi
45  parameter(tkelvi = 273.15d0)
46
47  !> Calories (1 cvar_al = xcal2j J)
48  double precision :: xcal2j
49  parameter(xcal2j = 4.1855d0)
50
51  !> Stephan constant for the radiative module \f$\sigma\f$
52  !> in \f$W.m^{-2}.K^{-4}\f$
53  double precision :: stephn
54  parameter(stephn = 5.6703d-8)
55
56  !> Perfect gas constant for air (mixture)
57  real(c_double), pointer, save :: rair
58
59  !> ratio gaz constant h2o/ dry air
60  real(c_double), pointer, save :: rvsra
61
62  !> latent heat of evaporation
63  real(c_double), pointer, save :: clatev
64
65  !> Boltzmann constant (\f$J.K^{-1}\f$)
66  double precision kboltz
67  parameter(kboltz = 1.38d-23)
68
69  !> Ideal gas constant (\f$J.mol^{-1}.K^{-1}\f$)
70  double precision cs_physical_constants_r
71  parameter(cs_physical_constants_r = 8.31446d0)
72
73  !> Gravity
74  real(c_double), pointer, save :: gx, gy, gz
75
76  !> Coriolis effects
77  integer(c_int), pointer, save :: icorio
78
79  !> Physical constants of the fluid
80  !> filling \ref xyzp0 indicator
81  integer(c_int), pointer, save :: ixyzp0
82
83  !> indicates if the isobaric specific heat \f$C_p\f$ is variable:
84  !>  - 0: constant, no property field is declared
85  !>  - 1: variable, \f$C_p\f$ is declared as a property field\n
86  !> When gas or coal combustion is activated, \ref icp is automatically set to 0
87  !> (constant \f$C_p\f$). With the electric module, it is automatically set to 1.
88  !> The user is not allowed to modify these default choices.\n
89  !> When \ref icp = 1 is specified, the code automatically modifies this value to
90  !> make \ref icp designate the effective index-number of the property "specific heat".
91  !> For each cell iel, the value of \f$C_p\f$ is then specified by the user in the
92  !> appropriate subroutine (\ref cs_user_physical_properties for the standard physics).\n
93  !> Useful if there is 1\f$\leqslant\f$N\f$\leqslant\f$\ref dimens::nscal "nscal" so that
94  !> iscsth(n)=1 (there is a scalar temperature) or with the compressible module
95  !> for non perfect gases.
96  integer(c_int), pointer, save :: icp
97
98  !> isochoric specific heat \f$ C_v \f$
99  integer(c_int), pointer, save :: icv
100
101  !> variable density field \f$ \rho \f$:
102  !>    - 1: true, its variation law be given either
103  !> in the GUI, or in the user subroutine
104  !> \ref cs_user_physical_properties .\n
105  !> See \subpage physical_properties for more informations.
106  !>    - 0: false, its value is the reference density
107  !> \ref ro0.
108  integer(c_int), pointer, save :: irovar
109
110  !> variable viscosity field \f$ \mu \f$:
111  !>    - 1: true, its variation law be given either
112  !> in the GUI, or in the user subroutine
113  !> \ref cs_user_physical_properties .\n
114  !> See \subpage physical_properties for more informations.
115  !>    - 0: false, its value is the reference molecular
116  !> dynamic viscosity \ref viscl0
117  integer(c_int), pointer, save :: ivivar
118
119  !> Sutherland law for laminar viscosity and thermal conductivity
120  !> Only useful in gas mix (igmix) specific physics
121  !> - 1: Sutherland law
122  !> - 0: low temperature law (linear except for helium)
123  integer(c_int), pointer, save :: ivsuth
124
125  !> reference density.\n
126  !>
127  !> Negative value: not initialized.
128  !> Its value is not used in gas or coal combustion modelling (it will be
129  !> calculated following the perfect gas law, with \f$P0\f$ and \f$T0\f$).
130  !> With the compressible module, it is also not used by the code,
131  !> but it may be (and often is) referenced by the user in user subroutines;
132  !> it is therefore better to specify its value.
133  !>
134  !> Always useful otherwise, even if a law defining the density is given by
135  !> the user subroutines \ref cs_user_physical_properties.
136  !> indeed, except with the
137  !> compressible module, CS  does not use
138  !> the total pressure \f$P\f$ when solving the Navier-Stokes equation, but a
139  !> reduced pressure .
140  !> \f$P^*=P-\rho_0\vect{g}.(\vect{x}-\vect{x}_0)+P^*_0-P_0\f$.
141  !> where
142  !> \f$\vect{x_0}\f$ is a reference point (see \ref xyzp0) and \f$P^*_0\f$ and \f$P_0\f$ are
143  !> reference values (see \ref pred0 and \ref p0). Hence, the term
144  !> \f$-\grad{P}+\rho\vect{g}\f$ in the equation is treated as
145  !> \f$-\grad{P^*}+(\rho-\rho_0)\vect{g}\f$. The closer \ref ro0 is to the value of \f$\rho\f$,
146  !> the more \f$P^*\f$ will tend to represent only the dynamic part of the pressure and
147  !> the faster and more precise its solution will be. Whatever the value of \ref ro0,
148  !> both \f$P\f$ and \f$P^*\f$ appear in the log and the post-processing outputs..
149  !> with the compressible module, the calculation is made directly on the total
150  !> pressure
151  real(c_double), pointer, save :: ro0
152
153  !> reference molecular dynamic viscosity.\n
154  !>
155  !> Negative value: not initialized.
156  !> Always useful, it is the used value unless the user specifies the
157  !> viscosity in the subroutine \ref cs_user_physical_properties
158  real(c_double), pointer, save :: viscl0
159
160  !> reference pressure for the total pressure.\n
161  !>
162  !> except with the compressible module, the total pressure \f$P\f$ is evaluated
163  !> from the reduced pressure \f$P^*\f$ so that \f$P\f$
164  !> is equal to \ref p0 at the reference position \f$\vect{x}_0\f$
165  !> (given by \ref xyzp0).
166  !> with the compressible module, the total pressure is solved directly.
167  !> always Useful
168  real(c_double), pointer, save :: p0
169
170  !> reference value for the reduced pressure \f$P^*\f$ (see \ref ro0).\n
171  !>
172  !> It is especially used to initialise the reduced pressure and as a reference
173  !> value for the outlet boundary conditions.
174  !> For an optimised precision in the resolution of \f$P^*\f$,
175  !>  it is wiser to keep \ref pred0 to 0.
176  !> With the compressible module, the "pressure" variable appearing in the
177  !> equations directly represents the total pressure.
178  !> It is therefore initialized to \ref p0 and not \ref pred0 (see \ref ro0).
179  !> Always useful, except with the compressible module
180  real(c_double), pointer, save :: pred0
181
182  !> coordinates of the reference point \f$\vect{x}_0\f$ for the total pressure.
183  !>  - When there are no Dirichlet conditions for the pressure (closed domain),
184  !> \ref xyzp0
185  !> does not need to be specified (unless the total pressure has a clear
186  !> physical meaning in the configuration treated).
187  !>  - When Dirichlet conditions on the pressure are specified but only through
188  !> stantard outlet conditions (as it is in most configurations),
189  !> \ref xyzp0 does not need to be specified by the user, since it will be
190  !> set to the coordinates of the reference outlet face
191  !> (i.e. the code will automatically select a
192  !> reference outlet boundary face and set \ref xyzp0 so that \f$P\f$ equals
193  !> \ref p0 at this face). Nonetheless, if \ref xyzp0 is specified by
194  !> the user, the calculation will remain correct.
195  !>  - When direct Dirichlet conditions are specified by the user (specific
196  !> value set on specific boundary faces), it is better to specify the
197  !> corresponding reference point (\em i.e. specify where the total pressure
198  !> is \ref p0). This way, the boundary conditions for the reduced pressure
199  !> will be close to \ref pred0, ensuring an optimal precision in the
200  !> resolution. If \ref xyzp0 is not specified, the reduced
201  !> pressure will be shifted, but the calculations will remain correct..
202  !>  - With the compressible module, the "pressure" variable appearing in the
203  !> equations directly represents the total pressure.
204  !> \ref xyzp0 is therefore not used..
205  !>
206  !> Always useful, except with the compressible module.
207  real(c_double), pointer, save :: xyzp0(:)
208
209  !> reference temperature.
210  !>
211  !> Useful for the specific physics gas or coal combustion (initialization
212  !> of the density), for the electricity modules to initialize the domain
213  !> temperature and for the compressible module (initializations).
214  !> It must be given in Kelvin.
215  real(c_double), pointer, save :: t0
216
217  !> Reference internal energy for the barotropic compressible module
218  double precision, save :: eint0
219
220  !> reference specific heat.
221  !>
222  !> Useful if there is 1 <= n <= nscaus,
223  !> so that \ref optcal::iscalt "iscalt" = n and \ref optcal::itherm "itherm" = 1
224  !> (there is a "temperature" scalar),
225  !> unless the user specifies the specific heat in the user subroutine
226  !> \ref cs_user_physical_properties (\ref cstphy::icp "icp" > 0) with the compressible module or
227  !>  coal combustion, \ref cp0 is also needed even when there is no user scalar.
228  !> \note None of the scalars from the specific physics is a temperature.
229  !> \note When using the Graphical Interface, \ref cp0 is also used to
230  !> calculate the diffusivity of the thermal scalars,
231  !> based on their conductivity; it is therefore needed, unless the
232  !> diffusivity is also specified in \ref cs_user_physical_properties.
233  real(c_double), pointer, save :: cp0
234
235  !> Reference isochoric specific heat.
236  !>
237  !> Useful for the compressible module (J/kg/K)
238  real(c_double), pointer, save :: cv0
239
240  !> Reference thermal conductivity.
241  !>
242  !> Useful with thermal model (W/m/K)
243  real(c_double), pointer, save :: lambda0
244
245  !> Molar mass of the perfect gas in \f$ kg/mol \f$
246  !> (if \ref cstphy::ieos "ieos"=1)
247  !>
248  !> Always useful
249  real(c_double), pointer, save :: xmasmr
250
251  !> Uniform variable thermodynamic pressure for the low-Mach algorithm
252  !>    - 1: true
253  !>    - 0: false
254  integer(c_int), pointer, save :: ipthrm
255
256  !> Thermodynamic pressure for the current time step
257  real(c_double), pointer, save :: pther
258
259  !> Thermodynamic pressure for the previous time step
260  real(c_double), pointer, save :: pthera
261
262  !>   pthermax: Thermodynamic maximum pressure for user clipping,
263  !>             used to model a venting effect
264  real(c_double), pointer, save :: pthermax
265
266  !> Leak surface
267  real(c_double), pointer, save :: sleak
268
269  !> Leak head loss (2.9 by default, from Idelcick)
270  real(c_double), pointer, save :: kleak
271
272  !> Initial reference density
273  real(c_double), pointer, save :: roref
274
275
276  !> \defgroup csttur Module for turbulence constants
277
278  !> \addtogroup csttur
279  !> \{
280
281  !> \f$ \kappa \f$ Karman constant. (= 0.42)
282  !> Useful if and only if \ref iturb >= 10.
283  !> (mixing length, \f$k-\varepsilon\f$, \f$R_{ij}-\varepsilon\f$,
284  !> LES, v2f or \f$k-\omega\f$)
285  double precision, save :: xkappa
286
287  !> constant of logarithmic law function:
288  !> \f$ \dfrac{1}{\kappa} \ln(y^+) + cstlog \f$
289  !>  (\f$ cstlog = 5.2 \f$)
290  !> constant of the logarithmic wall function.
291  !> Useful if and only if \ref iturb >= 10
292  !> (mixing length, \f$k-\varepsilon\f$, \f$R_{ij}-\varepsilon\f$,
293  !> LES, v2f or \f$k-\omega\f$)
294  double precision, save :: cstlog
295
296  !> limit value of \f$y^+\f$ for the viscous sublayer.
297  !> \ref ypluli depends on the chosen wall function: it is
298  !> initialized to 10.88 for the scalable wall function (\ref optcal::iwallf "iwallf"=4),
299  !> otherwise it is initialized to \f$1/\kappa\approx 2,38\f$.
300  !> In LES, \ref ypluli is taken by default to be 10.88.
301  !>
302  !> Always useful
303  real(c_double), pointer, save :: ypluli
304
305  !> Werner and Wengle coefficient
306  double precision, save :: apow
307
308  !> Werner and Wengle coefficient
309  double precision, save :: bpow
310
311  !> Werner and Wengle coefficient
312  double precision, save :: cpow
313
314  !> Werner and Wengle coefficient
315  double precision, save :: dpow
316
317  !> constant \f$C_\mu\f$ for all the RANS turbulence models
318  !> Warning, different values for the v2f model
319  !> Useful if and only if \ref iturb = 20, 21, 30, 31, 50, 51 or 60
320  !> (\f$k-\varepsilon\f$, \f$R_{ij}-\varepsilon\f$ or \f$k-\omega\f$)
321  real(c_double), pointer, save :: cmu
322
323  !> \f$ C_\mu^\frac{1}{4} \f$
324  real(c_double), pointer, save :: cmu025
325
326  !> constant \f$C_{\varepsilon 1}\f$ for all the RANS turbulence models except
327  !> for the v2f and the \f$k-\omega\f$ models.
328  !> Useful if and only if \ref iturb= 20,
329  !> 21, 30 or 31 (\f$k-\varepsilon\f$ or \f$R_{ij}-\varepsilon\f$)
330  double precision, save :: ce1
331
332  !> constant \f$C_{\varepsilon 2}\f$ for the \f$k-\varepsilon\f$ and
333  !> \f$R_{ij}-\varepsilon\f$ LRR models.
334  !> Useful if and only if \ref optcal::iturb "iturb"= 20, 21 or 30
335  !> (\f$k-\varepsilon\f$ or \f$R_{ij}-\varepsilon\f$ LRR)
336  double precision, save :: ce2
337
338  !> constant \f$C_{NL1}\f$ for the \f$k-\varepsilon\f$
339  !> model from Baglietto et al. (quadratric)
340  !> Useful if and only if \ref optcal::iturb "iturb"= 23
341  double precision, save :: cnl1
342
343  !> constant \f$C_{NL2}\f$ for the \f$k-\varepsilon\f$
344  !> model from Baglietto et al. (quadratric)
345  !> Useful if and only if \ref optcal::iturb "iturb"= 23
346  double precision, save :: cnl2
347
348  !> constant \f$C_{NL3}\f$ for the \f$k-\varepsilon\f$
349  !> model from Baglietto et al. (quadratric)
350  !> Useful if and only if \ref optcal::iturb "iturb"= 23
351  double precision, save :: cnl3
352
353  !> constant \f$C_{NL4}\f$ for the \f$k-\varepsilon\f$
354  !> model from Baglietto et al. (quadratric)
355  !> Useful if and only if \ref optcal::iturb "iturb"= 23
356  double precision, save :: cnl4
357
358  !> constant \f$C_{NL5}\f$ for the \f$k-\varepsilon\f$
359  !> model from Baglietto et al. (quadratric)
360  !> Useful if and only if \ref optcal::iturb "iturb"= 23
361  double precision, save :: cnl5
362
363  !> Coefficient of interfacial coefficient in k-eps,
364  !> used in Lagrange treatment
365  !>
366
367  !> constant \f$C_{\varepsilon 4}\f$ for the interfacial term
368  !> (Lagrangian module) in case of two-way coupling.
369  !> Useful in case of Lagrangian modelling,
370  !> in \f$k-\varepsilon\f$ and \f$R_{ij}-\varepsilon\f$ with two-way coupling.
371  double precision, save :: ce4
372
373  !> constant \f$C_1\f$ for the \f$R_{ij}-\varepsilon\f$ LRR model.
374  !> Useful if and only if \ref iturb=30
375  !> (\f$R_{ij}-\varepsilon\f$ LRR)
376  real(c_double), pointer, save :: crij1
377
378  !> constant \f$C_2\f$ for the \f$R_{ij}-\varepsilon\f$ LRR model.
379  !> Useful if and only if \ref iturb=30
380  !> (\f$R_{ij}-\varepsilon\f$ LRR)
381  real(c_double), pointer, save :: crij2
382
383  !> constant \f$C_3\f$ for the buoyant production term \f$R_{ij}-\varepsilon\f$
384  !>  models.
385  double precision, save :: crij3
386
387  !> constant \f$C_1^\prime\f$ for the \f$R_{ij}-\varepsilon\f$ LRR model,
388  !> corresponding to the wall echo terms.
389  !> Useful if and only if \ref iturb=30 and \ref optcal::irijec "irijec"=1
390  !> (\f$R_{ij}-\varepsilon\f$ LRR)
391  double precision, save :: crijp1
392
393  !> constant \f$C_2^\prime\f$ for the \f$R_{ij}-\varepsilon\f$ LRR model,
394  !> corresponding to the wall echo terms.
395  !> Useful if and only if \ref iturb=30 and \ref optcal::irijec "irijec"=1
396  !> (\f$R_{ij}-\varepsilon\f$ LRR)
397  double precision, save :: crijp2
398
399  !> constant \f$C_{\varepsilon 2}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
400  !> Useful if and only if \ref iturb=31
401  !> (\f$R_{ij}-\varepsilon\f$ SSG)
402  double precision, save :: cssge2
403
404  !> constant \f$C_{s1}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
405  !> Useful if and only if \ref iturb=31
406  !> (\f$R_{ij}-\varepsilon\f$ SSG)
407  double precision, save :: cssgs1
408
409  !> constant \f$C_{s2}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
410  !> Useful if and only if \ref iturb=31
411  !> (\f$R_{ij}-\varepsilon\f$ SSG)
412  double precision, save :: cssgs2
413
414  !> constant \f$C_{r1}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
415  !> Useful if and only if \ref iturb=31
416  !> (\f$R_{ij}-\varepsilon\f$ SSG)
417  double precision, save :: cssgr1
418
419  !> constant \f$C_{r2}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
420  !> Useful if and only if \ref iturb=31
421  !> (\f$R_{ij}-\varepsilon\f$ SSG)
422  double precision, save :: cssgr2
423
424  !> constant \f$C_{r3}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
425  !> Useful if and only if \ref iturb=31
426  !> (\f$R_{ij}-\varepsilon\f$ SSG)
427  double precision, save :: cssgr3
428
429  !> constant \f$C_{r4}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
430  !> Useful if and only if \ref iturb=31
431  !> (\f$R_{ij}-\varepsilon\f$ SSG)
432  double precision, save :: cssgr4
433
434
435  !> constant \f$C_{r1}\f$ for the \f$R_{ij}-\varepsilon\f$ SSG model.
436  !> Useful if and only if \ref iturb=31
437  !> (\f$R_{ij}-\varepsilon\f$ SSG)
438  double precision, save :: cssgr5
439
440  !> constant of the Rij-epsilon EBRSM
441  double precision, save :: cebms1
442
443  !> constant of the Rij-epsilon EBRSM
444  double precision, save :: cebms2
445
446  double precision, save :: cebmr1, cebmr2, cebmr3, cebmr4, cebmr5
447
448  !> constant \f$C_s\f$ for the \f$R_{ij}-\varepsilon\f$ models.
449  double precision, save :: csrij
450
451  !> constant of the Rij-epsilon EBRSM
452  double precision, save :: cebme2
453
454  !> constant of the Rij-epsilon EBRSM
455  double precision, save :: cebmmu
456
457  !> constant of the Rij-epsilon EBRSM
458  double precision, save :: xcl
459
460  !> constant in the expression of Ce1' for the Rij-epsilon EBRSM
461  double precision, save :: xa1
462
463  !> constant of the Rij-epsilon EBRSM
464  double precision, save :: xct
465
466  !> constant of the Rij-epsilon EBRSM
467  double precision, save :: xceta
468
469  !> specific constant of v2f "BL-v2k" (or phi-alpha)
470  double precision, save :: cpale1
471
472  !> specific constant of v2f "BL-v2k" (or phi-alpha)
473  double precision, save :: cpale2
474
475  !> specific constant of v2f "BL-v2k" (or phi-alpha)
476  double precision, save :: cpale3
477
478  !> specific constant of v2f "BL-v2k" (or phi-alpha)
479  double precision, save :: cpale4
480
481  !> specific constant of v2f "BL-v2k" (or phi-alpha)
482  double precision, save :: cpalc1
483
484  !> specific constant of v2f "BL-v2k" (or phi-alpha)
485  double precision, save :: cpalc2
486
487  !> specific constant of v2f "BL-v2k" (or phi-alpha)
488  double precision, save :: cpalct
489
490  !> specific constant of v2f "BL-v2k" (or phi-alpha)
491  double precision, save :: cpalcl
492
493  !> specific constant of v2f "BL-v2k" (or phi-alpha)
494  double precision, save :: cpalet
495
496  !> constant \f$\sigma_{k1}\f$ for the \f$k-\omega\f$ SST model.
497  !> Useful if and only if \ref iturb=60
498  double precision, save :: ckwsk1
499
500  !> constant \f$\sigma_{k2}\f$ for the \f$k-\omega\f$ SST model.
501  !> Useful if and only if \ref iturb=60
502  double precision, save :: ckwsk2
503
504  !> constant \f$\sigma_{\omega 1}\f$ for the \f$k-\omega\f$ SST model.
505  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
506  double precision, save :: ckwsw1
507
508  !> constant \f$\sigma_{\omega 2}\f$ for the \f$k-\omega\f$ SST model.
509  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
510  double precision, save :: ckwsw2
511
512  !> constant \f$\beta_1\f$ for the \f$k-\omega\f$ SST model.
513  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
514  double precision, save :: ckwbt1
515
516  !> constant \f$\beta_2\f$ for the \f$k-\omega\f$ SST model.
517  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
518  double precision, save :: ckwbt2
519
520  !> constant \f$ C_{DDES}\f$ for the hybrid \f$k-\omega\f$ SST model.
521  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST) and hybrid_turb=1
522  double precision, save :: cddes
523
524  !> constant \f$ C_{SAS}\f$ for the hybrid \f$k-\omega\f$ SST model.
525  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST) and hybrid_turb=3
526  double precision, save :: csas
527
528  !> constant \f$ C_{DDES}\f$ for the hybrid \f$k-\omega\f$ SST model.
529  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST) and hybrid_turb=3
530  double precision, save :: csas_eta2
531
532  !> \f$\frac{\beta_1}{C_\mu}-\frac{\kappa^2}{\sqrt{C_\mu}\sigma_{\omega 1}}\f$
533  !> constant \f$\gamma_1\f$ for the \f$k-\omega\f$ SST model.
534  !> Useful if and only if \ref iturb=60
535  !> (\f$k-\omega\f$ SST)
536  !> \warning: \f$\gamma_1\f$ is calculated before the call to
537  !> \ref usipsu. Hence, if \f$\beta_1\f$, \f$C_\mu\f$, \f$\kappa\f$
538  !> or \f$\sigma_{\omega 1}\f$ is modified in \ref usipsu,
539  !> \ref ckwgm1 must also be modified in accordance.
540  double precision, save :: ckwgm1
541
542  !> \f$\frac{\beta_2}{C_\mu}-\frac{\kappa^2}{\sqrt{C_\mu}\sigma_{\omega 2}}\f$
543  !> constant \f$\gamma_2\f$ for the \f$k-\omega\f$ SST model.
544  !> Useful if and only if \ref iturb=60
545  !> (\f$k-\omega\f$ SST)
546  !> \warning: \f$\gamma_2\f$ is calculated before the call to
547  !> \ref usipsu. Hence, if \f$\beta_2\f$, \f$C_\mu\f$, \f$\kappa\f$
548  !> or \f$\sigma_{\omega 2}\f$ is modified in \ref usipsu,
549  !> \ref ckwgm2 must also be modified in accordance.
550  double precision, save :: ckwgm2
551
552  !> specific constant of k-omega SST
553  !> constant \f$a_1\f$ for the \f$k-\omega\f$ SST model.
554  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
555  double precision, save :: ckwa1
556
557  !> constant \f$ c_1 \f$ for the \f$k-\omega\f$ SST model.
558  !> Useful if and only if \ref iturb=60 (\f$k-\omega\f$ SST)
559  !> specific constant of k-omega SST
560  double precision, save :: ckwc1
561
562  !> specific constant of Spalart-Allmaras
563  double precision, save :: csab1
564  !> specific constant of Spalart-Allmaras
565  double precision, save :: csab2
566
567  !> specific constant of Spalart-Allmaras
568  double precision, save :: csasig
569
570  !> specific constant of Spalart-Allmaras
571  double precision, save :: csav1
572  !> specific constant of Spalart-Allmaras
573  double precision, save :: csaw1
574  !> specific constant of Spalart-Allmaras
575  double precision, save :: csaw2
576  !> specific constant of Spalart-Allmaras
577  double precision, save :: csaw3
578
579  !> constant of the Spalart-Shur rotation/curvature correction
580  double precision, save :: cssr1
581  !> constant of the Spalart-Shur rotation/curvature correction
582  double precision, save :: cssr2
583  !> constant of the Spalart-Shur rotation/curvature correction
584  double precision, save :: cssr3
585
586  !> is a characteristic macroscopic
587  !> length of the domain, used for the initialization of the turbulence and
588  !> the potential clipping (with \ref optcal::iclkep "iclkep"=1)
589  !>  - Negative value: not initialized (the code then uses the cubic root of
590  !> the domain volume).
591  !>
592  !> Useful if and only if \ref iturb = 20, 21, 30, 31, 50 or 60 (RANS models)
593  real(c_double), pointer, save :: almax
594
595  !> the characteristic flow velocity,
596  !> used for the initialization of the turbulence.
597  !> Negative value: not initialized.
598  !>
599  !> Useful if and only if \ref iturb= 20, 21, 30, 31, 50 or 60 (RANS model)
600  !> and the turbulence is not initialized somewhere
601  !> else (restart file or subroutine \ref cs\_user\_initialization)
602  real(c_double), pointer, save :: uref
603
604  !> mixing length for the mixing length model
605  !>
606  !> Useful if and only if \ref iturb= 10 (mixing length)
607  real(c_double), pointer, save :: xlomlg
608
609  !> constant used in the definition of LES filtering diameter:
610  !> \f$ \delta = \text{xlesfl} . (\text{ales} . volume)^{\text{bles}}\f$
611  !> \ref xlesfl is a constant used to define, for
612  !> each cell \f$\omega_i\f$, the width of the (implicit) filter:
613  !> \f$\overline{\Delta}=xlesfl(ales*|\Omega_i|)^{bles}\f$\n
614  !> Useful if and only if \ref iturb = 40 or 41
615  real(c_double), pointer, save :: xlesfl
616
617  !> constant used to define, for each cell \f$Omega_i\f$,
618  !> the width of the (implicit) filter:
619  !>  - \f$\overline{\Delta}=xlesfl(ales*|Omega_i|)^{bles}\f$
620  !>
621  !> Useful if and only if \ref iturb = 40 or 41.
622  real(c_double), pointer, save :: ales
623
624  !> constant used to define, for each cell \f$Omega_i\f$,
625  !>
626  !> the width of the (implicit) filter:
627  !>  - \f$\overline{\Delta}=xlesfl(ales*|Omega_i|)^{bles}\f$
628  !>
629  !> Useful if and only if \ref iturb = 40 or 41
630  real(c_double), pointer, save :: bles
631
632  !> Smagorinsky constant used in the Smagorinsky model for LES.
633  !> The sub-grid scale viscosity is calculated by
634  !> \f$\displaystyle\mu_{sg}=
635  !> \rho C_{smago}^2\bar{\Delta}^2\sqrt{2\bar{S}_{ij}\bar{S}_{ij}}\f$
636  !> where \f$\bar{\Delta}\f$ is the width of the filter
637  !>  and \f$\bar{S}_{ij}\f$ the filtered strain rate.
638  !>
639  !> Useful if and only if \ref iturb = 40
640  !> \note In theory Smagorinsky constant is 0.18.
641  !> For a planar canal plan, 0.065 value is rather taken.
642  real(c_double), pointer, save :: csmago
643
644  !> ratio between
645  !> explicit and explicit filter width for a dynamic model
646  !> constant used to define, for each cell \f$\Omega_i\f$,
647  !> the width of the explicit filter used in the framework of
648  !> the LES dynamic model:
649  !> \f$\widetilde{\overline{\Delta}}=xlesfd\overline{\Delta}\f$.
650  !>
651  !> Useful if and only if \ref iturb = 41
652  real(c_double), pointer, save :: xlesfd
653
654  !> maximum allowed value for the variable \f$C\f$ appearing in the LES dynamic
655  !> model.
656  !> Any larger value yielded by the calculation
657  !> procedure of the dynamic model will be clipped to \f$ smagmx\f$.
658  !>
659  !> Useful if and only if \ref iturb = 41
660  real(c_double), pointer, save :: smagmx
661
662  !> minimum allowed value for the variable \f$C\f$ appearing in the LES dynamic
663  !> model.
664  !> Any smaller value yielded by the calculation
665  !> procedure of the dynamic model will be clipped to \f$ smagmn\f$.
666  !>
667  !> Useful if and only if \ref iturb = 41
668  real(c_double), pointer, save :: smagmn
669
670  !> van Driest constant appearing in the van Driest damping function
671  !> applied to the Smagorinsky constant:
672  !>  - \f$ (1-\exp^{(-y^+/cdries}) \f$.
673  !>
674  !> Useful if and only if \ref iturb = 40 or 41
675  real(c_double), pointer, save :: cdries
676
677  !> minimal control volume
678  double precision, save :: volmin
679  !> maximal control volume
680  double precision, save :: volmax
681  !> total domain volume
682  double precision, save :: voltot
683
684  !> constant \f$a_1\f$ for the v2f \f$\varphi\f$-model.
685  !> Useful if and only if \ref iturb=50
686  !> (v2f \f$\varphi\f$-model)
687  double precision, save :: cv2fa1
688
689  !> constant \f$C_{\varepsilon 2}\f$ for the v2f \f$\varphi\f$-model.
690  !> Useful if and only if \ref iturb=50
691  !> (v2f \f$\varphi\f$-model)
692  double precision, save :: cv2fe2
693
694  !> constant \f$C_1\f$ for the v2f \f$\varphi\f$-model.
695  !> Useful if and only if \ref iturb=50
696  !> (v2f \f$\varphi\f$-model)
697  double precision, save :: cv2fc1
698
699  !> constant \f$C_2\f$ for the v2f \f$\varphi\f$-model.
700  !> Useful if and only if \ref iturb=50
701  !> (v2f \f$\varphi\f$-model)
702  double precision, save :: cv2fc2
703
704  !> constant \f$C_T\f$ for the v2f \f$\varphi\f$-model.
705  !> Useful if and only if \ref iturb=50
706  !> (v2f \f$\varphi\f$-model)
707  double precision, save :: cv2fct
708
709  !> constant \f$C_L\f$ for the v2f \f$\varphi\f$-model.
710  !> Useful if and only if \ref iturb=50
711  !> (v2f \f$\varphi\f$-model)
712  double precision, save :: cv2fcl
713
714  !> constant \f$C_\eta\f$ for the v2f \f$\varphi\f$-model.
715  !> Useful if and only if \ref iturb=50
716  !> (v2f \f$\varphi\f$-model)
717  double precision, save :: cv2fet
718
719  !> constant of the WALE LES method
720  real(c_double), pointer, save :: cwale
721  !> coefficient of turbulent AFM flow model
722  double precision, save :: xiafm
723  !> coefficient of turbulent AFM flow model
724  double precision, save :: etaafm
725  !> coefficient of turbulent DFM flow model
726  double precision, save :: c1trit
727  !> coefficient of turbulent DFM flow model
728  double precision, save :: c2trit
729  !> coefficient of turbulent DFM flow model
730  double precision, save :: c3trit
731  !> coefficient of turbulent DFM flow model
732  double precision, save :: c4trit
733  !> constant of GGDH and AFM on the thermal scalar
734  double precision, save :: cthafm
735  !> constant of GGDH and AFM on the thermal scalar
736  double precision, save :: cthdfm
737  !> constant of EB-AFM and EB-DFM
738  double precision, save :: xclt
739  !> constant of EB-DFM
740  double precision, save :: rhebdfm
741  double precision, save :: cthebdfm
742
743  !> \}
744
745  !> \}
746
747  !=============================================================================
748
749  interface
750
751    !---------------------------------------------------------------------------
752
753    !> \cond DOXYGEN_SHOULD_SKIP_THIS
754
755    !---------------------------------------------------------------------------
756
757    ! Interface to C function retrieving pointers to members of the
758    ! global physical constants structure
759
760    subroutine cs_f_physical_constants_get_pointers(gx, gy, gz,              &
761                                                    icorio)                  &
762      bind(C, name='cs_f_physical_constants_get_pointers')
763      use, intrinsic :: iso_c_binding
764      implicit none
765      type(c_ptr), intent(out) :: gx, gy, gz, icorio
766    end subroutine cs_f_physical_constants_get_pointers
767
768    ! Interface to C function retrieving pointers to members of the
769    ! global fluid properties structure
770
771    subroutine cs_f_fluid_properties_get_pointers(ixyzp0,  &
772                                                  icp,     &
773                                                  icv,     &
774                                                  irovar,  &
775                                                  ivivar,  &
776                                                  ivsuth,  &
777                                                  ro0,     &
778                                                  viscl0,  &
779                                                  p0,      &
780                                                  pred0,   &
781                                                  xyzp0,   &
782                                                  t0,      &
783                                                  cp0,     &
784                                                  cv0,     &
785                                                  lambda0, &
786                                                  rair,    &
787                                                  rvsra,   &
788                                                  clatev,  &
789                                                  xmasmr,  &
790                                                  ipthrm,  &
791                                                  pther,   &
792                                                  pthera,  &
793                                                  pthermax,&
794                                                  sleak,   &
795                                                  kleak,   &
796                                                  roref)   &
797      bind(C, name='cs_f_fluid_properties_get_pointers')
798      use, intrinsic :: iso_c_binding
799      implicit none
800      type(c_ptr), intent(out) :: ixyzp0, icp, icv, irovar, ivivar, ivsuth
801      type(c_ptr), intent(out) :: ro0, viscl0, p0, pred0
802      type(c_ptr), intent(out) :: xyzp0, t0, cp0, cv0, lambda0
803      type(c_ptr), intent(out) :: rair, rvsra, clatev, xmasmr
804      type(c_ptr), intent(out) :: ipthrm
805      type(c_ptr), intent(out) :: pther, pthera, pthermax
806      type(c_ptr), intent(out) :: sleak, kleak, roref
807    end subroutine cs_f_fluid_properties_get_pointers
808
809    ! Interface to C function retrieving pointers to members of the
810    ! RANS turbulence model structure
811
812    subroutine cs_f_turb_reference_values(almax, uref, xlomlg) &
813      bind(C, name='cs_f_turb_reference_values')
814      use, intrinsic :: iso_c_binding
815      implicit none
816      type(c_ptr), intent(out) :: almax , uref  , xlomlg
817    end subroutine cs_f_turb_reference_values
818
819    ! Interface to C function retrieving pointers to members of the
820    ! global wall functions structure
821
822    subroutine cs_f_wall_reference_values(ypluli) &
823      bind(C, name='cs_f_wall_reference_values')
824      use, intrinsic :: iso_c_binding
825      implicit none
826      type(c_ptr), intent(out) :: ypluli
827    end subroutine cs_f_wall_reference_values
828
829    !---------------------------------------------------------------------------
830
831    ! Interface to C function completing the constants of the
832    ! turbulence model
833
834    subroutine cs_f_turb_complete_constants() &
835      bind(C, name='cs_turb_compute_constants')
836      use, intrinsic :: iso_c_binding
837    end subroutine cs_f_turb_complete_constants
838
839    !---------------------------------------------------------------------------
840
841    ! Interface to C function retrieving pointers to constants of the
842    ! turbulence model
843
844    subroutine cs_f_turb_model_constants_get_pointers(cmu, cmu025, crij1, crij2, &
845        csmago, xlesfd, smagmx, smagmn, cwale, xlesfl, ales, bles, cdries) &
846      bind(C, name='cs_f_turb_model_constants_get_pointers')
847      use, intrinsic :: iso_c_binding
848      implicit none
849      type(c_ptr), intent(out) :: cmu  , cmu025 , crij1 , crij2
850      type(c_ptr), intent(out) :: csmago, xlesfd, smagmx, smagmn
851      type(c_ptr), intent(out) :: cwale, xlesfl, ales, bles, cdries
852    end subroutine cs_f_turb_model_constants_get_pointers
853
854    !---------------------------------------------------------------------------
855
856    !> (DOXYGEN_SHOULD_SKIP_THIS) \endcond
857
858    !---------------------------------------------------------------------------
859
860  end interface
861
862  !=============================================================================
863
864contains
865
866  !=============================================================================
867
868  !> \brief Initialize Fortran physical constants API.
869  !> This maps Fortran pointers to global C structure members.
870
871  subroutine physical_constants_init
872
873    use, intrinsic :: iso_c_binding
874    implicit none
875
876    ! Local variables
877
878    type(c_ptr) :: c_gx, c_gy, c_gz, c_icorio
879
880    call cs_f_physical_constants_get_pointers(c_gx, c_gy, c_gz, c_icorio)
881
882    call c_f_pointer(c_gx, gx)
883    call c_f_pointer(c_gy, gy)
884    call c_f_pointer(c_gz, gz)
885    call c_f_pointer(c_icorio, icorio)
886
887  end subroutine physical_constants_init
888
889  !> \brief Initialize Fortran fluid properties API.
890  !> This maps Fortran pointers to global C structure members.
891
892  subroutine fluid_properties_init
893
894    use, intrinsic :: iso_c_binding
895    implicit none
896
897    ! Local variables
898
899    type(c_ptr) :: c_ixyzp0, c_icp, c_icv, c_irovar, c_ivivar
900    type(c_ptr) :: c_ivsuth, c_ro0, c_viscl0, c_p0
901    type(c_ptr) :: c_pred0, c_xyzp0, c_t0, c_cp0, c_cv0, c_lambda0
902    type(c_ptr) :: c_rair,c_rvsra,c_clatev, c_xmasmr
903    type(c_ptr) :: c_ipthrm
904    type(c_ptr) :: c_pther, c_pthera, c_pthermax
905    type(c_ptr) :: c_sleak, c_kleak, c_roref
906
907    call cs_f_fluid_properties_get_pointers(c_ixyzp0, c_icp, c_icv,         &
908                                            c_irovar, c_ivivar, c_ivsuth,   &
909                                            c_ro0, c_viscl0, c_p0, c_pred0, &
910                                            c_xyzp0, c_t0, c_cp0, c_cv0,    &
911                                            c_lambda0,                      &
912                                            c_rair,c_rvsra,c_clatev,        &
913                                            c_xmasmr,                       &
914                                            c_ipthrm, c_pther, c_pthera,    &
915                                            c_pthermax, c_sleak, c_kleak,   &
916                                            c_roref)
917
918
919    call c_f_pointer(c_ixyzp0, ixyzp0)
920    call c_f_pointer(c_icp, icp)
921    call c_f_pointer(c_icv, icv)
922    call c_f_pointer(c_irovar, irovar)
923    call c_f_pointer(c_ivivar, ivivar)
924    call c_f_pointer(c_ivsuth, ivsuth)
925    call c_f_pointer(c_ro0, ro0)
926    call c_f_pointer(c_viscl0, viscl0)
927    call c_f_pointer(c_p0, p0)
928    call c_f_pointer(c_pred0, pred0)
929    call c_f_pointer(c_xyzp0, xyzp0, [3])
930    call c_f_pointer(c_t0, t0)
931    call c_f_pointer(c_cp0, cp0)
932    call c_f_pointer(c_cv0, cv0)
933    call c_f_pointer(c_lambda0, lambda0)
934    call c_f_pointer(c_rair, rair)
935    call c_f_pointer(c_rvsra, rvsra)
936    call c_f_pointer(c_clatev, clatev)
937    call c_f_pointer(c_xmasmr, xmasmr)
938    call c_f_pointer(c_ipthrm, ipthrm)
939    call c_f_pointer(c_pther, pther)
940    call c_f_pointer(c_pthera, pthera)
941    call c_f_pointer(c_pthermax, pthermax)
942    call c_f_pointer(c_sleak, sleak)
943    call c_f_pointer(c_kleak, kleak)
944    call c_f_pointer(c_roref, roref)
945
946  end subroutine fluid_properties_init
947
948  !> \brief Initialize Fortran RANS turbulence model API.
949  !> This maps Fortran pointers to global C structure members.
950
951  subroutine turb_reference_values_init
952
953    use, intrinsic :: iso_c_binding
954    implicit none
955
956    ! Local variables
957
958    type(c_ptr) :: c_almax , c_uref  , c_xlomlg
959
960    call cs_f_turb_reference_values(c_almax, c_uref, c_xlomlg)
961
962    call c_f_pointer(c_almax , almax )
963    call c_f_pointer(c_uref  , uref  )
964    call c_f_pointer(c_xlomlg, xlomlg)
965
966  end subroutine turb_reference_values_init
967
968  !> \brief Initialize Fortran turbulence model constants.
969  !> This maps Fortran pointers to global C real numbers.
970
971  subroutine turb_model_constants_init
972
973    use, intrinsic :: iso_c_binding
974    implicit none
975
976    ! Local variables
977
978    type(c_ptr) :: c_cmu, c_cmu025, c_crij1, c_crij2
979    type(c_ptr) :: c_csmago, c_xlesfd, c_smagmx, c_smagmn, c_cwale
980    type(c_ptr) :: c_xlesfl, c_ales, c_bles, c_cdries
981
982    call cs_f_turb_model_constants_get_pointers(c_cmu, c_cmu025, c_crij1,    &
983                                                c_crij2, c_csmago, c_xlesfd, &
984                                                c_smagmx, c_smagmn, c_cwale, &
985                                                c_xlesfl, c_ales, c_bles,    &
986                                                c_cdries  )
987
988    call c_f_pointer(c_cmu    , cmu   )
989    call c_f_pointer(c_cmu025 , cmu025)
990    call c_f_pointer(c_crij1 , crij1)
991    call c_f_pointer(c_crij2 , crij2)
992    call c_f_pointer(c_csmago, csmago)
993    call c_f_pointer(c_xlesfd, xlesfd)
994    call c_f_pointer(c_smagmx, smagmx)
995    call c_f_pointer(c_smagmn, smagmn)
996    call c_f_pointer(c_cwale, cwale)
997    call c_f_pointer(c_xlesfl, xlesfl)
998    call c_f_pointer(c_ales  , ales  )
999    call c_f_pointer(c_bles  , bles  )
1000    call c_f_pointer(c_cdries, cdries)
1001
1002  end subroutine turb_model_constants_init
1003
1004  !=============================================================================
1005
1006end module cstphy
1007