source: LMDZ6/branches/blowing_snow/libf/phylmd/cdrag_mod.F90 @ 5456

Last change on this file since 5456 was 4481, checked in by evignon, 21 months ago

quelques petites modifs pour le prochain atelier tke (lundi 03/04/2023)

  • Property svn:executable set to *
File size: 13.8 KB
Line 
1!
2MODULE cdrag_mod
3!
4! This module contains some procedures for calculation of the cdrag
5! coefficients for turbulent diffusion at surface
6!
7  IMPLICIT NONE
8
9CONTAINS
10!
11!****************************************************************************************
12!
13!r original routine svn3623
14!
15 SUBROUTINE cdrag(knon,  nsrf,   &
16     speed, t1,    q1,    zgeop1, &
17     psol,  tsurf, qsurf, z0m, z0h,  &
18     ri_in, iri_in, &
19     cdm,  cdh,  zri,   pref)
20
21  USE dimphy
22  USE indice_sol_mod
23  USE print_control_mod, ONLY: lunout, prt_level
24  USE ioipsl_getin_p_mod, ONLY : getin_p
25  USE atke_turbulence_ini_mod, ONLY : ric, cinf, cepsilon, pr_slope, pr_asym, pr_neut
26
27  IMPLICIT NONE
28! ================================================================= c
29!
30! Objet : calcul des cdrags pour le moment (pcfm) et
31!         les flux de chaleur sensible et latente (pcfh) d'apr??s
32!         Louis 1982, Louis 1979, King et al 2001
33!         ou Zilitinkevich et al 2002  pour les cas stables, Louis 1979
34!         et 1982 pour les cas instables
35!
36! Modified history:
37!  writting on the 20/05/2016
38!  modified on the 13/12/2016 to be adapted to LMDZ6
39!
40! References:
41!   Louis, J. F., 1979: A parametric model of vertical eddy fluxes in the
42!     atmosphere. Boundary-Layer Meteorology. 01/1979; 17(2):187-202.
43!   Louis, J. F., Tiedtke, M. and Geleyn, J. F., 1982: `A short history of the
44!     operational PBL parametrization at ECMWF'. Workshop on boundary layer
45!     parametrization, November 1981, ECMWF, Reading, England.
46!     Page: 19. Equations in Table 1.
47!   Mascart P, Noilhan J, Giordani H 1995.A MODIFIED PARAMETERIZATION OF FLUX-PROFILE RELATIONSHIPS
48!    IN THE SURFACE LAYER USING DIFFERENT ROUGHNESS LENGTH VALUES FOR HEAT AND MOMENTUM
49!    Boundary-Layer Meteorology 72: 331-344
50!   Anton Beljaars. May 1992. The parametrization of the planetary boundary layer. 
51!     European Centre for Medium-Range Weather Forecasts.
52!     Equations: 110-113. Page 40.
53!   Miller,M.J., A.C.M.Beljaars, T.N.Palmer. 1992. The sensitivity of the ECMWF
54!     model to the parameterization of evaporation from the tropical oceans. J.
55!     Climate, 5:418-434.
56!   King J.C, Connolley, W.M ad Derbyshire S.H. 2001, Sensitivity of Modelled Antarctic climate
57!   to surface and boundary-layer flux parametrizations
58!   QJRMS, 127, pp 779-794
59!
60! ================================================================= c
61! ================================================================= c
62! On choisit le couple de fonctions de correction avec deux flags:
63! Un pour les cas instables, un autre pour les cas stables
64!
65! iflag_corr_insta:
66!                   1: Louis 1979 avec les modifications de Mascart 1995 (z0/= z0h)
67!                   2: Louis 1982
68!                   3: Laurent Li
69!
70! iflag_corr_sta:
71!                   1: Louis 1979 avec les modifications de Mascart 1995 (z0/= z0h)
72!                   2: Louis 1982
73!                   3: Laurent Li
74!                   4: King  2001 (SHARP)
75!                   5: MO 1st order theory (allow collapse of turbulence)
76!           
77!
78!*****************************************************************
79! Parametres d'entree
80!*****************************************************************
81 
82  INTEGER, INTENT(IN)                      :: knon, nsrf ! nombre de points de grille sur l'horizontal + type de surface
83  REAL, DIMENSION(klon), INTENT(IN)        :: speed ! module du vent au 1er niveau du modele
84  REAL, DIMENSION(klon), INTENT(IN)        :: zgeop1! geopotentiel au 1er niveau du modele
85  REAL, DIMENSION(klon), INTENT(IN)        :: tsurf ! Surface temperature (K)
86  REAL, DIMENSION(klon), INTENT(IN)        :: qsurf ! Surface humidity (Kg/Kg)
87  REAL, DIMENSION(klon), INTENT(IN)        :: z0m, z0h ! Rugosity at surface (m)
88  REAL, DIMENSION(klon), INTENT(IN)        :: ri_in ! Input Richardson 1st layer for first guess calculations of screen var.
89  INTEGER, INTENT(IN)                      :: iri_in! iflag to activate cdrag calculation using ri1
90  REAL, DIMENSION(klon), INTENT(IN)        :: t1  ! Temperature au premier niveau (K)
91  REAL, DIMENSION(klon), INTENT(IN)        :: q1  ! humidite specifique au premier niveau (kg/kg)
92  REAL, DIMENSION(klon), INTENT(IN)        :: psol ! pression au sol
93
94
95
96! Parametres de sortie
97!******************************************************************
98  REAL, DIMENSION(klon), INTENT(OUT)       :: cdm   ! Drag coefficient for momentum
99  REAL, DIMENSION(klon), INTENT(OUT)       :: cdh   ! Drag coefficient for heat flux
100  REAL, DIMENSION(klon), INTENT(OUT)       :: zri   ! Richardson number
101  REAL, DIMENSION(klon), INTENT(OUT)       :: pref  ! Pression au niveau zgeop/RG
102
103! Variables Locales
104!******************************************************************
105 
106
107  INCLUDE "YOMCST.h"
108  INCLUDE "YOETHF.h"
109  INCLUDE "clesphys.h"
110
111
112  REAL, PARAMETER       :: CKAP=0.40, CKAPT=0.42
113  REAL                   CEPDU2
114  REAL                   ALPHA
115  REAL                   CB,CC,CD,C2,C3
116  REAL                   MU, CM, CH, B, CMstar, CHstar
117  REAL                   PM, PH, BPRIME
118  REAL                   C
119  INTEGER                ng_q1                      ! Number of grids that q1 < 0.0
120  INTEGER                ng_qsurf                   ! Number of grids that qsurf < 0.0
121  INTEGER                i
122  REAL                   zdu2, ztsolv
123  REAL                   ztvd, zscf
124  REAL                   zucf, zcr
125  REAL                   friv, frih
126  REAL, DIMENSION(klon) :: FM, FH                   ! stability functions
127  REAL, DIMENSION(klon) :: cdmn, cdhn               ! Drag coefficient in neutral conditions
128  REAL zzzcd
129  REAL, DIMENSION(klon) :: sm, prandtl              ! Stability function from atke turbulence scheme
130  REAL                  ri0, ri1, cn                ! to have dimensionless term in sm and prandtl
131
132  LOGICAL, SAVE :: firstcall = .TRUE.
133  !$OMP THREADPRIVATE(firstcall)
134  INTEGER, SAVE :: iflag_corr_sta
135  !$OMP THREADPRIVATE(iflag_corr_sta)
136  INTEGER, SAVE :: iflag_corr_insta
137  !$OMP THREADPRIVATE(iflag_corr_insta)
138
139!===================================================================c
140! Valeurs numeriques des constantes
141!===================================================================c
142
143
144! Minimum du carre du vent
145
146 CEPDU2 = (0.1)**2
147
148! Louis 1982
149
150  CB=5.0
151  CC=5.0
152  CD=5.0
153
154
155! King 2001
156
157  C2=0.25
158  C3=0.0625
159
160 
161! Louis 1979
162
163  BPRIME=4.7
164  B=9.4
165 
166
167!MO
168
169  ALPHA=5.0
170
171! Consistent with atke scheme
172cn=(1./sqrt(cepsilon))**(2/3)
173ri0=2./rpi*(cinf - cn)*ric/cn
174ri1=-2./rpi * (pr_asym - pr_neut)/pr_slope
175
176
177! ================================================================= c
178! Tests avant de commencer
179! Fuxing WANG, 04/03/2015
180! To check if there are negative q1, qsurf values.
181!====================================================================c
182  ng_q1 = 0     ! Initialization
183  ng_qsurf = 0  ! Initialization
184  DO i = 1, knon
185     IF (q1(i).LT.0.0)     ng_q1 = ng_q1 + 1
186     IF (qsurf(i).LT.0.0)  ng_qsurf = ng_qsurf + 1
187  ENDDO
188  IF (ng_q1.GT.0 .and. prt_level > 5) THEN
189      WRITE(lunout,*)" *** Warning: Negative q1(humidity at 1st level) values in cdrag.F90 !"
190      WRITE(lunout,*)" The total number of the grids is: ", ng_q1
191      WRITE(lunout,*)" The negative q1 is set to zero "
192!      abort_message="voir ci-dessus"
193!      CALL abort_physic(modname,abort_message,1)
194  ENDIF
195  IF (ng_qsurf.GT.0 .and. prt_level > 5) THEN
196      WRITE(lunout,*)" *** Warning: Negative qsurf(humidity at surface) values in cdrag.F90 !"
197      WRITE(lunout,*)" The total number of the grids is: ", ng_qsurf
198      WRITE(lunout,*)" The negative qsurf is set to zero "
199!      abort_message="voir ci-dessus"
200!      CALL abort_physic(modname,abort_message,1)
201  ENDIF
202
203
204
205!=============================================================================c
206! Calcul du cdrag
207!=============================================================================c
208
209! On choisit les fonctions de stabilite utilisees au premier appel
210!**************************************************************************
211  IF (firstcall) THEN
212   iflag_corr_sta=2
213   iflag_corr_insta=2
214 
215   CALL getin_p('iflag_corr_sta',iflag_corr_sta)
216   CALL getin_p('iflag_corr_insta',iflag_corr_insta)
217
218   firstcall = .FALSE.
219 ENDIF
220
221!xxxxxxxxxxxxxxxxxxxxxxx
222  DO i = 1, knon        ! Boucle sur l'horizontal
223!xxxxxxxxxxxxxxxxxxxxxxx
224
225
226! calculs preliminaires:
227!***********************
228     
229
230     zdu2 = MAX(CEPDU2, speed(i)**2)
231     pref(i) = EXP(LOG(psol(i)) - zgeop1(i)/(RD*t1(i)* &
232                 (1.+ RETV * max(q1(i),0.0))))           ! negative q1 set to zero
233     ztsolv = tsurf(i) * (1.0+RETV*max(qsurf(i),0.0))    ! negative qsurf set to zero
234     ztvd = (t1(i)+zgeop1(i)/RCPD/(1.+RVTMP2*max(q1(i),0.0))) &
235          *(1.+RETV*max(q1(i),0.0))                      ! negative q1 set to zero
236     zri(i) = zgeop1(i)*(ztvd-ztsolv)/(zdu2*ztvd)
237     IF (iri_in.EQ.1) THEN
238       zri(i) = ri_in(i)
239     ENDIF
240
241
242! Coefficients CD neutres : k^2/ln(z/z0) et k^2/(ln(z/z0)*ln(z/z0h)):
243!********************************************************************
244
245     zzzcd=CKAP/LOG(1.+zgeop1(i)/(RG*z0m(i)))
246     cdmn(i) = zzzcd*zzzcd
247     cdhn(i) = zzzcd*(CKAP/LOG(1.+zgeop1(i)/(RG*z0h(i))))
248
249
250! Calcul des fonctions de stabilit?? FMs, FHs, FMi, FHi :
251!*******************************************************
252
253
254
255!''''''''''''''
256! Cas instables
257!''''''''''''''
258
259 IF (zri(i) .LT. 0.) THEN   
260
261
262        SELECT CASE (iflag_corr_insta)
263   
264        CASE (1) ! Louis 1979 + Mascart 1995
265
266           MU=LOG(MAX(z0m(i)/z0h(i),0.01))
267           CMstar=6.8741+2.6933*MU-0.3601*(MU**2)+0.0154*(MU**3)
268           PM=0.5233-0.0815*MU+0.0135*(MU**2)-0.001*(MU**3)
269           CHstar=3.2165+4.3431*MU+0.536*(MU**2)-0.0781*(MU**3)
270           PH=0.5802-0.1571*MU+0.0327*(MU**2)-0.0026*(MU**3)
271           CH=CHstar*B*CKAP/LOG(z0m(i)+zgeop1(i)/(RG*z0m(i))) &
272            & * CKAPT/LOG(z0h(i)+zgeop1(i)/(RG*z0h(i)))       &
273            & * ((zgeop1(i)/(RG*z0h(i)))**PH)
274           CM=CMstar*B*CKAP/LOG(z0m(i)+zgeop1(i)/(RG*z0m(i))) &
275            & *CKAP/LOG(z0m(i)+zgeop1(i)/(RG*z0m(i)))         &
276            & * ((zgeop1(i)/(RG*z0m(i)))**PM)
277         
278
279
280
281           FM(i)=1.-B*zri(i)/(1.+CM*SQRT(ABS(zri(i))))
282           FH(i)=1.-B*zri(i)/(1.+CH*SQRT(ABS(zri(i))))
283
284        CASE (2) ! Louis 1982
285
286           zucf = 1./(1.+3.0*CB*CC*cdmn(i)*SQRT(ABS(zri(i)) &
287                *(1.0+zgeop1(i)/(RG*z0m(i)))))
288           FM(i) = AMAX1((1.-2.0*CB*zri(i)*zucf),f_ri_cd_min)
289           FH(i) = AMAX1((1.-3.0*CB*zri(i)*zucf),f_ri_cd_min)
290
291
292        CASE (3) ! Laurent Li
293
294           
295           FM(i) = MAX(SQRT(1.0-18.0*zri(i)),f_ri_cd_min)
296           FH(i) = MAX(SQRT(1.0-18.0*zri(i)),f_ri_cd_min)
297
298         CASE (6) ! Consistent with turbulence scheme (in stationary case) derived in atke (2023)
299         
300           sm(i) = 2./rpi * (cinf - cn) * atan(-zri(i)/ri0) + cn
301           prandtl(i) = -2./rpi * (pr_asym - pr_neut) * atan(zri(i)/ri1) + pr_neut
302           FM(i) = MAX((sm(i)**(3/2) * sqrt(cepsilon) * (1 - zri(i) / prandtl(i))**(1/2)),f_ri_cd_min)
303           FH(i) = MAX((FM(i) / prandtl(i)), f_ri_cd_min)
304
305         CASE default ! Louis 1982
306
307           zucf = 1./(1.+3.0*CB*CC*cdmn(i)*SQRT(ABS(zri(i)) &
308                *(1.0+zgeop1(i)/(RG*z0m(i)))))
309           FM(i) = AMAX1((1.-2.0*CB*zri(i)*zucf),f_ri_cd_min)
310           FH(i) = AMAX1((1.-3.0*CB*zri(i)*zucf),f_ri_cd_min)
311
312
313         END SELECT
314
315
316
317! Calcul des drags
318
319
320       cdm(i)=cdmn(i)*FM(i)
321       cdh(i)=f_cdrag_ter*cdhn(i)*FH(i)
322
323
324! Traitement particulier des cas oceaniques
325! on applique Miller et al 1992 en l'absence de gustiness
326
327  IF (nsrf == is_oce) THEN
328!        cdh(i)=f_cdrag_oce*cdhn(i)*FH(i)
329
330        IF(iflag_gusts==0) THEN
331           zcr = (0.0016/(cdmn(i)*SQRT(zdu2)))*ABS(ztvd-ztsolv)**(1./3.)
332           cdh(i) =f_cdrag_oce* cdhn(i)*(1.0+zcr**1.25)**(1./1.25)
333        ENDIF
334
335
336        cdm(i)=MIN(cdm(i),cdmmax)
337        cdh(i)=MIN(cdh(i),cdhmax)
338
339  END IF
340
341
342
343 ELSE                         
344
345!'''''''''''''''
346! Cas stables :
347!'''''''''''''''
348        zri(i) = MIN(20.,zri(i))
349
350       SELECT CASE (iflag_corr_sta)
351   
352        CASE (1) ! Louis 1979 + Mascart 1995
353   
354           FM(i)=MAX(1./((1+BPRIME*zri(i))**2),f_ri_cd_min)
355           FH(i)=FM(i)
356
357
358        CASE (2) ! Louis 1982
359           
360           zscf = SQRT(1.+CD*ABS(zri(i)))
361           FM(i)= AMAX1(1. / (1.+2.*CB*zri(i)/zscf), f_ri_cd_min)
362           FH(i)= AMAX1(1./ (1.+3.*CB*zri(i)*zscf), f_ri_cd_min )
363
364
365        CASE (3) ! Laurent Li
366   
367           FM(i)=MAX(1.0 / (1.0+10.0*zri(i)*(1+8.0*zri(i))),f_ri_cd_min)
368           FH(i)=FM(i)
369
370
371        CASE (4)  ! King 2001
372         
373           if (zri(i) .LT. C2/2.) then
374           FM(i)=MAX((1.-zri(i)/C2)**2,f_ri_cd_min)
375           FH(i)=  FM(i)
376
377
378           else
379           FM(i)=MAX(C3*((C2/zri(i))**2),f_ri_cd_min)
380           FH(i)= FM(i)
381           endif
382
383
384        CASE (5) ! MO
385   
386          if (zri(i) .LT. 1./alpha) then
387
388           FM(i)=MAX((1.-alpha*zri(i))**2,f_ri_cd_min)
389           FH(i)=FM(i)
390
391           else
392
393
394           FM(i)=MAX(1E-7,f_ri_cd_min)
395           FH(i)=FM(i)
396
397          endif
398
399        CASE (6) ! Consistent with turbulence scheme (in stationary case) derived in atke (2023)
400          sm(i) = MAX(0.,cn*(1.-zri(i)/ric))
401          prandtl(i) = pr_neut + zri(i) * pr_slope
402          FM(i) = MAX((sm(i)**(3/2) * sqrt(cepsilon) * (1 - zri(i) / prandtl(i))**(1/2)),f_ri_cd_min)
403          FH(i) = MAX((FM(i) / prandtl(i)), f_ri_cd_min)
404
405        CASE default ! Louis 1982
406
407           zscf = SQRT(1.+CD*ABS(zri(i)))
408           FM(i)= AMAX1(1. / (1.+2.*CB*zri(i)/zscf), f_ri_cd_min)
409           FH(i)= AMAX1(1./ (1.+3.*CB*zri(i)*zscf), f_ri_cd_min )
410
411
412
413   END SELECT
414
415! Calcul des drags
416
417
418       cdm(i)=cdmn(i)*FM(i)
419       cdh(i)=f_cdrag_ter*cdhn(i)*FH(i)
420
421       IF(nsrf.EQ.is_oce) THEN
422
423        cdh(i)=f_cdrag_oce*cdhn(i)*FH(i)
424        cdm(i)=MIN(cdm(i),cdmmax)
425        cdh(i)=MIN(cdh(i),cdhmax)
426
427      ENDIF
428
429
430 ENDIF
431
432!xxxxxxxxxxx
433  END DO   !  Fin de la boucle sur l'horizontal
434!xxxxxxxxxxx
435! ================================================================= c
436     
437END SUBROUTINE cdrag
438
439END MODULE cdrag_mod
Note: See TracBrowser for help on using the repository browser.