Newer
Older
!MNH_LIC Copyright 1994-2020 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
!-----------------------------------------------------------------
! ####################
MODULE MODE_PRANDTL
USE PARKIND1, ONLY : JPRB
USE YOMHOOK , ONLY : LHOOK, DR_HOOK
! ####################
!
!* modification 08/2010 V. Masson smoothing of the discontinuity in functions
! used for implicitation of exchange coefficients
! 05/2020 V. Masson and C. Lac : bug in D_PHI3DTDZ2_O_DDTDZ
!
USE MODD_CTURB, ONLY : XCTV, XCSHF, XCTD, XPHI_LIM, XCPR3, XCPR4, XCPR5
USE MODD_PARAMETERS, ONLY : JPVEXT_TURB
!
IMPLICIT NONE
!----------------------------------------------------------------------------
CONTAINS
!----------------------------------------------------------------------------
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
SUBROUTINE PRANDTL(KKA,KKU,KKL,KRR,KRRI,OTURB_DIAG, &
HTURBDIM, &
TPFILE, &
PDXX,PDYY,PDZZ,PDZX,PDZY, &
PTHVREF,PLOCPEXNM,PATHETA,PAMOIST, &
PLM,PLEPS,PTKEM,PTHLM,PRM,PSVM,PSRCM, &
PREDTH1,PREDR1, &
PRED2TH3, PRED2R3, PRED2THR3, &
PREDS1,PRED2THS3, PRED2RS3, &
PBLL_O_E, &
PETHETA, PEMOIST )
! ###########################################################
!
!
!!**** *PRANDTL* - routine to compute the Prandtl turbulent numbers
!!
!! PURPOSE
!! -------
! The purpose of this routine is to compute the Redelsperger
! numbers and then get the turbulent Prandtl and Schmidt numbers:
! * for the heat fluxes - PHI3 = 1/ Prandtl
! * for the moisture fluxes - PSI3 = 1/ Schmidt
!
!!** METHOD
!! ------
!! The following steps are performed:
!!
!! 1 - default values of 1 are taken for phi3 and psi3 and different masks
!! are defined depending on the presence of turbulence, stratification and
!! humidity. The 1D Redelsperger numbers are computed
!! * ZREDTH1 : (g / THVREF ) (LT**2 / TKE ) ETHETA (D Theta / Dz)
!! * ZREDR1 : (g / THVREF ) (LT**2 / TKE ) EMOIST (D TW / Dz)
!! 2 - 3D Redelsperger numbers are computed only for turbulent
!! grid points where ZREDTH1 or ZREDR1 are > 0.
!! 3 - PHI3 is computed only for turbulent grid points where ZREDTH1 > 0
!! (turbulent thermally stratified points)
!! 4 - PSI3 is computed only for turbulent grid points where ZREDR1 > 0
!! (turbulent moist points)
!!
!!
!! EXTERNAL
!! --------
!! FUNCTIONs ETHETA and EMOIST :
!! allows to compute the coefficients
!! for the turbulent correlation between any variable
!! and the virtual potential temperature, of its correlations
!! with the conservative potential temperature and the humidity
!! conservative variable:
!! ------- ------- -------
!! A' Thv' = ETHETA A' Thl' + EMOIST A' Rnp'
!!
!! GX_M_M, GY_M_M, GZ_M_M : Cartesian gradient operators
!! MZM : Shuman function (mean operator in the z direction)
!! Module MODI_ETHETA : interface module for ETHETA
!! Module MODI_EMOIST : interface module for EMOIST
!! Module MODI_SHUMAN : interface module for Shuman operators
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!! Module MODD_CST : contains physical constants
!! XG : gravity constant
!!
!! Module MODD_CTURB: contains the set of constants for
!! the turbulence scheme
!! XCTV,XCPR2 : constants for the turbulent prandtl numbers
!! XTKEMIN : minimum value allowed for the TKE
!!
!! Module MODD_PARAMETERS
!! JPVEXT_TURB : number of vertical marginal points
!!
!! REFERENCE
!! ---------
!! Book 2 of documentation (routine PRANDTL)
!! Book 1 of documentation (Chapter: Turbulence)
!!
!! AUTHOR
!! ------
!! Joan Cuxart * INM and Meteo-France *
!!
!! MODIFICATIONS
!! -------------
!! Original 18/10/94
!! Modifications: Feb 14, 1995 (J.Cuxart and J.Stein)
!! Doctorization and Optimization
!! Modifications: March 21, 1995 (J.M. Carriere)
!! Introduction of cloud water
!! Modifications: March 21, 1995 (J. Cuxart and J.Stein)
!! Phi3 and Psi3 at w point + cleaning
!! Modifications: July 2, 1995 (J.Cuxart and Ph.Bougeault)
!! change the value of Phi3 and Psi3 if negative
!! Modifications: Sept 20, 1995 (J. Stein, J. Cuxart, J.L. Redelsperger)
!! remove the Where + use REDTH1+REDR1 for the tests
!! Modifications: October 10, 1995 (J. Cuxart and J.Stein)
!! Psi3 for tPREDS1he scalar variables
!! Modifications: February 27, 1996 (J.Stein) optimization
!! Modifications: June 15, 1996 (P.Jabouille) return to the previous
!! computation of Phi3 and Psi3
!! Modifications: October 10, 1996 (J. Stein) change the temporal
!! discretization
!! Modifications: May 23, 1997 (J. Stein) bug in 3D Redels number at ground
!! with orography
!! Modifications: Feb 20, 1998 (J. Stein) bug in all the 3D cases due to
!! the use of ZW1 instead of ZW2
!! Feb 20, 2003 (JP Pinty) Add PFRAC_ICE
!! July 2005 (Tomas, Masson) implicitation of PHI3 and PSI3
!! October 2009 (G. Tanguy) add ILENCH=LEN(YCOMMENT) after
!! change of YCOMMENT
!! 2012-02 Y. Seity, add possibility to run with reversed
!! vertical levels
!! Modifications: July 2015 (Wim de Rooy) LHARAT (Racmo turbulence) switch
!! 2017-09 J.Escobar, use epsilon XMNH_TINY_12 for R*4
!! Philippe Wautelet: 05/2016-04/2018: new data structures and calls for I/O
!! JL Redelsperger 03/2021 : adding Ocean case for temperature only
!! --------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
USE PARKIND1, ONLY : JPRB
USE YOMHOOK , ONLY : LHOOK, DR_HOOK
!
USE MODD_CST
USE MODD_CONF
USE MODD_CTURB
USE MODD_DYN_n, ONLY: LOCEAN
USE MODD_FIELD, ONLY: TFIELDDATA, TYPEREAL
USE MODD_IO, ONLY: TFILEDATA
USE MODD_PARAMETERS
!
USE MODI_GRADIENT_M
USE MODE_EMOIST, ONLY: EMOIST
USE MODE_ETHETA, ONLY: ETHETA
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
USE MODI_SHUMAN, ONLY: MZM
USE MODE_IO_FIELD_WRITE, ONLY: IO_FIELD_WRITE
!
IMPLICIT NONE
!
!
!* 0.1 declarations of arguments
!
INTEGER, INTENT(IN) :: KKA !near ground array index
INTEGER, INTENT(IN) :: KKU !uppest atmosphere array index
INTEGER, INTENT(IN) :: KKL !vert. levels type 1=MNH -1=ARO
INTEGER, INTENT(IN) :: KRR ! number of moist var.
INTEGER, INTENT(IN) :: KRRI ! number of ice var.
!
LOGICAL, INTENT(IN) :: OTURB_DIAG ! switch to write some
! diagnostic fields in the syncronous FM-file
CHARACTER(LEN=4), INTENT(IN) :: HTURBDIM ! Kind of turbulence param.
TYPE(TFILEDATA), INTENT(IN) :: TPFILE ! Output file
REAL, DIMENSION(:,:,:), INTENT(IN) :: PDXX,PDYY,PDZZ,PDZX,PDZY
! metric coefficients
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTHVREF ! Virtual Potential Temp.
! of the reference state
REAL, DIMENSION(:,:,:), INTENT(IN) :: PLOCPEXNM ! Lv(T)/Cp/Exner at t-1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PATHETA ! coefficients between
REAL, DIMENSION(:,:,:), INTENT(IN) :: PAMOIST ! s and Thetal and Rnp
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PLM ! Turbulent Mixing length
REAL, DIMENSION(:,:,:), INTENT(IN) :: PLEPS ! Dissipative length
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTHLM,PTKEM! Conservative Potential
! Temperature and TKE at t-1
REAL, DIMENSION(:,:,:,:), INTENT(IN) :: PRM ! Mixing ratios at t-1
! with PRM(:,:,:,1) = cons.
! mixing ratio
REAL, DIMENSION(:,:,:,:), INTENT(IN) :: PSVM ! Scalars at t-1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PSRCM
! s'r'c/2Sigma_s2 at t-1 multiplied by Lambda_3
!
!
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PREDTH1 ! Redelsperger number R_theta
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PREDR1 ! Redelsperger number R_q
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PRED2TH3 ! Redelsperger number R*2_theta
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PRED2R3 ! Redelsperger number R*2_q
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PRED2THR3! Redelsperger number R*2_thq
REAL, DIMENSION(:,:,:,:), INTENT(OUT):: PREDS1 ! Redelsperger number R_s
REAL, DIMENSION(:,:,:,:), INTENT(OUT):: PRED2THS3! Redelsperger number R*2_thsv
REAL, DIMENSION(:,:,:,:), INTENT(OUT):: PRED2RS3 ! Redelsperger number R*2_qsv
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PBLL_O_E! beta*Lk*Leps/tke
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PETHETA ! coefficient E_theta
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PEMOIST ! coefficient E_moist
!
!
! 0.2 declaration of local variables
!
REAL, DIMENSION(SIZE(PTHLM,1),SIZE(PTHLM,2),SIZE(PTHLM,3)) :: &
ZW1, ZW2, ZW3
! working variables
!
INTEGER :: IKB ! vertical index value for the first inner mass point
INTEGER :: IKE ! vertical index value for the last inner mass point
INTEGER :: IRESP ! Return code of FM routines
INTEGER :: ILENG ! Length of the data field in LFIFM file
INTEGER :: IGRID ! C-grid indicator in LFIFM file
INTEGER :: ILENCH ! Length of comment string in LFIFM file
CHARACTER (LEN=100) :: YCOMMENT ! comment string in LFIFM file
CHARACTER (LEN=16) :: YRECFM ! Name of the desired field in LFIFM file
INTEGER:: ISV ! number of scalar variables
INTEGER:: JSV ! loop index for the scalar variables
INTEGER :: JLOOP
REAL :: ZMINVAL
TYPE(TFIELDDATA) :: TZFIELD
! ---------------------------------------------------------------------------
!
!* 1. DEFAULT VALUES, 1D REDELSPERGER NUMBERS
! ----------------------------------------
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('PRANDTL',0,ZHOOK_HANDLE)
IF (LHARAT) THEN
PREDTH1(:,:,:)=0.
PREDR1(:,:,:)=0.
PRED2TH3(:,:,:)=0.
PRED2R3(:,:,:)=0.
PRED2THR3(:,:,:)=0.
PREDS1(:,:,:,:)=0.
PRED2THS3(:,:,:,:)=0.
PRED2RS3(:,:,:,:)=0.
PBLL_O_E(:,:,:)=0.
ENDIF
!
IKB = KKA+JPVEXT_TURB*KKL
IKE = KKU-JPVEXT_TURB*KKL
ILENG=SIZE(PTHLM,1)*SIZE(PTHLM,2)*SIZE(PTHLM,3)
ISV =SIZE(PSVM,4)
!
PETHETA(:,:,:) = MZM(ETHETA(KRR,KRRI,PTHLM,PRM,PLOCPEXNM,PATHETA,PSRCM), KKA, KKU, KKL)
PEMOIST(:,:,:) = MZM(EMOIST(KRR,KRRI,PTHLM,PRM,PLOCPEXNM,PAMOIST,PSRCM), KKA, KKU, KKL)
PETHETA(:,:,KKA) = 2.*PETHETA(:,:,IKB) - PETHETA(:,:,IKB+KKL)
PEMOIST(:,:,KKA) = 2.*PEMOIST(:,:,IKB) - PEMOIST(:,:,IKB+KKL)
!
!---------------------------------------------------------------------------
IF (.NOT. LHARAT) THEN
!
! 1.3 1D Redelsperger numbers
!
PBLL_O_E(:,:,:) = MZM(XG / PTHVREF(:,:,:) * PLM(:,:,:) * PLEPS(:,:,:) / PTKEM(:,:,:), KKA, KKU, KKL)
IF (KRR /= 0) THEN ! moist case
PREDTH1(:,:,:)= XCTV*PBLL_O_E(:,:,:) * PETHETA(:,:,:) * &

