source: LMDZ6/branches/Amaury_dev/libf/phylmd/conccm.F90 @ 5218

Last change on this file since 5218 was 5153, checked in by abarral, 5 months ago

Revert FCTTRE to INCLUDE to assess impact of inlining

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 28.7 KB
Line 
1! $Header$
2
3SUBROUTINE conccm(dtime, paprs, pplay, t, q, conv_q, d_t, d_q, rain, snow, &
4        kbascm, ktopcm)
5
6  USE dimphy
7  USE lmdz_yoethf
8  USE lmdz_yomcst
9
10  IMPLICIT NONE
11  ! ======================================================================
12  ! Auteur(s): Z.X. Li (LMD/CNRS) date: le 14 mars 1996
13  ! Objet: Schema simple (avec flux de masse) pour la convection
14  ! (schema standard du modele NCAR CCM2)
15  ! ======================================================================
16
17  ! Entree:
18  REAL dtime ! pas d'integration
19  REAL paprs(klon, klev + 1) ! pression inter-couche (Pa)
20  REAL pplay(klon, klev) ! pression au milieu de couche (Pa)
21  REAL t(klon, klev) ! temperature (K)
22  REAL q(klon, klev) ! humidite specifique (g/g)
23  REAL conv_q(klon, klev) ! taux de convergence humidite (g/g/s)
24  ! Sortie:
25  REAL d_t(klon, klev) ! incrementation temperature
26  REAL d_q(klon, klev) ! incrementation vapeur
27  REAL rain(klon) ! pluie (mm/s)
28  REAL snow(klon) ! neige (mm/s)
29  INTEGER kbascm(klon) ! niveau du bas de convection
30  INTEGER ktopcm(klon) ! niveau du haut de convection
31
32  REAL pt(klon, klev)
33  REAL pq(klon, klev)
34  REAL pres(klon, klev)
35  REAL dp(klon, klev)
36  REAL zgeom(klon, klev)
37  REAL cmfprs(klon)
38  REAL cmfprt(klon)
39  INTEGER ntop(klon)
40  INTEGER nbas(klon)
41  INTEGER i, k
42  REAL zlvdcp, zlsdcp, zdelta, zz, za, zb
43
44  LOGICAL usekuo ! utiliser convection profonde (schema Kuo)
45  PARAMETER (usekuo = .TRUE.)
46
47  REAL d_t_bis(klon, klev)
48  REAL d_q_bis(klon, klev)
49  REAL rain_bis(klon)
50  REAL snow_bis(klon)
51  INTEGER ibas_bis(klon)
52  INTEGER itop_bis(klon)
53  REAL d_ql_bis(klon, klev)
54  REAL rneb_bis(klon, klev)
55
56  ! initialiser les variables de sortie (pour securite)
57  DO i = 1, klon
58    rain(i) = 0.0
59    snow(i) = 0.0
60    kbascm(i) = 0
61    ktopcm(i) = 0
62  END DO
63  DO k = 1, klev
64    DO i = 1, klon
65      d_t(i, k) = 0.0
66      d_q(i, k) = 0.0
67    END DO
68  END DO
69
70  ! preparer les variables d'entree (attention: l'ordre des niveaux
71  ! verticaux augmente du haut vers le bas)
72  DO k = 1, klev
73    DO i = 1, klon
74      pt(i, k) = t(i, klev - k + 1)
75      pq(i, k) = q(i, klev - k + 1)
76      pres(i, k) = pplay(i, klev - k + 1)
77      dp(i, k) = paprs(i, klev + 1 - k) - paprs(i, klev + 1 - k + 1)
78    END DO
79  END DO
80  DO i = 1, klon
81    zgeom(i, klev) = rd * t(i, 1) / (0.5 * (paprs(i, 1) + pplay(i, &
82            1))) * (paprs(i, 1) - pplay(i, 1))
83  END DO
84  DO k = 2, klev
85    DO i = 1, klon
86      zgeom(i, klev + 1 - k) = zgeom(i, klev + 1 - k + 1) + rd * 0.5 * (t(i, k - 1) + t(i, k)) / &
87              paprs(i, k) * (pplay(i, k - 1) - pplay(i, k))
88    END DO
89  END DO
90
91  CALL cmfmca(dtime, pres, dp, zgeom, pt, pq, cmfprt, cmfprs, ntop, nbas)
92
93  DO k = 1, klev
94    DO i = 1, klon
95      d_q(i, klev + 1 - k) = pq(i, k) - q(i, klev + 1 - k)
96      d_t(i, klev + 1 - k) = pt(i, k) - t(i, klev + 1 - k)
97    END DO
98  END DO
99
100  DO i = 1, klon
101    rain(i) = cmfprt(i) * rhoh2o
102    snow(i) = cmfprs(i) * rhoh2o
103    kbascm(i) = klev + 1 - nbas(i)
104    ktopcm(i) = klev + 1 - ntop(i)
105  END DO
106
107  IF (usekuo) THEN
108    CALL conkuo(dtime, paprs, pplay, t, q, conv_q, d_t_bis, d_q_bis, &
109            d_ql_bis, rneb_bis, rain_bis, snow_bis, ibas_bis, itop_bis)
110    DO k = 1, klev
111      DO i = 1, klon
112        d_t(i, k) = d_t(i, k) + d_t_bis(i, k)
113        d_q(i, k) = d_q(i, k) + d_q_bis(i, k)
114      END DO
115    END DO
116    DO i = 1, klon
117      rain(i) = rain(i) + rain_bis(i)
118      snow(i) = snow(i) + snow_bis(i)
119      kbascm(i) = min(kbascm(i), ibas_bis(i))
120      ktopcm(i) = max(ktopcm(i), itop_bis(i))
121    END DO
122    DO k = 1, klev ! eau liquide convective est
123      DO i = 1, klon ! dispersee dans l'air
124        zlvdcp = rlvtt / rcpd / (1.0 + rvtmp2 * q(i, k))
125        zlsdcp = rlstt / rcpd / (1.0 + rvtmp2 * q(i, k))
126        zdelta = max(0., sign(1., rtt - t(i, k)))
127        zz = d_ql_bis(i, k) ! re-evap. de l'eau liquide
128        zb = max(0.0, zz)
129        za = -max(0.0, zz) * (zlvdcp * (1. - zdelta) + zlsdcp * zdelta)
130        d_t(i, k) = d_t(i, k) + za
131        d_q(i, k) = d_q(i, k) + zb
132      END DO
133    END DO
134  END IF
135
136END SUBROUTINE conccm
137SUBROUTINE cmfmca(deltat, p, dp, gz, tb, shb, cmfprt, cmfprs, cnt, cnb)
138  USE dimphy
139  USE lmdz_yoethf
140
141  USE lmdz_yomcst
142
143  IMPLICIT NONE
144 INCLUDE "FCTTRE.h"
145  ! -----------------------------------------------------------------------
146  ! Moist convective mass flux procedure:
147  ! If stratification is unstable to nonentraining parcel ascent,
148  ! complete an adjustment making use of a simple cloud model
149
150  ! Code generalized to allow specification of parcel ("updraft")
151  ! properties, as well as convective transport of an arbitrary
152  ! number of passive constituents (see cmrb array).
153  ! ----------------------------Code History-------------------------------
154  ! Original version:  J. J. Hack, March 22, 1990
155  ! Standardized:      J. Rosinski, June 1992
156  ! Reviewed:          J. Hack, G. Taylor, August 1992
157  ! Adaptation au LMD: Z.X. Li, mars 1996 (reference: Hack 1994,
158  ! J. Geophys. Res. vol 99, D3, 5551-5568). J'ai
159  ! introduit les constantes et les fonctions thermo-
160  ! dynamiques du Centre Europeen. J'ai elimine le
161  ! re-indicage du code en esperant que cela pourra
162  ! simplifier la lecture et la comprehension.
163  ! -----------------------------------------------------------------------
164  INTEGER pcnst ! nombre de traceurs passifs
165  PARAMETER (pcnst = 1)
166  ! ------------------------------Arguments--------------------------------
167  ! Input arguments
168
169  REAL deltat ! time step (seconds)
170  REAL p(klon, klev) ! pressure
171  REAL dp(klon, klev) ! delta-p
172  REAL gz(klon, klev) ! geopotential (a partir du sol)
173
174  REAL thtap(klon) ! PBL perturbation theta
175  REAL shp(klon) ! PBL perturbation specific humidity
176  REAL pblh(klon) ! PBL height (provided by PBL routine)
177  REAL cmrp(klon, pcnst) ! constituent perturbations in PBL
178
179  ! Updated arguments:
180
181  REAL tb(klon, klev) ! temperature (t bar)
182  REAL shb(klon, klev) ! specific humidity (sh bar)
183  REAL cmrb(klon, klev, pcnst) ! constituent mixing ratios (cmr bar)
184
185  ! Output arguments
186
187  REAL cmfdt(klon, klev) ! dT/dt due to moist convection
188  REAL cmfdq(klon, klev) ! dq/dt due to moist convection
189  REAL cmfmc(klon, klev) ! moist convection cloud mass flux
190  REAL cmfdqr(klon, klev) ! dq/dt due to convective rainout
191  REAL cmfsl(klon, klev) ! convective lw static energy flux
192  REAL cmflq(klon, klev) ! convective total water flux
193  REAL cmfprt(klon) ! convective precipitation rate
194  REAL cmfprs(klon) ! convective snowfall rate
195  REAL qc(klon, klev) ! dq/dt due to rainout terms
196  INTEGER cnt(klon) ! top level of convective activity
197  INTEGER cnb(klon) ! bottom level of convective activity
198  ! ------------------------------Parameters-------------------------------
199  REAL c0 ! rain water autoconversion coefficient
200  PARAMETER (c0 = 1.0E-4)
201  REAL dzmin ! minimum convective depth for precipitation
202  PARAMETER (dzmin = 0.0)
203  REAL betamn ! minimum overshoot parameter
204  PARAMETER (betamn = 0.10)
205  REAL cmftau ! characteristic adjustment time scale
206  PARAMETER (cmftau = 3600.)
207  INTEGER limcnv ! top interface level limit for convection
208  PARAMETER (limcnv = 1)
209  REAL tpmax ! maximum acceptable t perturbation (degrees C)
210  PARAMETER (tpmax = 1.50)
211  REAL shpmax ! maximum acceptable q perturbation (g/g)
212  PARAMETER (shpmax = 1.50E-3)
213  REAL tiny ! arbitrary small num used in transport estimates
214  PARAMETER (tiny = 1.0E-36)
215  REAL eps ! convergence criteria (machine dependent)
216  PARAMETER (eps = 1.0E-13)
217  REAL tmelt ! freezing point of water(req'd for rain vs snow)
218  PARAMETER (tmelt = 273.15)
219  REAL ssfac ! supersaturation bound (detrained air)
220  PARAMETER (ssfac = 1.001)
221
222  ! ---------------------------Local workspace-----------------------------
223  REAL gam(klon, klev) ! L/cp (d(qsat)/dT)
224  REAL sb(klon, klev) ! dry static energy (s bar)
225  REAL hb(klon, klev) ! moist static energy (h bar)
226  REAL shbs(klon, klev) ! sat. specific humidity (sh bar star)
227  REAL hbs(klon, klev) ! sat. moist static energy (h bar star)
228  REAL shbh(klon, klev + 1) ! specific humidity on interfaces
229  REAL sbh(klon, klev + 1) ! s bar on interfaces
230  REAL hbh(klon, klev + 1) ! h bar on interfaces
231  REAL cmrh(klon, klev + 1) ! interface constituent mixing ratio
232  REAL prec(klon) ! instantaneous total precipitation
233  REAL dzcld(klon) ! depth of convective layer (m)
234  REAL beta(klon) ! overshoot parameter (fraction)
235  REAL betamx ! local maximum on overshoot
236  REAL eta(klon) ! convective mass flux (kg/m^2 s)
237  REAL etagdt ! eta*grav*deltat
238  REAL cldwtr(klon) ! cloud water (mass)
239  REAL rnwtr(klon) ! rain water  (mass)
240  REAL sc(klon) ! dry static energy   ("in-cloud")
241  REAL shc(klon) ! specific humidity   ("in-cloud")
242  REAL hc(klon) ! moist static energy ("in-cloud")
243  REAL cmrc(klon) ! constituent mix rat ("in-cloud")
244  REAL dq1(klon) ! shb  convective change (lower lvl)
245  REAL dq2(klon) ! shb  convective change (mid level)
246  REAL dq3(klon) ! shb  convective change (upper lvl)
247  REAL ds1(klon) ! sb   convective change (lower lvl)
248  REAL ds2(klon) ! sb   convective change (mid level)
249  REAL ds3(klon) ! sb   convective change (upper lvl)
250  REAL dcmr1(klon) ! cmrb convective change (lower lvl)
251  REAL dcmr2(klon) ! cmrb convective change (mid level)
252  REAL dcmr3(klon) ! cmrb convective change (upper lvl)
253  REAL flotab(klon) ! hc - hbs (mesure d'instabilite)
254  LOGICAL ldcum(klon) ! .TRUE. si la convection existe
255  LOGICAL etagt0 ! true if eta > 0.0
256  REAL dt ! current 2 delta-t (model time step)
257  REAL cats ! modified characteristic adj. time
258  REAL rdt ! 1./dt
259  REAL qprime ! modified specific humidity pert.
260  REAL tprime ! modified thermal perturbation
261  REAL pblhgt ! bounded pbl height (max[pblh,1m])
262  REAL fac1 ! intermediate scratch variable
263  REAL shprme ! intermediate specific humidity pert.
264  REAL qsattp ! saturation mixing ratio for
265  !  thermally perturbed PBL parcels
266  REAL dz ! local layer depth
267  REAL b1 ! bouyancy measure in detrainment lvl
268  REAL b2 ! bouyancy measure in condensation lvl
269  REAL g ! bounded vertical gradient of hb
270  REAL tmass ! total mass available for convective exchange
271  REAL denom ! intermediate scratch variable
272  REAL qtest1 ! used in negative q test (middle lvl)
273  REAL qtest2 ! used in negative q test (lower lvl)
274  REAL fslkp ! flux lw static energy (bot interface)
275  REAL fslkm ! flux lw static energy (top interface)
276  REAL fqlkp ! flux total water (bottom interface)
277  REAL fqlkm ! flux total water (top interface)
278  REAL botflx ! bottom constituent mixing ratio flux
279  REAL topflx ! top constituent mixing ratio flux
280  REAL efac1 ! ratio cmrb to convectively induced change (bot lvl)
281  REAL efac2 ! ratio cmrb to convectively induced change (mid lvl)
282  REAL efac3 ! ratio cmrb to convectively induced change (top lvl)
283
284  INTEGER i, k ! indices horizontal et vertical
285  INTEGER km1 ! k-1 (index offset)
286  INTEGER kp1 ! k+1 (index offset)
287  INTEGER m ! constituent index
288  INTEGER ktp ! temporary index used to track top
289  INTEGER is ! nombre de points a ajuster
290
291  REAL tmp1, tmp2, tmp3, tmp4
292  REAL zx_t, zx_p, zx_q, zx_qs, zx_gam
293  REAL zcor, zdelta, zcvm5
294
295  REAL qhalf, sh1, sh2, shbs1, shbs2
296  qhalf(sh1, sh2, shbs1, shbs2) = min(max(sh1, sh2), &
297          (shbs2 * sh1 + shbs1 * sh2) / (shbs1 + shbs2))
298
299  ! -----------------------------------------------------------------------
300  ! pas de traceur pour l'instant
301  DO m = 1, pcnst
302    DO k = 1, klev
303      DO i = 1, klon
304        cmrb(i, k, m) = 0.0
305      END DO
306    END DO
307  END DO
308
309  ! Les perturbations de la couche limite sont zero pour l'instant
310
311  DO m = 1, pcnst
312    DO i = 1, klon
313      cmrp(i, m) = 0.0
314    END DO
315  END DO
316  DO i = 1, klon
317    thtap(i) = 0.0
318    shp(i) = 0.0
319    pblh(i) = 1.0
320  END DO
321
322  ! Ensure that characteristic adjustment time scale (cmftau) assumed
323  ! in estimate of eta isn't smaller than model time scale (deltat)
324
325  dt = deltat
326  cats = max(dt, cmftau)
327  rdt = 1.0 / dt
328
329  ! Compute sb,hb,shbs,hbs
330
331  DO k = 1, klev
332    DO i = 1, klon
333      zx_t = tb(i, k)
334      zx_p = p(i, k)
335      zx_q = shb(i, k)
336      zdelta = max(0., sign(1., rtt - zx_t))
337      zcvm5 = r5les * rlvtt * (1. - zdelta) + r5ies * rlstt * zdelta
338      zcvm5 = zcvm5 / rcpd / (1.0 + rvtmp2 * zx_q)
339      zx_qs = r2es * foeew(zx_t, zdelta) / zx_p
340      zx_qs = min(0.5, zx_qs)
341      zcor = 1. / (1. - retv * zx_qs)
342      zx_qs = zx_qs * zcor
343      zx_gam = foede(zx_t, zdelta, zcvm5, zx_qs, zcor)
344      shbs(i, k) = zx_qs
345      gam(i, k) = zx_gam
346    END DO
347  END DO
348
349  DO k = limcnv, klev
350    DO i = 1, klon
351      sb(i, k) = rcpd * tb(i, k) + gz(i, k)
352      hb(i, k) = sb(i, k) + rlvtt * shb(i, k)
353      hbs(i, k) = sb(i, k) + rlvtt * shbs(i, k)
354    END DO
355  END DO
356
357  ! Compute sbh, shbh
358
359  DO k = limcnv + 1, klev
360    km1 = k - 1
361    DO i = 1, klon
362      sbh(i, k) = 0.5 * (sb(i, km1) + sb(i, k))
363      shbh(i, k) = qhalf(shb(i, km1), shb(i, k), shbs(i, km1), shbs(i, k))
364      hbh(i, k) = sbh(i, k) + rlvtt * shbh(i, k)
365    END DO
366  END DO
367
368  ! Specify properties at top of model (not used, but filling anyway)
369
370  DO i = 1, klon
371    sbh(i, limcnv) = sb(i, limcnv)
372    shbh(i, limcnv) = shb(i, limcnv)
373    hbh(i, limcnv) = hb(i, limcnv)
374  END DO
375
376  ! Zero vertically independent control, tendency & diagnostic arrays
377
378  DO i = 1, klon
379    prec(i) = 0.0
380    dzcld(i) = 0.0
381    cnb(i) = 0
382    cnt(i) = klev + 1
383  END DO
384
385  DO k = 1, klev
386    DO i = 1, klon
387      cmfdt(i, k) = 0.
388      cmfdq(i, k) = 0.
389      cmfdqr(i, k) = 0.
390      cmfmc(i, k) = 0.
391      cmfsl(i, k) = 0.
392      cmflq(i, k) = 0.
393    END DO
394  END DO
395
396  ! Begin moist convective mass flux adjustment procedure.
397  ! Formalism ensures that negative cloud liquid water can never occur
398
399  DO k = klev - 1, limcnv + 1, -1
400    km1 = k - 1
401    kp1 = k + 1
402    DO i = 1, klon
403      eta(i) = 0.0
404      beta(i) = 0.0
405      ds1(i) = 0.0
406      ds2(i) = 0.0
407      ds3(i) = 0.0
408      dq1(i) = 0.0
409      dq2(i) = 0.0
410      dq3(i) = 0.0
411
412      ! Specification of "cloud base" conditions
413
414      qprime = 0.0
415      tprime = 0.0
416
417      ! Assign tprime within the PBL to be proportional to the quantity
418      ! thtap (which will be bounded by tpmax), passed to this routine by
419      ! the PBL routine.  Don't allow perturbation to produce a dry
420      ! adiabatically unstable parcel.  Assign qprime within the PBL to be
421      ! an appropriately modified value of the quantity shp (which will be
422      ! bounded by shpmax) passed to this routine by the PBL routine.  The
423      ! quantity qprime should be less than the local saturation value
424      ! (qsattp=qsat[t+tprime,p]).  In both cases, thtap and shp are
425      ! linearly reduced toward zero as the PBL top is approached.
426
427      pblhgt = max(pblh(i), 1.0)
428      IF (gz(i, kp1) / rg<=pblhgt .AND. dzcld(i)==0.0) THEN
429        fac1 = max(0.0, 1.0 - gz(i, kp1) / rg / pblhgt)
430        tprime = min(thtap(i), tpmax) * fac1
431        qsattp = shbs(i, kp1) + rcpd / rlvtt * gam(i, kp1) * tprime
432        shprme = min(min(shp(i), shpmax) * fac1, max(qsattp - shb(i, kp1), 0.0))
433        qprime = max(qprime, shprme)
434      ELSE
435        tprime = 0.0
436        qprime = 0.0
437      END IF
438
439      ! Specify "updraft" (in-cloud) thermodynamic properties
440
441      sc(i) = sb(i, kp1) + rcpd * tprime
442      shc(i) = shb(i, kp1) + qprime
443      hc(i) = sc(i) + rlvtt * shc(i)
444      flotab(i) = hc(i) - hbs(i, k)
445      dz = dp(i, k) * rd * tb(i, k) / rg / p(i, k)
446      IF (flotab(i)>0.0) THEN
447        dzcld(i) = dzcld(i) + dz
448      ELSE
449        dzcld(i) = 0.0
450      END IF
451    END DO
452
453    ! Check on moist convective instability
454
455    is = 0
456    DO i = 1, klon
457      IF (flotab(i)>0.0) THEN
458        ldcum(i) = .TRUE.
459        is = is + 1
460      ELSE
461        ldcum(i) = .FALSE.
462      END IF
463    END DO
464
465    IF (is==0) THEN
466      DO i = 1, klon
467        dzcld(i) = 0.0
468      END DO
469      GO TO 70
470    END IF
471
472    ! Current level just below top level => no overshoot
473
474    IF (k<=limcnv + 1) THEN
475      DO i = 1, klon
476        IF (ldcum(i)) THEN
477          cldwtr(i) = sb(i, k) - sc(i) + flotab(i) / (1.0 + gam(i, k))
478          cldwtr(i) = max(0.0, cldwtr(i))
479          beta(i) = 0.0
480        END IF
481      END DO
482      GO TO 20
483    END IF
484
485    ! First guess at overshoot parameter using crude buoyancy closure
486    ! 10% overshoot assumed as a minimum and 1-c0*dz maximum to start
487    ! If pre-existing supersaturation in detrainment layer, beta=0
488    ! cldwtr is temporarily equal to RLVTT*l (l=> liquid water)
489
490    DO i = 1, klon
491      IF (ldcum(i)) THEN
492        cldwtr(i) = sb(i, k) - sc(i) + flotab(i) / (1.0 + gam(i, k))
493        cldwtr(i) = max(0.0, cldwtr(i))
494        betamx = 1.0 - c0 * max(0.0, (dzcld(i) - dzmin))
495        b1 = (hc(i) - hbs(i, km1)) * dp(i, km1)
496        b2 = (hc(i) - hbs(i, k)) * dp(i, k)
497        beta(i) = max(betamn, min(betamx, 1.0 + b1 / b2))
498        IF (hbs(i, km1)<=hb(i, km1)) beta(i) = 0.0
499      END IF
500    END DO
501
502    ! Bound maximum beta to ensure physically realistic solutions
503
504    ! First check constrains beta so that eta remains positive
505    ! (assuming that eta is already positive for beta equal zero)
506    ! La premiere contrainte de beta est que le flux eta doit etre positif.
507
508    DO i = 1, klon
509      IF (ldcum(i)) THEN
510        tmp1 = (1.0 + gam(i, k)) * (sc(i) - sbh(i, kp1) + cldwtr(i)) - &
511                (hbh(i, kp1) - hc(i)) * dp(i, k) / dp(i, kp1)
512        tmp2 = (1.0 + gam(i, k)) * (sc(i) - sbh(i, k))
513        IF ((beta(i) * tmp2 - tmp1)>0.0) THEN
514          betamx = 0.99 * (tmp1 / tmp2)
515          beta(i) = max(0.0, min(betamx, beta(i)))
516        END IF
517
518        ! Second check involves supersaturation of "detrainment layer"
519        ! small amount of supersaturation acceptable (by ssfac factor)
520        ! La 2e contrainte est que la convection ne doit pas sursaturer
521        ! la "detrainment layer", Neanmoins, une petite sursaturation
522        ! est acceptee (facteur ssfac).
523
524        IF (hb(i, km1)<hbs(i, km1)) THEN
525          tmp1 = (1.0 + gam(i, k)) * (sc(i) - sbh(i, kp1) + cldwtr(i)) - &
526                  (hbh(i, kp1) - hc(i)) * dp(i, k) / dp(i, kp1)
527          tmp1 = tmp1 / dp(i, k)
528          tmp2 = gam(i, km1) * (sbh(i, k) - sc(i) + cldwtr(i)) - hbh(i, k) + hc(i) - &
529                  sc(i) + sbh(i, k)
530          tmp3 = (1.0 + gam(i, k)) * (sc(i) - sbh(i, k)) / dp(i, k)
531          tmp4 = (dt / cats) * (hc(i) - hbs(i, k)) * tmp2 / (dp(i, km1) * (hbs(i, km1) - hb(i, &
532                  km1))) + tmp3
533          IF ((beta(i) * tmp4 - tmp1)>0.0) THEN
534            betamx = ssfac * (tmp1 / tmp4)
535            beta(i) = max(0.0, min(betamx, beta(i)))
536          END IF
537        ELSE
538          beta(i) = 0.0
539        END IF
540
541        ! Third check to avoid introducing 2 delta x thermodynamic
542        ! noise in the vertical ... constrain adjusted h (or theta e)
543        ! so that the adjustment doesn't contribute to "kinks" in h
544
545        g = min(0.0, hb(i, k) - hb(i, km1))
546        tmp3 = (hb(i, k) - hb(i, km1) - g) * (cats / dt) / (hc(i) - hbs(i, k))
547        tmp1 = (1.0 + gam(i, k)) * (sc(i) - sbh(i, kp1) + cldwtr(i)) - &
548                (hbh(i, kp1) - hc(i)) * dp(i, k) / dp(i, kp1)
549        tmp1 = tmp1 / dp(i, k)
550        tmp1 = tmp3 * tmp1 + (hc(i) - hbh(i, kp1)) / dp(i, k)
551        tmp2 = tmp3 * (1.0 + gam(i, k)) * (sc(i) - sbh(i, k)) / dp(i, k) + &
552                (hc(i) - hbh(i, k) - cldwtr(i)) * (1.0 / dp(i, k) + 1.0 / dp(i, kp1))
553        IF ((beta(i) * tmp2 - tmp1)>0.0) THEN
554          betamx = 0.0
555          IF (tmp2/=0.0) betamx = tmp1 / tmp2
556          beta(i) = max(0.0, min(betamx, beta(i)))
557        END IF
558      END IF
559    END DO
560
561    ! Calculate mass flux required for stabilization.
562
563    ! Ensure that the convective mass flux, eta, is positive by
564    ! setting negative values of eta to zero..
565    ! Ensure that estimated mass flux cannot move more than the
566    ! minimum of total mass contained in either layer k or layer k+1.
567    ! Also test for other pathological cases that result in non-
568    ! physical states and adjust eta accordingly.
569
570    20  CONTINUE
571    DO i = 1, klon
572      IF (ldcum(i)) THEN
573        beta(i) = max(0.0, beta(i))
574        tmp1 = hc(i) - hbs(i, k)
575        tmp2 = ((1.0 + gam(i, k)) * (sc(i) - sbh(i, kp1) + cldwtr(i)) - beta(i) * (1.0 + gam(&
576                i, k)) * (sc(i) - sbh(i, k))) / dp(i, k) - (hbh(i, kp1) - hc(i)) / dp(i, kp1)
577        eta(i) = tmp1 / (tmp2 * rg * cats)
578        tmass = min(dp(i, k), dp(i, kp1)) / rg
579        IF (eta(i)>tmass * rdt .OR. eta(i)<=0.0) eta(i) = 0.0
580
581        ! Check on negative q in top layer (bound beta)
582
583        IF (shc(i) - shbh(i, k)<0.0 .AND. beta(i) * eta(i)/=0.0) THEN
584          denom = eta(i) * rg * dt * (shc(i) - shbh(i, k)) / dp(i, km1)
585          beta(i) = max(0.0, min(-0.999 * shb(i, km1) / denom, beta(i)))
586        END IF
587
588        ! Check on negative q in middle layer (zero eta)
589
590        qtest1 = shb(i, k) + eta(i) * rg * dt * ((shc(i) - shbh(i, &
591                kp1)) - (1.0 - beta(i)) * cldwtr(i) / rlvtt - beta(i) * (shc(i) - shbh(i, &
592                k))) / dp(i, k)
593        IF (qtest1<=0.0) eta(i) = 0.0
594
595        ! Check on negative q in lower layer (bound eta)
596
597        fac1 = -(shbh(i, kp1) - shc(i)) / dp(i, kp1)
598        qtest2 = shb(i, kp1) - eta(i) * rg * dt * fac1
599        IF (qtest2<0.0) THEN
600          eta(i) = 0.99 * shb(i, kp1) / (rg * dt * fac1)
601        END IF
602      END IF
603    END DO
604
605
606    ! Calculate cloud water, rain water, and thermodynamic changes
607
608    DO i = 1, klon
609      IF (ldcum(i)) THEN
610        etagdt = eta(i) * rg * dt
611        cldwtr(i) = etagdt * cldwtr(i) / rlvtt / rg
612        rnwtr(i) = (1.0 - beta(i)) * cldwtr(i)
613        ds1(i) = etagdt * (sbh(i, kp1) - sc(i)) / dp(i, kp1)
614        dq1(i) = etagdt * (shbh(i, kp1) - shc(i)) / dp(i, kp1)
615        ds2(i) = (etagdt * (sc(i) - sbh(i, kp1)) + rlvtt * rg * cldwtr(i) - beta(i) * etagdt &
616                * (sc(i) - sbh(i, k))) / dp(i, k)
617        dq2(i) = (etagdt * (shc(i) - shbh(i, kp1)) - rg * rnwtr(i) - beta(i) * etagdt * (shc &
618                (i) - shbh(i, k))) / dp(i, k)
619        ds3(i) = beta(i) * (etagdt * (sc(i) - sbh(i, k)) - rlvtt * rg * cldwtr(i)) / dp(i, &
620                km1)
621        dq3(i) = beta(i) * etagdt * (shc(i) - shbh(i, k)) / dp(i, km1)
622
623        ! Isolate convective fluxes for later diagnostics
624
625        fslkp = eta(i) * (sc(i) - sbh(i, kp1))
626        fslkm = beta(i) * (eta(i) * (sc(i) - sbh(i, k)) - rlvtt * cldwtr(i) * rdt)
627        fqlkp = eta(i) * (shc(i) - shbh(i, kp1))
628        fqlkm = beta(i) * eta(i) * (shc(i) - shbh(i, k))
629
630
631        ! Update thermodynamic profile (update sb, hb, & hbs later)
632
633        tb(i, kp1) = tb(i, kp1) + ds1(i) / rcpd
634        tb(i, k) = tb(i, k) + ds2(i) / rcpd
635        tb(i, km1) = tb(i, km1) + ds3(i) / rcpd
636        shb(i, kp1) = shb(i, kp1) + dq1(i)
637        shb(i, k) = shb(i, k) + dq2(i)
638        shb(i, km1) = shb(i, km1) + dq3(i)
639        prec(i) = prec(i) + rnwtr(i) / rhoh2o
640
641        ! Update diagnostic information for final budget
642        ! Tracking temperature & specific humidity tendencies,
643        ! rainout term, convective mass flux, convective liquid
644        ! water static energy flux, and convective total water flux
645
646        cmfdt(i, kp1) = cmfdt(i, kp1) + ds1(i) / rcpd * rdt
647        cmfdt(i, k) = cmfdt(i, k) + ds2(i) / rcpd * rdt
648        cmfdt(i, km1) = cmfdt(i, km1) + ds3(i) / rcpd * rdt
649        cmfdq(i, kp1) = cmfdq(i, kp1) + dq1(i) * rdt
650        cmfdq(i, k) = cmfdq(i, k) + dq2(i) * rdt
651        cmfdq(i, km1) = cmfdq(i, km1) + dq3(i) * rdt
652        cmfdqr(i, k) = cmfdqr(i, k) + (rg * rnwtr(i) / dp(i, k)) * rdt
653        cmfmc(i, kp1) = cmfmc(i, kp1) + eta(i)
654        cmfmc(i, k) = cmfmc(i, k) + beta(i) * eta(i)
655        cmfsl(i, kp1) = cmfsl(i, kp1) + fslkp
656        cmfsl(i, k) = cmfsl(i, k) + fslkm
657        cmflq(i, kp1) = cmflq(i, kp1) + rlvtt * fqlkp
658        cmflq(i, k) = cmflq(i, k) + rlvtt * fqlkm
659        qc(i, k) = (rg * rnwtr(i) / dp(i, k)) * rdt
660      END IF
661    END DO
662
663    ! Next, convectively modify passive constituents
664
665    DO m = 1, pcnst
666      DO i = 1, klon
667        IF (ldcum(i)) THEN
668
669          ! If any of the reported values of the constituent is negative in
670          ! the three adjacent levels, nothing will be done to the profile
671
672          IF ((cmrb(i, kp1, m)<0.0) .OR. (cmrb(i, k, m)<0.0) .OR. (cmrb(i, km1, &
673                  m)<0.0)) GO TO 40
674
675          ! Specify constituent interface values (linear interpolation)
676
677          cmrh(i, k) = 0.5 * (cmrb(i, km1, m) + cmrb(i, k, m))
678          cmrh(i, kp1) = 0.5 * (cmrb(i, k, m) + cmrb(i, kp1, m))
679
680          ! Specify perturbation properties of constituents in PBL
681
682          pblhgt = max(pblh(i), 1.0)
683          IF (gz(i, kp1) / rg<=pblhgt .AND. dzcld(i)==0.) THEN
684            fac1 = max(0.0, 1.0 - gz(i, kp1) / rg / pblhgt)
685            cmrc(i) = cmrb(i, kp1, m) + cmrp(i, m) * fac1
686          ELSE
687            cmrc(i) = cmrb(i, kp1, m)
688          END IF
689
690          ! Determine fluxes, flux divergence => changes due to convection
691          ! Logic must be included to avoid producing negative values. A bit
692          ! messy since there are no a priori assumptions about profiles.
693          ! Tendency is modified (reduced) when pending disaster detected.
694
695          etagdt = eta(i) * rg * dt
696          botflx = etagdt * (cmrc(i) - cmrh(i, kp1))
697          topflx = beta(i) * etagdt * (cmrc(i) - cmrh(i, k))
698          dcmr1(i) = -botflx / dp(i, kp1)
699          efac1 = 1.0
700          efac2 = 1.0
701          efac3 = 1.0
702
703          IF (cmrb(i, kp1, m) + dcmr1(i)<0.0) THEN
704            efac1 = max(tiny, abs(cmrb(i, kp1, m) / dcmr1(i)) - eps)
705          END IF
706
707          IF (efac1==tiny .OR. efac1>1.0) efac1 = 0.0
708          dcmr1(i) = -efac1 * botflx / dp(i, kp1)
709          dcmr2(i) = (efac1 * botflx - topflx) / dp(i, k)
710
711          IF (cmrb(i, k, m) + dcmr2(i)<0.0) THEN
712            efac2 = max(tiny, abs(cmrb(i, k, m) / dcmr2(i)) - eps)
713          END IF
714
715          IF (efac2==tiny .OR. efac2>1.0) efac2 = 0.0
716          dcmr2(i) = (efac1 * botflx - efac2 * topflx) / dp(i, k)
717          dcmr3(i) = efac2 * topflx / dp(i, km1)
718
719          IF (cmrb(i, km1, m) + dcmr3(i)<0.0) THEN
720            efac3 = max(tiny, abs(cmrb(i, km1, m) / dcmr3(i)) - eps)
721          END IF
722
723          IF (efac3==tiny .OR. efac3>1.0) efac3 = 0.0
724          efac3 = min(efac2, efac3)
725          dcmr2(i) = (efac1 * botflx - efac3 * topflx) / dp(i, k)
726          dcmr3(i) = efac3 * topflx / dp(i, km1)
727
728          cmrb(i, kp1, m) = cmrb(i, kp1, m) + dcmr1(i)
729          cmrb(i, k, m) = cmrb(i, k, m) + dcmr2(i)
730          cmrb(i, km1, m) = cmrb(i, km1, m) + dcmr3(i)
731        END IF
732      40    END DO
733    END DO ! end of m=1,pcnst loop
734
735    IF (k==limcnv + 1) GO TO 60 ! on ne pourra plus glisser
736
737    ! Dans la procedure de glissage ascendant, les variables thermo-
738    ! dynamiques des couches k et km1 servent au calcul des couches
739    ! superieures. Elles ont donc besoin d'une mise-a-jour.
740
741    DO i = 1, klon
742      IF (ldcum(i)) THEN
743        zx_t = tb(i, k)
744        zx_p = p(i, k)
745        zx_q = shb(i, k)
746        zdelta = max(0., sign(1., rtt - zx_t))
747        zcvm5 = r5les * rlvtt * (1. - zdelta) + r5ies * rlstt * zdelta
748        zcvm5 = zcvm5 / rcpd / (1.0 + rvtmp2 * zx_q)
749        zx_qs = r2es * foeew(zx_t, zdelta) / zx_p
750        zx_qs = min(0.5, zx_qs)
751        zcor = 1. / (1. - retv * zx_qs)
752        zx_qs = zx_qs * zcor
753        zx_gam = foede(zx_t, zdelta, zcvm5, zx_qs, zcor)
754        shbs(i, k) = zx_qs
755        gam(i, k) = zx_gam
756
757        zx_t = tb(i, km1)
758        zx_p = p(i, km1)
759        zx_q = shb(i, km1)
760        zdelta = max(0., sign(1., rtt - zx_t))
761        zcvm5 = r5les * rlvtt * (1. - zdelta) + r5ies * rlstt * zdelta
762        zcvm5 = zcvm5 / rcpd / (1.0 + rvtmp2 * zx_q)
763        zx_qs = r2es * foeew(zx_t, zdelta) / zx_p
764        zx_qs = min(0.5, zx_qs)
765        zcor = 1. / (1. - retv * zx_qs)
766        zx_qs = zx_qs * zcor
767        zx_gam = foede(zx_t, zdelta, zcvm5, zx_qs, zcor)
768        shbs(i, km1) = zx_qs
769        gam(i, km1) = zx_gam
770
771        sb(i, k) = sb(i, k) + ds2(i)
772        sb(i, km1) = sb(i, km1) + ds3(i)
773        hb(i, k) = sb(i, k) + rlvtt * shb(i, k)
774        hb(i, km1) = sb(i, km1) + rlvtt * shb(i, km1)
775        hbs(i, k) = sb(i, k) + rlvtt * shbs(i, k)
776        hbs(i, km1) = sb(i, km1) + rlvtt * shbs(i, km1)
777
778        sbh(i, k) = 0.5 * (sb(i, k) + sb(i, km1))
779        shbh(i, k) = qhalf(shb(i, km1), shb(i, k), shbs(i, km1), shbs(i, k))
780        hbh(i, k) = sbh(i, k) + rlvtt * shbh(i, k)
781        sbh(i, km1) = 0.5 * (sb(i, km1) + sb(i, k - 2))
782        shbh(i, km1) = qhalf(shb(i, k - 2), shb(i, km1), shbs(i, k - 2), &
783                shbs(i, km1))
784        hbh(i, km1) = sbh(i, km1) + rlvtt * shbh(i, km1)
785      END IF
786    END DO
787
788    ! Ensure that dzcld is reset if convective mass flux zero
789    ! specify the current vertical extent of the convective activity
790    ! top of convective layer determined by size of overshoot param.
791
792    60  CONTINUE
793    DO i = 1, klon
794      etagt0 = eta(i) > 0.0
795      IF (.NOT. etagt0) dzcld(i) = 0.0
796      IF (etagt0 .AND. beta(i)>betamn) THEN
797        ktp = km1
798      ELSE
799        ktp = k
800      END IF
801      IF (etagt0) THEN
802        cnt(i) = min(cnt(i), ktp)
803        cnb(i) = max(cnb(i), k)
804      END IF
805    END DO
806  70 END DO ! end of k loop
807
808  ! determine whether precipitation, prec, is frozen (snow) or not
809
810  DO i = 1, klon
811    IF (tb(i, klev)<tmelt .AND. tb(i, klev - 1)<tmelt) THEN
812      cmfprs(i) = prec(i) * rdt
813    ELSE
814      cmfprt(i) = prec(i) * rdt
815    END IF
816  END DO
817
818  RETURN ! we're all done ... return to calling procedure
819END SUBROUTINE cmfmca
Note: See TracBrowser for help on using the repository browser.