source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phylmd/cv_driver.F90 @ 5434

Last change on this file since 5434 was 2481, checked in by fhourdin, 9 years ago

Introduction d'une dependance epmax=f(Cape) sur proposition de Camille Risi

  • 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: 25.3 KB
Line 
1
2! $Header$
3
4SUBROUTINE cv_driver(len, nd, ndp1, ntra, iflag_con, t1, q1, qs1, u1, v1, &
5    tra1, p1, ph1, iflag1, ft1, fq1, fu1, fv1, ftra1, precip1, vprecip1, &
6    cbmf1, sig1, w01, icb1, inb1, delt, ma1, upwd1, dnwd1, dnwd01, qcondc1, &
7    wd1, cape1, da1, phi1, mp1, phi21, d1a1, dam1, sij1, clw1, elij1, & !
8                                                                        ! RomP
9    evap1, ep1, epmlmmm1, eplamm1, & ! RomP
10    wdtraina1, wdtrainm1, & ! RomP
11    epmax_diag1) ! epmax_cape
12
13  USE dimphy
14  IMPLICIT NONE
15
16  ! .............................START PROLOGUE............................
17
18
19  ! All argument names (except len,nd,ntra,nloc,delt and the flags) have a
20  ! "1" appended.
21  ! The "1" is removed for the corresponding compressed (local) variables.
22
23  ! PARAMETERS:
24  ! Name            Type         Usage            Description
25  ! ----------      ----------     -------  ----------------------------
26
27  ! len           Integer        Input        first (i) dimension
28  ! nd            Integer        Input        vertical (k) dimension
29  ! ndp1          Integer        Input        nd + 1
30  ! ntra          Integer        Input        number of tracors
31  ! iflag_con     Integer        Input        version of convect (3/4)
32  ! t1            Real           Input        temperature
33  ! q1            Real           Input        specific hum
34  ! qs1           Real           Input        sat specific hum
35  ! u1            Real           Input        u-wind
36  ! v1            Real           Input        v-wind
37  ! tra1          Real           Input        tracors
38  ! p1            Real           Input        full level pressure
39  ! ph1           Real           Input        half level pressure
40  ! iflag1        Integer        Output       flag for Emanuel conditions
41  ! ft1           Real           Output       temp tend
42  ! fq1           Real           Output       spec hum tend
43  ! fu1           Real           Output       u-wind tend
44  ! fv1           Real           Output       v-wind tend
45  ! ftra1         Real           Output       tracor tend
46  ! precip1       Real           Output       precipitation
47  ! VPrecip1      Real           Output       vertical profile of
48  ! precipitations
49  ! cbmf1         Real           Output       cloud base mass flux
50  ! sig1          Real           In/Out       section adiabatic updraft
51  ! w01           Real           In/Out       vertical velocity within adiab
52  ! updraft
53  ! delt          Real           Input        time step
54  ! Ma1           Real           Output       mass flux adiabatic updraft
55  ! upwd1         Real           Output       total upward mass flux
56  ! (adiab+mixed)
57  ! dnwd1         Real           Output       saturated downward mass flux
58  ! (mixed)
59  ! dnwd01        Real           Output       unsaturated downward mass flux
60  ! qcondc1       Real           Output       in-cld mixing ratio of
61  ! condensed water
62  ! wd1           Real           Output       downdraft velocity scale for
63  ! sfc fluxes
64  ! cape1         Real           Output       CAPE
65
66  ! wdtrainA1     Real           Output   precipitation detrained from
67  ! adiabatic draught;
68  ! used in tracer transport (cvltr)
69  ! wdtrainM1     Real           Output   precipitation detrained from mixed
70  ! draughts;
71  ! used in tracer transport (cvltr)
72  ! da1           Real           Output   used in tracer transport (cvltr)
73  ! phi1          Real           Output   used in tracer transport (cvltr)
74  ! mp1           Real           Output   used in tracer transport (cvltr)
75
76  ! phi21         Real           Output   used in tracer transport (cvltr)
77
78  ! d1a1          Real           Output   used in tracer transport (cvltr)
79  ! dam1          Real           Output   used in tracer transport (cvltr)
80
81  ! evap1         Real           Output
82  ! ep1           Real           Output
83  ! sij1        Real           Output
84  ! elij1         Real           Output
85
86  ! S. Bony, Mar 2002:
87  ! * Several modules corresponding to different physical processes
88  ! * Several versions of convect may be used:
89  ! - iflag_con=3: version lmd  (previously named convect3)
90  ! - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
91  ! + tard:     - iflag_con=5: version lmd with ice (previously named convectg)
92  ! S. Bony, Oct 2002:
93  ! * Vectorization of convect3 (ie version lmd)
94
95  ! ..............................END PROLOGUE.............................
96
97
98  ! Input
99  INTEGER len
100  INTEGER nd
101  INTEGER ndp1
102  INTEGER noff
103  INTEGER iflag_con
104  INTEGER ntra
105  REAL delt
106  REAL t1(len, nd)
107  REAL q1(len, nd)
108  REAL qs1(len, nd)
109  REAL u1(len, nd)
110  REAL v1(len, nd)
111  REAL tra1(len, nd, ntra)
112  REAL p1(len, nd)
113  REAL ph1(len, ndp1)
114
115  ! Output
116  INTEGER iflag1(len)
117  REAL ft1(len, nd)
118  REAL fq1(len, nd)
119  REAL fu1(len, nd)
120  REAL fv1(len, nd)
121  REAL ftra1(len, nd, ntra)
122  REAL precip1(len)
123  REAL cbmf1(len)
124  REAL sig1(klon, klev)
125  REAL w01(klon, klev)
126  REAL vprecip1(len, nd+1)
127  REAL evap1(len, nd) !RomP
128  REAL ep1(len, nd) !RomP
129  REAL ma1(len, nd)
130  REAL upwd1(len, nd)
131  REAL dnwd1(len, nd)
132  REAL dnwd01(len, nd)
133
134  REAL qcondc1(len, nd) ! cld
135  REAL wd1(len) ! gust
136  REAL cape1(len)
137
138  ! RomP >>>
139  REAL wdtraina1(len, nd), wdtrainm1(len, nd)
140  REAL sij1(len, nd, nd), elij1(len, nd, nd)
141  REAL da1(len, nd), phi1(len, nd, nd), mp1(len, nd)
142
143  REAL phi21(len, nd, nd)
144  REAL d1a1(len, nd), dam1(len, nd)
145  REAL epmlmmm1(len, nd, nd), eplamm1(len, nd)
146  ! RomP <<<
147  REAL epmax_diag1 (len) ! epmax_cape     
148
149  ! -------------------------------------------------------------------
150  ! Original Prologue by Kerry Emanuel.
151  ! -------------------------------------------------------------------
152  ! --- ARGUMENTS
153  ! -------------------------------------------------------------------
154  ! --- On input:
155
156  ! t:   Array of absolute temperature (K) of dimension ND, with first
157  ! index corresponding to lowest model level. Note that this array
158  ! will be altered by the subroutine if dry convective adjustment
159  ! occurs and if IPBL is not equal to 0.
160
161  ! q:   Array of specific humidity (gm/gm) of dimension ND, with first
162  ! index corresponding to lowest model level. Must be defined
163  ! at same grid levels as T. Note that this array will be altered
164  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
165
166  ! qs:  Array of saturation specific humidity of dimension ND, with first
167  ! index corresponding to lowest model level. Must be defined
168  ! at same grid levels as T. Note that this array will be altered
169  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
170
171  ! u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
172  ! index corresponding with the lowest model level. Defined at
173  ! same levels as T. Note that this array will be altered if
174  ! dry convective adjustment occurs and if IPBL is not equal to 0.
175
176  ! v:   Same as u but for meridional velocity.
177
178  ! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
179  ! where NTRA is the number of different tracers. If no
180  ! convective tracer transport is needed, define a dummy
181  ! input array of dimension (ND,1). Tracers are defined at
182  ! same vertical levels as T. Note that this array will be altered
183  ! if dry convective adjustment occurs and if IPBL is not equal to 0.
184
185  ! p:   Array of pressure (mb) of dimension ND, with first
186  ! index corresponding to lowest model level. Must be defined
187  ! at same grid levels as T.
188
189  ! ph:  Array of pressure (mb) of dimension ND+1, with first index
190  ! corresponding to lowest level. These pressures are defined at
191  ! levels intermediate between those of P, T, Q and QS. The first
192  ! value of PH should be greater than (i.e. at a lower level than)
193  ! the first value of the array P.
194
195  ! nl:  The maximum number of levels to which convection can penetrate, plus
196  ! 1.
197  ! NL MUST be less than or equal to ND-1.
198
199  ! delt: The model time step (sec) between calls to CONVECT
200
201  ! ----------------------------------------------------------------------------
202  ! ---   On Output:
203
204  ! iflag: An output integer whose value denotes the following:
205  ! VALUE   INTERPRETATION
206  ! -----   --------------
207  ! 0     Moist convection occurs.
208  ! 1     Moist convection occurs, but a CFL condition
209  ! on the subsidence warming is violated. This
210  ! does not cause the scheme to terminate.
211  ! 2     Moist convection, but no precip because ep(inb) lt 0.0001
212  ! 3     No moist convection because new cbmf is 0 and old cbmf is 0.
213  ! 4     No moist convection; atmosphere is not
214  ! unstable
215  ! 6     No moist convection because ihmin le minorig.
216  ! 7     No moist convection because unreasonable
217  ! parcel level temperature or specific humidity.
218  ! 8     No moist convection: lifted condensation
219  ! level is above the 200 mb level.
220  ! 9     No moist convection: cloud base is higher
221  ! then the level NL-1.
222
223  ! ft:   Array of temperature tendency (K/s) of dimension ND, defined at
224  ! same
225  ! grid levels as T, Q, QS and P.
226
227  ! fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
228  ! defined at same grid levels as T, Q, QS and P.
229
230  ! fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
231  ! defined at same grid levels as T.
232
233  ! fv:   Same as FU, but for forcing of meridional velocity.
234
235  ! ftra: Array of forcing of tracer content, in tracer mixing ratio per
236  ! second, defined at same levels as T. Dimensioned (ND,NTRA).
237
238  ! precip: Scalar convective precipitation rate (mm/day).
239
240  ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
241
242  ! wd:   A convective downdraft velocity scale. For use in surface
243  ! flux parameterizations. See convect.ps file for details.
244
245  ! tprime: A convective downdraft temperature perturbation scale (K).
246  ! For use in surface flux parameterizations. See convect.ps
247  ! file for details.
248
249  ! qprime: A convective downdraft specific humidity
250  ! perturbation scale (gm/gm).
251  ! For use in surface flux parameterizations. See convect.ps
252  ! file for details.
253
254  ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
255  ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
256  ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"
257  ! by the calling program between calls to CONVECT.
258
259  ! det:   Array of detrainment mass flux of dimension ND.
260
261  ! -------------------------------------------------------------------
262
263  ! Local arrays
264
265
266  INTEGER i, k, n, il, j
267  INTEGER icbmax
268  INTEGER nk1(klon)
269  INTEGER icb1(klon)
270  INTEGER inb1(klon)
271  INTEGER icbs1(klon)
272
273  REAL plcl1(klon)
274  REAL tnk1(klon)
275  REAL qnk1(klon)
276  REAL gznk1(klon)
277  REAL pnk1(klon)
278  REAL qsnk1(klon)
279  REAL pbase1(klon)
280  REAL buoybase1(klon)
281
282  REAL lv1(klon, klev)
283  REAL cpn1(klon, klev)
284  REAL tv1(klon, klev)
285  REAL gz1(klon, klev)
286  REAL hm1(klon, klev)
287  REAL h1(klon, klev)
288  REAL tp1(klon, klev)
289  REAL tvp1(klon, klev)
290  REAL clw1(klon, klev)
291  REAL th1(klon, klev)
292
293  INTEGER ncum
294
295  ! (local) compressed fields:
296
297  ! ym      integer nloc
298  ! ym      parameter (nloc=klon) ! pour l'instant
299#define nloc klon
300  INTEGER idcum(nloc)
301  INTEGER iflag(nloc), nk(nloc), icb(nloc)
302  INTEGER nent(nloc, klev)
303  INTEGER icbs(nloc)
304  INTEGER inb(nloc), inbis(nloc)
305
306  REAL cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
307  REAL t(nloc, klev), q(nloc, klev), qs(nloc, klev)
308  REAL u(nloc, klev), v(nloc, klev)
309  REAL gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)
310  REAL p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)
311  REAL clw(nloc, klev)
312  REAL dph(nloc, klev)
313  REAL pbase(nloc), buoybase(nloc), th(nloc, klev)
314  REAL tvp(nloc, klev)
315  REAL sig(nloc, klev), w0(nloc, klev)
316  REAL hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)
317  REAL frac(nloc), buoy(nloc, klev)
318  REAL cape(nloc)
319  REAL m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)
320  REAL uent(nloc, klev, klev), vent(nloc, klev, klev)
321  REAL ments(nloc, klev, klev), qents(nloc, klev, klev)
322  REAL sij(nloc, klev, klev), elij(nloc, klev, klev)
323  REAL qp(nloc, klev), up(nloc, klev), vp(nloc, klev)
324  REAL wt(nloc, klev), water(nloc, klev), evap(nloc, klev)
325  REAL b(nloc, klev), ft(nloc, klev), fq(nloc, klev)
326  REAL fu(nloc, klev), fv(nloc, klev)
327  REAL upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)
328  REAL ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)
329  REAL tps(nloc, klev), qprime(nloc), tprime(nloc)
330  REAL precip(nloc)
331  REAL vprecip(nloc, klev+1)
332  REAL tra(nloc, klev, ntra), trap(nloc, klev, ntra)
333  REAL ftra(nloc, klev, ntra), traent(nloc, klev, klev, ntra)
334  REAL qcondc(nloc, klev) ! cld
335  REAL wd(nloc) ! gust
336
337  ! RomP >>>
338  REAL da(nloc, klev), phi(nloc, klev, klev), mp(nloc, klev)
339  REAL epmlmmm(nloc, klev, klev), eplamm(nloc, klev)
340  REAL phi2(nloc, klev, klev)
341  REAL d1a(nloc, klev), dam(nloc, klev)
342  REAL wdtraina(nloc, klev), wdtrainm(nloc, klev)
343  REAL sigd(nloc)
344  ! RomP <<<
345  REAL epmax_diag(nloc) ! epmax_cape
346
347  nent(:, :) = 0
348  ! -------------------------------------------------------------------
349  ! --- SET CONSTANTS AND PARAMETERS
350  ! -------------------------------------------------------------------
351  ! print *, '-> cv_driver'      !jyg
352  ! -- set simulation flags:
353  ! (common cvflag)
354
355  CALL cv_flag(0)
356
357  ! -- set thermodynamical constants:
358  ! (common cvthermo)
359
360  CALL cv_thermo(iflag_con)
361
362  ! -- set convect parameters
363
364  ! includes microphysical parameters and parameters that
365  ! control the rate of approach to quasi-equilibrium)
366  ! (common cvparam)
367
368
369  IF (iflag_con==30) THEN
370    CALL cv30_param(nd, delt)
371  END IF
372
373  IF (iflag_con==4) THEN
374    CALL cv_param(nd)
375  END IF
376
377  ! ---------------------------------------------------------------------
378  ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
379  ! ---------------------------------------------------------------------
380
381  ft1(:, :) = 0.0
382  fq1(:, :) = 0.0
383  fu1(:, :) = 0.0
384  fv1(:, :) = 0.0
385  tvp1(:, :) = 0.0
386  tp1(:, :) = 0.0
387  clw1(:, :) = 0.0
388  ! ym
389  clw(:, :) = 0.0
390  gz1(:, :) = 0.
391  vprecip1(:, :) = 0.
392  ma1(:, :) = 0.0
393  upwd1(:, :) = 0.0
394  dnwd1(:, :) = 0.0
395  dnwd01(:, :) = 0.0
396  qcondc1(:, :) = 0.0
397
398  ftra1(:, :, :) = 0.0
399
400  elij1(:, :, :) = 0.0
401  sij1(:, :, :) = 0.0
402
403  precip1(:) = 0.0
404  iflag1(:) = 0
405  wd1(:) = 0.0
406  cape1(:) = 0.0
407  epmax_diag1(:) = 0.0 ! epmax_cape
408
409
410  IF (iflag_con==30) THEN
411    DO il = 1, len
412      sig1(il, nd) = sig1(il, nd) + 1.
413      sig1(il, nd) = amin1(sig1(il,nd), 12.1)
414    END DO
415  END IF
416
417  ! RomP >>>
418  wdtraina1(:, :) = 0.
419  wdtrainm1(:, :) = 0.
420  da1(:, :) = 0.
421  phi1(:, :, :) = 0.
422  epmlmmm1(:, :, :) = 0.
423  eplamm1(:, :) = 0.
424  mp1(:, :) = 0.
425  evap1(:, :) = 0.
426  ep1(:, :) = 0.
427  sij1(:, :, :) = 0.
428  elij1(:, :, :) = 0.
429  phi21(:, :, :) = 0.
430  d1a1(:, :) = 0.
431  dam1(:, :) = 0.
432  ! RomP <<<
433
434  ! --------------------------------------------------------------------
435  ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
436  ! --------------------------------------------------------------------
437
438  IF (iflag_con==30) THEN
439
440    ! print*,'Emanuel version 30 '
441    CALL cv30_prelim(len, nd, ndp1, t1, q1, p1, ph1 & ! nd->na
442      , lv1, cpn1, tv1, gz1, h1, hm1, th1)
443  END IF
444
445  IF (iflag_con==4) THEN
446    CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, &
447      hm1)
448  END IF
449
450  ! --------------------------------------------------------------------
451  ! --- CONVECTIVE FEED
452  ! --------------------------------------------------------------------
453
454  IF (iflag_con==30) THEN
455    CALL cv30_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1 & !
456                                                             ! nd->na
457      , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
458  END IF
459
460  IF (iflag_con==4) THEN
461    CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1, nk1, icb1, icbmax, &
462      iflag1, tnk1, qnk1, gznk1, plcl1)
463  END IF
464
465  ! --------------------------------------------------------------------
466  ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
467  ! (up through ICB for convect4, up through ICB+1 for convect3)
468  ! Calculates the lifted parcel virtual temperature at nk, the
469  ! actual temperature, and the adiabatic liquid water content.
470  ! --------------------------------------------------------------------
471
472  IF (iflag_con==30) THEN
473    CALL cv30_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1 & ! nd->na
474      , tp1, tvp1, clw1, icbs1)
475  END IF
476
477  IF (iflag_con==4) THEN
478    CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, tp1, &
479      tvp1, clw1)
480  END IF
481
482  ! -------------------------------------------------------------------
483  ! --- TRIGGERING
484  ! -------------------------------------------------------------------
485
486  IF (iflag_con==30) THEN
487    CALL cv30_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1 & !
488                                                                 ! nd->na
489      , pbase1, buoybase1, iflag1, sig1, w01)
490  END IF
491
492  IF (iflag_con==4) THEN
493    CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
494  END IF
495
496  ! =====================================================================
497  ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
498  ! =====================================================================
499
500  ncum = 0
501  DO i = 1, len
502    IF (iflag1(i)==0) THEN
503      ncum = ncum + 1
504      idcum(ncum) = i
505    END IF
506  END DO
507
508  ! print*,'cv_driver : klon, ncum = ',len,ncum
509
510  IF (ncum>0) THEN
511
512    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
513    ! --- COMPRESS THE FIELDS
514    ! (-> vectorization over convective gridpoints)
515    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
516
517    IF (iflag_con==30) THEN
518      CALL cv30_compress(len, nloc, ncum, nd, ntra, iflag1, nk1, icb1, icbs1, &
519        plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, &
520        gz1, th1, tra1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, &
521        w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, &
522        q, qs, u, v, gz, th, tra, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, &
523        w0)
524    END IF
525
526    IF (iflag_con==4) THEN
527      CALL cv_compress(len, nloc, ncum, nd, iflag1, nk1, icb1, cbmf1, plcl1, &
528        tnk1, qnk1, gznk1, t1, q1, qs1, u1, v1, gz1, h1, lv1, cpn1, p1, ph1, &
529        tv1, tp1, tvp1, clw1, iflag, nk, icb, cbmf, plcl, tnk, qnk, gznk, t, &
530        q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw, dph)
531    END IF
532
533    ! -------------------------------------------------------------------
534    ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
535    ! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
536    ! ---   &
537    ! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
538    ! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
539    ! ---   &
540    ! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
541    ! -------------------------------------------------------------------
542
543    IF (iflag_con==30) THEN
544      CALL cv30_undilute2(nloc, ncum, nd, icb, icbs, nk & !na->nd
545        , tnk, qnk, gznk, t, q, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, &
546        inb, tp, tvp, clw, hp, ep, sigp, buoy)
547    END IF
548
549    IF (iflag_con==4) THEN
550      CALL cv_undilute2(nloc, ncum, nd, icb, nk, tnk, qnk, gznk, t, q, qs, &
551        gz, p, dph, h, tv, lv, inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
552    END IF
553
554    ! -------------------------------------------------------------------
555    ! --- CLOSURE
556    ! -------------------------------------------------------------------
557
558    IF (iflag_con==30) THEN
559      CALL cv30_closure(nloc, ncum, nd, icb, inb & ! na->nd
560        , pbase, p, ph, tv, buoy, sig, w0, cape, m)
561
562      ! epmax_cape
563      call cv30_epmax_fn_cape(nloc,ncum,nd &
564                ,cape,ep,hp,icb,inb,clw,nk,t,h,lv &
565                ,epmax_diag)
566        ! on écrase ep et recalcule hp
567    END IF
568
569    IF (iflag_con==4) THEN
570      CALL cv_closure(nloc, ncum, nd, nk, icb, tv, tvp, p, ph, dph, plcl, &
571        cpn, iflag, cbmf)
572    END IF
573   
574
575    ! -------------------------------------------------------------------
576    ! --- MIXING
577    ! -------------------------------------------------------------------
578
579    IF (iflag_con==30) THEN
580      CALL cv30_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb & !
581                                                                ! na->nd
582        , ph, t, q, qs, u, v, tra, h, lv, qnk, hp, tv, tvp, ep, clw, m, sig, &
583        ment, qent, uent, vent, sij, elij, ments, qents, traent)
584    END IF
585
586    IF (iflag_con==4) THEN
587      CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis, ph, t, q, qs, u, v, &
588        h, lv, qnk, hp, tv, tvp, ep, clw, cbmf, m, ment, qent, uent, vent, &
589        nent, sij, elij)
590    END IF
591
592    ! -------------------------------------------------------------------
593    ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
594    ! -------------------------------------------------------------------
595
596    IF (iflag_con==30) THEN
597      ! RomP >>>
598      CALL cv30_unsat(nloc, ncum, nd, nd, ntra, icb, inb & ! na->nd
599        , t, q, qs, gz, u, v, tra, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, &
600        ment, elij, delt, plcl, mp, qp, up, vp, trap, wt, water, evap, b, &
601        wdtraina, wdtrainm)
602      ! RomP <<<
603    END IF
604
605    IF (iflag_con==4) THEN
606      CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph, h, lv, &
607        ep, sigp, clw, m, ment, elij, iflag, mp, qp, up, vp, wt, water, evap)
608    END IF
609
610    ! -------------------------------------------------------------------
611    ! --- YIELD
612    ! (tendencies, precipitation, variables of interface with other
613    ! processes, etc)
614    ! -------------------------------------------------------------------
615
616    IF (iflag_con==30) THEN
617      CALL cv30_yield(nloc, ncum, nd, nd, ntra & ! na->nd
618        , icb, inb, delt, t, q, u, v, tra, gz, p, ph, h, hp, lv, cpn, th, ep, &
619        clw, m, tp, mp, qp, up, vp, trap, wt, water, evap, b, ment, qent, &
620        uent, vent, nent, elij, traent, sig, tv, tvp, iflag, precip, vprecip, &
621        ft, fq, fu, fv, ftra, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, &
622        wd)
623    END IF
624
625    IF (iflag_con==4) THEN
626      CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt, t, q, u, v, gz, p, &
627        ph, h, hp, lv, cpn, ep, clw, frac, m, mp, qp, up, vp, wt, water, &
628        evap, ment, qent, uent, vent, nent, elij, tv, tvp, iflag, wd, qprime, &
629        tprime, precip, cbmf, ft, fq, fu, fv, ma, qcondc)
630    END IF
631
632    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
633    ! --- passive tracers
634    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
635
636    IF (iflag_con==30) THEN
637      ! RomP >>>
638      CALL cv30_tracer(nloc, len, ncum, nd, nd, ment, sij, da, phi, phi2, &
639        d1a, dam, ep, vprecip, elij, clw, epmlmmm, eplamm, icb, inb)
640      ! RomP <<<
641    END IF
642
643    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
644    ! --- UNCOMPRESS THE FIELDS
645    ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
646    ! set iflag1 =42 for non convective points
647    DO i = 1, len
648      iflag1(i) = 42
649    END DO
650
651    IF (iflag_con==30) THEN
652      CALL cv30_uncompress(nloc, len, ncum, nd, ntra, idcum, iflag, precip, &
653        vprecip, evap, ep, sig, w0 & !RomP
654        , ft, fq, fu, fv, ftra, inb, ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
655        da, phi, mp, phi2, d1a, dam, sij & !RomP
656        , elij, clw, epmlmmm, eplamm & !RomP
657        , wdtraina, wdtrainm,epmax_diag &     !RomP
658        , iflag1, precip1, vprecip1, evap1, ep1, sig1, w01 & !RomP
659        , ft1, fq1, fu1, fv1, ftra1, inb1, ma1, upwd1, dnwd1, dnwd01, &
660        qcondc1, wd1, cape1, da1, phi1, mp1, phi21, d1a1, dam1, sij1 & !RomP
661        , elij1, clw1, epmlmmm1, eplamm1 & !RomP
662        , wdtraina1, wdtrainm1,epmax_diag1) !RomP
663    END IF
664
665    IF (iflag_con==4) THEN
666      CALL cv_uncompress(nloc, len, ncum, nd, idcum, iflag, precip, cbmf, ft, &
667        fq, fu, fv, ma, qcondc, iflag1, precip1, cbmf1, ft1, fq1, fu1, fv1, &
668        ma1, qcondc1)
669    END IF
670
671  END IF ! ncum>0
672
673  ! print *, 'fin cv_driver ->'      !jyg
674  RETURN
675END SUBROUTINE cv_driver
676
677! ==================================================================
678SUBROUTINE cv_flag(iflag_ice_thermo)
679  IMPLICIT NONE
680
681  ! Argument : iflag_ice_thermo : ice thermodynamics is taken into account if
682  ! iflag_ice_thermo >=1
683  INTEGER iflag_ice_thermo
684
685  include "cvflag.h"
686
687  ! -- si .TRUE., on rend la gravite plus explicite et eventuellement
688  ! differente de 10.0 dans convect3:
689  cvflag_grav = .TRUE.
690  cvflag_ice = iflag_ice_thermo >= 1
691
692  RETURN
693END SUBROUTINE cv_flag
694
695! ==================================================================
696SUBROUTINE cv_thermo(iflag_con)
697  IMPLICIT NONE
698
699  ! -------------------------------------------------------------
700  ! Set thermodynamical constants for convectL
701  ! -------------------------------------------------------------
702
703  include "YOMCST.h"
704  include "cvthermo.h"
705
706  INTEGER iflag_con
707
708
709  ! original set from convect:
710  IF (iflag_con==4) THEN
711    cpd = 1005.7
712    cpv = 1870.0
713    cl = 4190.0
714    rrv = 461.5
715    rrd = 287.04
716    lv0 = 2.501E6
717    g = 9.8
718    t0 = 273.15
719    grav = g
720  ELSE
721
722    ! constants consistent with LMDZ:
723    cpd = rcpd
724    cpv = rcpv
725    cl = rcw
726    ci = rcs
727    rrv = rv
728    rrd = rd
729    lv0 = rlvtt
730    lf0 = rlstt - rlvtt
731    g = rg ! not used in convect3
732    ! ori      t0  = RTT
733    t0 = 273.15 ! convect3 (RTT=273.16)
734    ! maf       grav= 10.    ! implicitely or explicitely used in convect3
735    grav = g ! implicitely or explicitely used in convect3
736  END IF
737
738  rowl = 1000.0 !(a quelle variable de YOMCST cela correspond-il?)
739
740  clmcpv = cl - cpv
741  clmcpd = cl - cpd
742  clmci = cl - ci
743  cpdmcp = cpd - cpv
744  cpvmcpd = cpv - cpd
745  cpvmcl = cl - cpv ! for convect3
746  eps = rrd/rrv
747  epsi = 1.0/eps
748  epsim1 = epsi - 1.0
749  ! ginv=1.0/g
750  ginv = 1.0/grav
751  hrd = 0.5*rrd
752
753  RETURN
754END SUBROUTINE cv_thermo
Note: See TracBrowser for help on using the repository browser.