RODIER Quentin
committed
& GZ_M_W(KKA, KKU, KKL,PTHLM,PDZZ)
PREDR1(:,:,:) = XCTV*PBLL_O_E(:,:,:) * PEMOIST(:,:,:) * &

RODIER Quentin
committed
& GZ_M_W(KKA, KKU, KKL,PRM(:,:,:,1),PDZZ)
ELSE ! dry case

RODIER Quentin
committed
PREDTH1(:,:,:)= XCTV*PBLL_O_E(:,:,:) * GZ_M_W(KKA, KKU, KKL,PTHLM,PDZZ)
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
PREDR1(:,:,:) = 0.
END IF
!
! 3. Limits on 1D Redelperger numbers
! --------------------------------
!
ZMINVAL = (1.-1./XPHI_LIM)
!
ZW1 = 1.
ZW2 = 1.
!
WHERE (PREDTH1+PREDR1<-ZMINVAL)
ZW1 = (-ZMINVAL) / (PREDTH1+PREDR1)
END WHERE
!
WHERE (PREDTH1<-ZMINVAL)
ZW2 = (-ZMINVAL) / (PREDTH1)
END WHERE
ZW2 = MIN(ZW1,ZW2)
!
ZW1 = 1.
WHERE (PREDR1<-ZMINVAL)
ZW1 = (-ZMINVAL) / (PREDR1)
END WHERE
ZW1 = MIN(ZW2,ZW1)
!
!
! 3. Modification of Mixing length and dissipative length
! ----------------------------------------------------
!
PBLL_O_E(:,:,:) = PBLL_O_E(:,:,:) * ZW1(:,:,:)
PREDTH1 (:,:,:) = PREDTH1 (:,:,:) * ZW1(:,:,:)
PREDR1 (:,:,:) = PREDR1 (:,:,:) * ZW1(:,:,:)
!
! 4. Threshold for very small (in absolute value) Redelperger numbers
! ----------------------------------------------------------------
!
ZW2=SIGN(1.,PREDTH1(:,:,:))
PREDTH1(:,:,:)= ZW2(:,:,:) * MAX(1.E-30, ZW2(:,:,:)*PREDTH1(:,:,:))
!
IF (KRR /= 0) THEN ! dry case
ZW2=SIGN(1.,PREDR1(:,:,:))
PREDR1(:,:,:)= ZW2(:,:,:) * MAX(1.E-30, ZW2(:,:,:)*PREDR1(:,:,:))
END IF
!
!
!---------------------------------------------------------------------------
!
! For the scalar variables
DO JSV=1,ISV

