source: LMDZ6/trunk/libf/phylmd/cv_driver_mod.F90 @ 6087

Last change on this file since 6087 was 6058, checked in by fhourdin, 6 weeks ago

Travail pour la replayisation de la convection

Reunion de tous les anciens common devenus modules, dans lmdz_cv_ini.
Déplacement de presque toutes les routines d'initialisation dans lmdz_cv_ini.
Encapsulage de certains sous-programmes dans des modules.
Suppression de programmes inutilisés (cv3_crit et cv3_incp)
Reste :

  • à sortir des routines d'initialisation "_pre" de cv_driver et

cva_driver

  • à passer le variables argunement en intent(in/out/inout).

La convergence numérique a été testée pour
iflag_con=3/30/4
en 3D parallèle.
La compilation de la version isotopique fonctionne.

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