RODIER Quentin
committed
PREDS1(:,:,:,JSV)=XCTV*PBLL_O_E(:,:,:)*GZ_M_W(KKA, KKU, KKL,PSVM(:,:,:,JSV),PDZZ)
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
END DO
!
DO JSV=1,ISV
ZW2=SIGN(1.,PREDS1(:,:,:,JSV))
PREDS1(:,:,:,JSV)= ZW2(:,:,:) * MAX(1.E-30, ZW2(:,:,:)*PREDS1(:,:,:,JSV))
END DO
!
!---------------------------------------------------------------------------
!
!* 2. 3D REDELSPERGER NUMBERS
! ------------------------
!
IF(HTURBDIM=='1DIM') THEN ! 1D case
!
!
PRED2TH3(:,:,:) = PREDTH1(:,:,:)**2
!
PRED2R3(:,:,:) = PREDR1(:,:,:) **2
!
PRED2THR3(:,:,:) = PREDTH1(:,:,:) * PREDR1(:,:,:)
!
ELSE IF (L2D) THEN ! 3D case in a 2D model
!
IF (KRR /= 0) THEN ! moist 3D case
PRED2TH3(:,:,:)= PREDTH1(:,:,:)**2+(XCTV*PBLL_O_E(:,:,:)*PETHETA(:,:,:) )**2 * &
MZM(GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2TH3(:,:,IKB)=PRED2TH3(:,:,IKB+KKL)
!
PRED2R3(:,:,:)= PREDR1(:,:,:)**2 + (XCTV*PBLL_O_E(:,:,:)*PEMOIST(:,:,:))**2 * &
MZM(GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2R3(:,:,IKB)=PRED2R3(:,:,IKB+KKL)
!
PRED2THR3(:,:,:)= PREDR1(:,:,:) * PREDTH1(:,:,:) + XCTV**2*PBLL_O_E(:,:,:)**2 * &
PEMOIST(:,:,:) * PETHETA(:,:,:) * &
MZM(GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL), KKA, KKU, KKL)
PRED2THR3(:,:,IKB)=PRED2THR3(:,:,IKB+KKL)
!
ELSE ! dry 3D case in a 2D model
PRED2TH3(:,:,:) = PREDTH1(:,:,:)**2 + XCTV**2*PBLL_O_E(:,:,:)**2 * &
MZM(GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2TH3(:,:,IKB)=PRED2TH3(:,:,IKB+KKL)
!
PRED2R3(:,:,:) = 0.
!
PRED2THR3(:,:,:) = 0.
!
END IF
!
ELSE ! 3D case in a 3D model
!
IF (KRR /= 0) THEN ! moist 3D case
PRED2TH3(:,:,:)= PREDTH1(:,:,:)**2 + ( XCTV*PBLL_O_E(:,:,:)*PETHETA(:,:,:) )**2 * &
MZM(GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2 &
+ GY_M_M(PTHLM,PDYY,PDZZ,PDZY, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2TH3(:,:,IKB)=PRED2TH3(:,:,IKB+KKL)
!
PRED2R3(:,:,:)= PREDR1(:,:,:)**2 + (XCTV*PBLL_O_E(:,:,:)*PEMOIST(:,:,:))**2 * &
MZM(GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2 + &
GY_M_M(PRM(:,:,:,1),PDYY,PDZZ,PDZY, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2R3(:,:,IKB)=PRED2R3(:,:,IKB+KKL)
!
PRED2THR3(:,:,:)= PREDR1(:,:,:) * PREDTH1(:,:,:) + XCTV**2*PBLL_O_E(:,:,:)**2 * &
PEMOIST(:,:,:) * PETHETA(:,:,:) * &
MZM(GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL)+ &
GY_M_M(PRM(:,:,:,1),PDYY,PDZZ,PDZY, KKA, KKU, KKL)* &
GY_M_M(PTHLM,PDYY,PDZZ,PDZY, KKA, KKU, KKL), KKA, KKU, KKL)
PRED2THR3(:,:,IKB)=PRED2THR3(:,:,IKB+KKL)
!
ELSE ! dry 3D case in a 3D model
PRED2TH3(:,:,:) = PREDTH1(:,:,:)**2 + XCTV**2*PBLL_O_E(:,:,:)**2 * &
MZM(GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL)**2 &
+ GY_M_M(PTHLM,PDYY,PDZZ,PDZY, KKA, KKU, KKL)**2, KKA, KKU, KKL)
PRED2TH3(:,:,IKB)=PRED2TH3(:,:,IKB+KKL)
!
PRED2R3(:,:,:) = 0.
!
PRED2THR3(:,:,:) = 0.
!
END IF
!
END IF ! end of the if structure on the turbulence dimensionnality
!
!
!---------------------------------------------------------------------------
!
! 5. Prandtl numbers for scalars
! ---------------------------
DO JSV=1,ISV
!
IF(HTURBDIM=='1DIM') THEN
! 1D case
PRED2THS3(:,:,:,JSV) = PREDS1(:,:,:,JSV) * PREDTH1(:,:,:)
IF (KRR /= 0) THEN
PRED2RS3(:,:,:,JSV) = PREDR1(:,:,:) *PREDS1(:,:,:,JSV)
ELSE
PRED2RS3(:,:,:,JSV) = 0.
END IF
!
ELSE IF (L2D) THEN ! 3D case in a 2D model
!
IF (KRR /= 0) THEN
ZW1 = MZM((XG / PTHVREF * PLM * PLEPS / PTKEM)**2, KKA, KKU, KKL) *PETHETA
ELSE
ZW1 = MZM((XG / PTHVREF * PLM * PLEPS / PTKEM)**2, KKA, KKU, KKL)
END IF
PRED2THS3(:,:,:,JSV) = PREDTH1(:,:,:) * PREDS1(:,:,:,JSV) + &
ZW1* &
MZM(GX_M_M(PSVM(:,:,:,JSV),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL), &
KKA, KKU, KKL)
!
IF (KRR /= 0) THEN
PRED2RS3(:,:,:,JSV) = PREDR1(:,:,:) * PREDS1(:,:,:,JSV) + &
ZW1 * PEMOIST * &
MZM(GX_M_M(PSVM(:,:,:,JSV),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL), &
KKA, KKU, KKL)
ELSE
PRED2RS3(:,:,:,JSV) = 0.
END IF
!
ELSE ! 3D case in a 3D model
!
IF (KRR /= 0) THEN
ZW1 = MZM((XG / PTHVREF * PLM * PLEPS / PTKEM)**2, KKA, KKU, KKL) *PETHETA
ELSE
ZW1 = MZM((XG / PTHVREF * PLM * PLEPS / PTKEM)**2, KKA, KKU, KKL)
END IF
PRED2THS3(:,:,:,JSV) = PREDTH1(:,:,:) * PREDS1(:,:,:,JSV) + &
ZW1* &
MZM(GX_M_M(PSVM(:,:,:,JSV),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PTHLM,PDXX,PDZZ,PDZX, KKA, KKU, KKL) &
+GY_M_M(PSVM(:,:,:,JSV),PDYY,PDZZ,PDZY, KKA, KKU, KKL)* &
GY_M_M(PTHLM,PDYY,PDZZ,PDZY, KKA, KKU, KKL), &
KKA, KKU, KKL)
!
IF (KRR /= 0) THEN
PRED2RS3(:,:,:,JSV) = PREDR1(:,:,:) * PREDS1(:,:,:,JSV) + &
ZW1 * PEMOIST * &
MZM(GX_M_M(PSVM(:,:,:,JSV),PDXX,PDZZ,PDZX, KKA, KKU, KKL)* &
GX_M_M(PRM(:,:,:,1),PDXX,PDZZ,PDZX, KKA, KKU, KKL) &
+GY_M_M(PSVM(:,:,:,JSV),PDYY,PDZZ,PDZY, KKA, KKU, KKL)* &
GY_M_M(PRM(:,:,:,1),PDYY,PDZZ,PDZY, KKA, KKU, KKL), &
KKA, KKU, KKL)
ELSE
PRED2RS3(:,:,:,JSV) = 0.
END IF
!
END IF ! end of HTURBDIM if-block
!
END DO
!
!---------------------------------------------------------------------------
!
!* 6. SAVES THE REDELSPERGER NUMBERS
! ------------------------------
!
IF ( OTURB_DIAG .AND. TPFILE%LOPENED ) THEN
!
! stores the RED_TH1
TZFIELD%CMNHNAME = 'RED_TH1'
TZFIELD%CSTDNAME = ''
TZFIELD%CLONGNAME = 'RED_TH1'
TZFIELD%CUNITS = '1'
TZFIELD%CDIR = 'XY'
TZFIELD%CCOMMENT = 'X_Y_Z_RED_TH1'
TZFIELD%NGRID = 4
TZFIELD%NTYPE = TYPEREAL
TZFIELD%NDIMS = 3
TZFIELD%LTIMEDEP = .TRUE.
CALL IO_Field_write(TPFILE,TZFIELD,PREDTH1)
!
! stores the RED_R1
TZFIELD%CMNHNAME = 'RED_R1'
TZFIELD%CSTDNAME = ''
TZFIELD%CLONGNAME = 'RED_R1'
TZFIELD%CUNITS = '1'
TZFIELD%CDIR = 'XY'
TZFIELD%CCOMMENT = 'X_Y_Z_RED_R1'
TZFIELD%NGRID = 4
TZFIELD%NTYPE = TYPEREAL
TZFIELD%NDIMS = 3
TZFIELD%LTIMEDEP = .TRUE.
CALL IO_Field_write(TPFILE,TZFIELD,PREDR1)
!
! stores the RED2_TH3
TZFIELD%CMNHNAME = 'RED2_TH3'
TZFIELD%CSTDNAME = ''
TZFIELD%CLONGNAME = 'RED2_TH3'
TZFIELD%CUNITS = '1'
TZFIELD%CDIR = 'XY'
TZFIELD%CCOMMENT = 'X_Y_Z_RED2_TH3'
TZFIELD%NGRID = 4
TZFIELD%NTYPE = TYPEREAL
TZFIELD%NDIMS = 3
TZFIELD%LTIMEDEP = .TRUE.
CALL IO_Field_write(TPFILE,TZFIELD,PRED2TH3)
!
! stores the RED2_R3
TZFIELD%CMNHNAME = 'RED2_R3'
TZFIELD%CSTDNAME = ''
TZFIELD%CLONGNAME = 'RED2_R3'
TZFIELD%CUNITS = '1'
TZFIELD%CDIR = 'XY'
TZFIELD%CCOMMENT = 'X_Y_Z_RED2_R3'
TZFIELD%NGRID = 4
TZFIELD%NTYPE = TYPEREAL
TZFIELD%NDIMS = 3
TZFIELD%LTIMEDEP = .TRUE.
CALL IO_Field_write(TPFILE,TZFIELD,PRED2R3)
!
! stores the RED2_THR3
TZFIELD%CMNHNAME = 'RED2_THR3'
TZFIELD%CSTDNAME = ''
TZFIELD%CLONGNAME = 'RED2_THR3'
TZFIELD%CUNITS = '1'
TZFIELD%CDIR = 'XY'
TZFIELD%CCOMMENT = 'X_Y_Z_RED2_THR3'
TZFIELD%NGRID = 4
TZFIELD%NTYPE = TYPEREAL
TZFIELD%NDIMS = 3
TZFIELD%LTIMEDEP = .TRUE.
CALL IO_Field_write(TPFILE,TZFIELD,PRED2THR3)
!
END IF
!
!---------------------------------------------------------------------------
ENDIF ! (Done only if LHARAT is FALSE)
!
IF (LHOOK) CALL DR_HOOK('PRANDTL',1,ZHOOK_HANDLE)
END SUBROUTINE PRANDTL
!
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
SUBROUTINE SMOOTH_TURB_FUNCT(PPHI3,PF_LIM,PF)
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPHI3 ! Phi3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PF_LIM ! Value of F when Phi3 is
! ! larger than Phi_lim
REAL, DIMENSION(:,:,:), INTENT(INOUT) :: PF ! function F to smooth
!
REAL, DIMENSION(SIZE(PF,1),SIZE(PF,2),SIZE(PF,3)) :: ZCOEF
!
!* adds a artificial correction to smooth the function near the discontinuity
! point at Phi3 = Phi_lim
! This smoothing is applied between 0.9*phi_lim (=2.7) and Phi_lim (=3)
! Note that in the Boundary layer, phi is usually between 0.8 and 1
!
!
ZCOEF = MAX(MIN(( 10.*(1.-PPHI3/XPHI_LIM)) ,1.), 0.)
!
PF(:,:,:) = ZCOEF(:,:,:) * PF &
+ (1.-ZCOEF(:,:,:)) * PF_LIM
!
END SUBROUTINE SMOOTH_TURB_FUNCT
!----------------------------------------------------------------------------
FUNCTION PHI3(PREDTH1,PREDR1,PRED2TH3,PRED2R3,PRED2THR3,HTURBDIM,OUSERV)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2TH3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2R3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2THR3
CHARACTER(len=4), INTENT(IN) :: HTURBDIM ! 1DIM or 3DIM turb. scheme
LOGICAL, INTENT(IN) :: OUSERV ! flag to use vapor
REAL, DIMENSION(SIZE(PREDTH1,1),SIZE(PREDTH1,2),SIZE(PREDTH1,3)) :: PHI3
!
REAL, DIMENSION(SIZE(PREDTH1,1),SIZE(PREDTH1,2),SIZE(PREDTH1,3)) :: ZW1, ZW2
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:PHI3',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PREDTH1,3)-JPVEXT_TURB
!
IF (HTURBDIM=='3DIM') THEN
!* 3DIM case
IF (OUSERV) THEN
ZW1(:,:,:) = 1. + 1.5* (PREDTH1(:,:,:)+PREDR1(:,:,:)) + &
( 0.5 * (PREDTH1(:,:,:)**2+PREDR1(:,:,:)**2) &
+ PREDTH1(:,:,:) * PREDR1(:,:,:) &
)
ZW2(:,:,:) = 0.5 * (PRED2TH3(:,:,:)-PRED2R3(:,:,:))
PHI3(:,:,:)= 1. - &
( ( (1.+PREDR1(:,:,:)) * &
(PRED2THR3(:,:,:) + PRED2TH3(:,:,:)) / PREDTH1(:,:,:) &
) + ZW2(:,:,:) &
) / ZW1(:,:,:)
ELSE
ZW1(:,:,:) = 1. + 1.5* PREDTH1(:,:,:) + &
0.5* PREDTH1(:,:,:)**2
ZW2(:,:,:) = 0.5* PRED2TH3(:,:,:)
PHI3(:,:,:)= 1. - &
(PRED2TH3(:,:,:) / PREDTH1(:,:,:) + ZW2(:,:,:)) / ZW1(:,:,:)
END IF
WHERE( PHI3 <= 0. .OR. PHI3 > XPHI_LIM )
PHI3 = XPHI_LIM
END WHERE
ELSE
!* 1DIM case
IF (OUSERV) THEN
PHI3(:,:,:)= 1./(1.+PREDTH1(:,:,:)+PREDR1(:,:,:))
ELSE
PHI3(:,:,:)= 1./(1.+PREDTH1(:,:,:))
END IF
END IF
!
PHI3(:,:,IKB-1)=PHI3(:,:,IKB)
PHI3(:,:,IKE+1)=PHI3(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:PHI3',1,ZHOOK_HANDLE)
END FUNCTION PHI3
!----------------------------------------------------------------------------
FUNCTION PSI_SV(PREDTH1,PREDR1,PREDS1,PRED2THS,PRED2RS,PPHI3,PPSI3)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:,:), INTENT(IN) :: PREDS1
REAL, DIMENSION(:,:,:,:), INTENT(IN) :: PRED2THS
REAL, DIMENSION(:,:,:,:), INTENT(IN) :: PRED2RS
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPHI3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPSI3
REAL, DIMENSION(SIZE(PRED2THS,1),SIZE(PRED2THS,2),SIZE(PRED2THS,3),SIZE(PRED2THS,4)) :: PSI_SV
!
INTEGER :: IKB, IKE
INTEGER :: JSV
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:PSI_SV',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PREDTH1,3)-JPVEXT_TURB
!
DO JSV=1,SIZE(PSI_SV,4)
PSI_SV(:,:,:,JSV) = ( 1. &
- (XCPR3+XCPR5) * (PRED2THS(:,:,:,JSV)/PREDS1(:,:,:,JSV)-PREDTH1) &
- (XCPR4+XCPR5) * (PRED2RS (:,:,:,JSV)/PREDS1(:,:,:,JSV)-PREDR1 ) &
- XCPR3 * PREDTH1 * PPHI3 - XCPR4 * PREDR1 * PPSI3 &
) / ( 1. + XCPR5 * ( PREDTH1 + PREDR1 ) )
! control of the PSI_SV positivity
WHERE ( (PSI_SV(:,:,:,JSV) <=0.).AND. (PREDTH1+PREDR1) <= 0. )
PSI_SV(:,:,:,JSV)=XPHI_LIM
END WHERE
PSI_SV(:,:,:,JSV) = MAX( 1.E-4, MIN(XPHI_LIM,PSI_SV(:,:,:,JSV)) )
!
PSI_SV(:,:,IKB-1,JSV)=PSI_SV(:,:,IKB,JSV)
PSI_SV(:,:,IKE+1,JSV)=PSI_SV(:,:,IKE,JSV)
END DO
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:PSI_SV',1,ZHOOK_HANDLE)
END FUNCTION PSI_SV
!----------------------------------------------------------------------------
FUNCTION D_PHI3DTDZ_O_DDTDZ(PPHI3,PREDTH1,PREDR1,PRED2TH3,PRED2THR3,HTURBDIM,OUSERV)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPHI3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2TH3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2THR3
CHARACTER(len=4), INTENT(IN) :: HTURBDIM ! 1DIM or 3DIM turb. scheme
LOGICAL, INTENT(IN) :: OUSERV ! flag to use vapor
REAL, DIMENSION(SIZE(PREDTH1,1),SIZE(PREDTH1,2),SIZE(PREDTH1,3)) :: D_PHI3DTDZ_O_DDTDZ
INTEGER :: IKB, IKE,JL,JK
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DTDZ_O_DDTDZ',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PREDTH1,3)-JPVEXT_TURB
!
IF (HTURBDIM=='3DIM') THEN
!* 3DIM case
IF (OUSERV) THEN
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:) &
* (1. - PREDTH1(:,:,:) * (3./2.+PREDTH1+PREDR1) &
/((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))) &
+ (1.+PREDR1)*(PRED2THR3+PRED2TH3) &
/ (PREDTH1*(1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1))) &
- (1./2.*PREDTH1+PREDR1 * (1.+PREDTH1+PREDR1)) &
/ ((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))
ELSEWHERE
D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:)
ENDWHERE
!
ELSE
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:) &
* (1. - PREDTH1(:,:,:) * (3./2.+PREDTH1) &
/((1.+PREDTH1)*(1.+1./2.*PREDTH1))) &
+ PRED2TH3 / (PREDTH1*(1.+PREDTH1)*(1.+1./2.*PREDTH1)) &
- 1./2.*PREDTH1 / ((1.+PREDTH1)*(1.+1./2.*PREDTH1))
ELSEWHERE
D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:)
ENDWHERE
!
END IF
ELSE
!* 1DIM case
! WHERE (PPHI3(:,:,:)/=XPHI_LIM)
! D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:) &
! * (1. - PREDTH1(:,:,:)*PPHI3(:,:,:))
! ELSEWHERE
! D_PHI3DTDZ_O_DDTDZ(:,:,:) = PPHI3(:,:,:)
! ENDWHERE
DO JL=1,SIZE(PPHI3,1)
DO JK=1,SIZE(PPHI3,3)
IF ( ABS(PPHI3(JL,1,JK)-XPHI_LIM) < 1.E-12 ) THEN
D_PHI3DTDZ_O_DDTDZ(JL,1,JK)=PPHI3(JL,1,JK)*&
& (1. - PREDTH1(JL,1,JK)*PPHI3(JL,1,JK))
ELSE
D_PHI3DTDZ_O_DDTDZ(JL,1,JK)=PPHI3(JL,1,JK)
ENDIF
ENDDO
ENDDO
END IF
!
!
D_PHI3DTDZ_O_DDTDZ(:,:,IKB-1)=D_PHI3DTDZ_O_DDTDZ(:,:,IKB)
D_PHI3DTDZ_O_DDTDZ(:,:,IKE+1)=D_PHI3DTDZ_O_DDTDZ(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DTDZ_O_DDTDZ',1,ZHOOK_HANDLE)
END FUNCTION D_PHI3DTDZ_O_DDTDZ
!----------------------------------------------------------------------------
FUNCTION D_PHI3DRDZ_O_DDRDZ(PPHI3,PREDTH1,PREDR1,PRED2TH3,PRED2THR3,HTURBDIM,OUSERV)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPHI3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2TH3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2THR3
CHARACTER(len=4), INTENT(IN) :: HTURBDIM ! 1DIM or 3DIM turb. scheme
LOGICAL, INTENT(IN) :: OUSERV ! flag to use vapor
REAL, DIMENSION(SIZE(PREDTH1,1),SIZE(PREDTH1,2),SIZE(PREDTH1,3)) :: D_PHI3DRDZ_O_DDRDZ
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DRDZ_O_DDRDZ',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PREDTH1,3)-JPVEXT_TURB
!
!
IF (HTURBDIM=='3DIM') THEN
!* 3DIM case
IF (OUSERV) THEN
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DRDZ_O_DDRDZ(:,:,:) = &
PPHI3(:,:,:) * (1.-PREDR1(:,:,:)*(3./2.+PREDTH1+PREDR1) &
/ ((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))) &
- PREDR1(:,:,:) * (PRED2THR3+PRED2TH3) / (PREDTH1 &
* (1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1))) &
+ PREDR1(:,:,:) * (1./2.+PREDTH1+PREDR1) &
/ ((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))
ELSEWHERE
D_PHI3DRDZ_O_DDRDZ(:,:,:) = PPHI3(:,:,:)
END WHERE
ELSE
D_PHI3DRDZ_O_DDRDZ(:,:,:) = PPHI3(:,:,:)
END IF
ELSE
!* 1DIM case
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DRDZ_O_DDRDZ(:,:,:) = PPHI3(:,:,:) &
* (1. - PREDR1(:,:,:)*PPHI3(:,:,:))
ELSEWHERE
D_PHI3DRDZ_O_DDRDZ(:,:,:) = PPHI3(:,:,:)
END WHERE
END IF
!
!
D_PHI3DRDZ_O_DDRDZ(:,:,IKB-1)=D_PHI3DRDZ_O_DDRDZ(:,:,IKB)
D_PHI3DRDZ_O_DDRDZ(:,:,IKE+1)=D_PHI3DRDZ_O_DDRDZ(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DRDZ_O_DDRDZ',1,ZHOOK_HANDLE)
END FUNCTION D_PHI3DRDZ_O_DDRDZ
!----------------------------------------------------------------------------
FUNCTION D_PHI3DTDZ2_O_DDTDZ(PPHI3,PREDTH1,PREDR1,PRED2TH3,PRED2THR3,PDTDZ,HTURBDIM,OUSERV)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PPHI3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2TH3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PRED2THR3
REAL, DIMENSION(:,:,:), INTENT(IN) :: PDTDZ
CHARACTER(len=4), INTENT(IN) :: HTURBDIM ! 1DIM or 3DIM turb. scheme
LOGICAL, INTENT(IN) :: OUSERV ! flag to use vapor
REAL, DIMENSION(SIZE(PREDTH1,1),SIZE(PREDTH1,2),SIZE(PREDTH1,3)) :: D_PHI3DTDZ2_O_DDTDZ
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DTDZ2_O_DDTDZ',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PREDTH1,3)-JPVEXT_TURB
!
!
IF (HTURBDIM=='3DIM') THEN
!* 3DIM case
IF (OUSERV) THEN
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:) &
* PDTDZ(:,:,:)*(2.-PREDTH1(:,:,:)*(3./2.+PREDTH1+PREDR1) &
/((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))) &
+ (1.+PREDR1)*(PRED2THR3+PRED2TH3) &
/ (PREDTH1*(1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1))) &
- (1./2.*PREDTH1+PREDR1 * (1.+PREDTH1+PREDR1)) &
/ ((1.+PREDTH1+PREDR1)*(1.+1./2.*(PREDTH1+PREDR1)))
ELSEWHERE
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:) * 2. * PDTDZ(:,:,:)
ENDWHERE
!
ELSE
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:) &
* PDTDZ(:,:,:)*(2.-PREDTH1(:,:,:)*(3./2.+PREDTH1) &
/((1.+PREDTH1)*(1.+1./2.*PREDTH1))) &
+ PRED2TH3 / (PREDTH1*(1.+PREDTH1)*(1.+1./2.*PREDTH1)) &
- 1./2.*PREDTH1 / ((1.+PREDTH1)*(1.+1./2.*PREDTH1))
ELSEWHERE
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:) * 2. * PDTDZ(:,:,:)
ENDWHERE
END IF
ELSE
!* 1DIM case
WHERE (PPHI3(:,:,:)/=XPHI_LIM)
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:)*PDTDZ(:,:,:) &
* (2. - PREDTH1(:,:,:)*PPHI3(:,:,:))
ELSEWHERE
D_PHI3DTDZ2_O_DDTDZ(:,:,:) = PPHI3(:,:,:) * 2. * PDTDZ(:,:,:)
END WHERE
END IF
!
!
D_PHI3DTDZ2_O_DDTDZ(:,:,IKB-1)=D_PHI3DTDZ2_O_DDTDZ(:,:,IKB)
D_PHI3DTDZ2_O_DDTDZ(:,:,IKE+1)=D_PHI3DTDZ2_O_DDTDZ(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_PHI3DTDZ2_O_DDTDZ',1,ZHOOK_HANDLE)
END FUNCTION D_PHI3DTDZ2_O_DDTDZ
!----------------------------------------------------------------------------
FUNCTION M3_WTH_WTH2(PREDTH1,PREDR1,PD,PBLL_O_E,PETHETA)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PBLL_O_E
REAL, DIMENSION(:,:,:), INTENT(IN) :: PETHETA
REAL, DIMENSION(SIZE(PD,1),SIZE(PD,2),SIZE(PD,3)) :: M3_WTH_WTH2
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_WTH2',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PD,3)-JPVEXT_TURB
M3_WTH_WTH2(:,:,:) = XCSHF*PBLL_O_E*PETHETA*0.5/XCTD &
* (1.+0.5*PREDTH1+PREDR1) / PD
M3_WTH_WTH2(:,:,IKB-1)=M3_WTH_WTH2(:,:,IKB)
M3_WTH_WTH2(:,:,IKE+1)=M3_WTH_WTH2(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_WTH2',1,ZHOOK_HANDLE)
END FUNCTION M3_WTH_WTH2
!----------------------------------------------------------------------------
FUNCTION D_M3_WTH_WTH2_O_DDTDZ(PM3_WTH_WTH2,PREDTH1,PREDR1,PD,PBLL_O_E,PETHETA)
REAL, DIMENSION(:,:,:), INTENT(IN) :: PM3_WTH_WTH2
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PBLL_O_E
REAL, DIMENSION(:,:,:), INTENT(IN) :: PETHETA
REAL, DIMENSION(SIZE(PD,1),SIZE(PD,2),SIZE(PD,3)) :: D_M3_WTH_WTH2_O_DDTDZ
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_M3_WTH_WTH2_O_DDTDZ',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PD,3)-JPVEXT_TURB
D_M3_WTH_WTH2_O_DDTDZ(:,:,:) = ( 0.5*XCSHF*PBLL_O_E*PETHETA*0.5/XCTD/PD &
- PM3_WTH_WTH2/PD*(1.5+PREDTH1+PREDR1) )&
* PBLL_O_E * PETHETA * XCTV
!
D_M3_WTH_WTH2_O_DDTDZ(:,:,IKB-1)=D_M3_WTH_WTH2_O_DDTDZ(:,:,IKB)
D_M3_WTH_WTH2_O_DDTDZ(:,:,IKE+1)=D_M3_WTH_WTH2_O_DDTDZ(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_M3_WTH_WTH2_O_DDTDZ',1,ZHOOK_HANDLE)
END FUNCTION D_M3_WTH_WTH2_O_DDTDZ
!----------------------------------------------------------------------------
FUNCTION M3_WTH_W2TH(KKA,KKU,KKL,PREDTH1,PREDR1,PD,PKEFF,PTKE)
INTEGER, INTENT(IN) :: KKA
INTEGER, INTENT(IN) :: KKU
INTEGER, INTENT(IN) :: KKL
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PKEFF
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTKE
REAL, DIMENSION(SIZE(PD,1),SIZE(PD,2),SIZE(PD,3)) :: M3_WTH_W2TH
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_W2TH',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PD,3)-JPVEXT_TURB
M3_WTH_W2TH(:,:,:) = XCSHF*PKEFF*1.5/MZM(PTKE, KKA, KKU, KKL) &
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
* (1. - 0.5*PREDR1*(1.+PREDR1)/PD ) / (1.+PREDTH1)
!
M3_WTH_W2TH(:,:,IKB-1)=M3_WTH_W2TH(:,:,IKB)
M3_WTH_W2TH(:,:,IKE+1)=M3_WTH_W2TH(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_W2TH',1,ZHOOK_HANDLE)
END FUNCTION M3_WTH_W2TH
!----------------------------------------------------------------------------
FUNCTION D_M3_WTH_W2TH_O_DDTDZ(KKA,KKU,KKL,PREDTH1,PREDR1,PD,PBLL_O_E,PETHETA,PKEFF,PTKE)
INTEGER, INTENT(IN) :: KKA
INTEGER, INTENT(IN) :: KKU
INTEGER, INTENT(IN) :: KKL
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PBLL_O_E
REAL, DIMENSION(:,:,:), INTENT(IN) :: PETHETA
REAL, DIMENSION(:,:,:), INTENT(IN) :: PKEFF
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTKE
REAL, DIMENSION(SIZE(PD,1),SIZE(PD,2),SIZE(PD,3)) :: D_M3_WTH_W2TH_O_DDTDZ
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_M3_WTH_W2TH_O_DDTDZ',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PD,3)-JPVEXT_TURB
D_M3_WTH_W2TH_O_DDTDZ(:,:,:) = &
- XCSHF*PKEFF*1.5/MZM(PTKE, KKA, KKU, KKL)/(1.+PREDTH1)**2*XCTV*PBLL_O_E*PETHETA &
* (1. - 0.5*PREDR1*(1.+PREDR1)/PD*( 1.+(1.+PREDTH1)*(1.5+PREDR1+PREDTH1)/PD) )
!
D_M3_WTH_W2TH_O_DDTDZ(:,:,IKB-1)=D_M3_WTH_W2TH_O_DDTDZ(:,:,IKB)
D_M3_WTH_W2TH_O_DDTDZ(:,:,IKE+1)=D_M3_WTH_W2TH_O_DDTDZ(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:D_M3_WTH_W2TH_O_DDTDZ',1,ZHOOK_HANDLE)
END FUNCTION D_M3_WTH_W2TH_O_DDTDZ
!----------------------------------------------------------------------------
FUNCTION M3_WTH_W2R(KKA,KKU,KKL,PD,PKEFF,PTKE,PBLL_O_E,PEMOIST,PDTDZ)
INTEGER, INTENT(IN) :: KKA
INTEGER, INTENT(IN) :: KKU
INTEGER, INTENT(IN) :: KKL
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PKEFF
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTKE
REAL, DIMENSION(:,:,:), INTENT(IN) :: PBLL_O_E
REAL, DIMENSION(:,:,:), INTENT(IN) :: PEMOIST
REAL, DIMENSION(:,:,:), INTENT(IN) :: PDTDZ
REAL, DIMENSION(SIZE(PD,1),SIZE(PD,2),SIZE(PD,3)) :: M3_WTH_W2R
INTEGER :: IKB, IKE
!
REAL(KIND=JPRB) :: ZHOOK_HANDLE
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_W2R',0,ZHOOK_HANDLE)
IKB = 1+JPVEXT_TURB
IKE = SIZE(PD,3)-JPVEXT_TURB
M3_WTH_W2R(:,:,:) = - XCSHF*PKEFF*0.75*XCTV*PBLL_O_E/MZM(PTKE, KKA, KKU, KKL)*PEMOIST*PDTDZ/PD
!
M3_WTH_W2R(:,:,IKB-1)=M3_WTH_W2R(:,:,IKB)
M3_WTH_W2R(:,:,IKE+1)=M3_WTH_W2R(:,:,IKE)
!
IF (LHOOK) CALL DR_HOOK('MODE_PRANDTL:M3_WTH_W2R',1,ZHOOK_HANDLE)
END FUNCTION M3_WTH_W2R
!----------------------------------------------------------------------------
FUNCTION D_M3_WTH_W2R_O_DDTDZ(KKA,KKU,KKL,PREDTH1,PREDR1,PD,PKEFF,PTKE,PBLL_O_E,PEMOIST)
INTEGER, INTENT(IN) :: KKA
INTEGER, INTENT(IN) :: KKU
INTEGER, INTENT(IN) :: KKL
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDTH1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PREDR1
REAL, DIMENSION(:,:,:), INTENT(IN) :: PD
REAL, DIMENSION(:,:,:), INTENT(IN) :: PKEFF
REAL, DIMENSION(:,:,:), INTENT(IN) :: PTKE
REAL, DIMENSION(:,:,:), INTENT(IN) :: PBLL_O_E