source: LMDZ5/branches/testing/libf/phylmd/cva_driver.F90 @ 5446

Last change on this file since 5446 was 2870, checked in by Laurent Fairhead, 8 years ago

Merged trunk changes r2842:2865 into testing branch

  • 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: 49.6 KB
Line 
1
2! $Id: cva_driver.F90 2870 2017-05-04 07:31:05Z fhourdin $
3
4SUBROUTINE cva_driver(len, nd, ndp1, ntra, nloc, k_upper, &
5                      iflag_con, iflag_mix, iflag_ice_thermo, iflag_clos, ok_conserv_q, &
6!!                      delt, t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, &  ! jyg
7                      delt, comp_threshold, &                                      ! jyg
8                      t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, &          ! jyg
9                      u1, v1, tra1, &
10                      p1, ph1, &
11                      Ale1, Alp1, omega1, &
12                      sig1feed1, sig2feed1, wght1, &
13                      iflag1, ft1, fq1, fu1, fv1, ftra1, &
14                      precip1, kbas1, ktop1, &
15                      cbmf1, plcl1, plfc1, wbeff1, &
16                      sig1, w01, & !input/output
17                      ptop21, sigd1, &
18                      ma1, mip1, Vprecip1, Vprecipi1, upwd1, dnwd1, dnwd01, &      ! jyg
19                      qcondc1, wd1, &
20                      cape1, cin1, tvp1, &
21                      ftd1, fqd1, &
22                      Plim11, Plim21, asupmax1, supmax01, asupmaxmin1, &
23                      lalim_conv1, &
24!!                      da1,phi1,mp1,phi21,d1a1,dam1,sigij1,clw1, &        ! RomP
25!!                      elij1,evap1,ep1,epmlmMm1,eplaMm1, &                ! RomP
26                      da1, phi1, mp1, phi21, d1a1, dam1, sigij1, wghti1, & ! RomP, RL
27                      clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, &        ! RomP, RL
28                      wdtrainA1, wdtrainM1, qtc1, sigt1, tau_cld_cv, &
29                      coefw_cld_cv, &                                      ! RomP, AJ
30                      epmax_diag1)  ! epmax_cape
31! **************************************************************
32! *
33! CV_DRIVER                                                   *
34! *
35! *
36! written by   : Sandrine Bony-Lena , 17/05/2003, 11.19.41    *
37! modified by :                                               *
38! **************************************************************
39! **************************************************************
40
41  USE print_control_mod, ONLY: prt_level, lunout
42  IMPLICIT NONE
43
44! .............................START PROLOGUE............................
45
46
47! All argument names (except len,nd,ntra,nloc,delt and the flags) have a "1" appended.
48! The "1" is removed for the corresponding compressed variables.
49! PARAMETERS:
50! Name            Type         Usage            Description
51! ----------      ----------     -------  ----------------------------
52
53! len           Integer        Input        first (i) dimension
54! nd            Integer        Input        vertical (k) dimension
55! ndp1          Integer        Input        nd + 1
56! ntra          Integer        Input        number of tracors
57! nloc          Integer        Input        dimension of arrays for compressed fields
58! k_upper       Integer        Input        upmost level for vertical loops
59! iflag_con     Integer        Input        version of convect (3/4)
60! iflag_mix     Integer        Input        version of mixing  (0/1/2)
61! iflag_ice_thermo Integer        Input        accounting for ice thermodynamics (0/1)
62! iflag_clos    Integer        Input        version of closure (0/1)
63! tau_cld_cv    Real           Input        characteristic time of dissipation of mixing fluxes
64! coefw_cld_cv  Real           Input        coefficient for updraft velocity in convection
65! ok_conserv_q  Logical        Input        when true corrections for water conservation are swtiched on
66! delt          Real           Input        time step
67! comp_threshold Real           Input       threshold on the fraction of convective points below which
68!                                            fields  are compressed
69! t1            Real           Input        temperature (sat draught envt)
70! q1            Real           Input        specific hum (sat draught envt)
71! qs1           Real           Input        sat specific hum (sat draught envt)
72! t1_wake       Real           Input        temperature (unsat draught envt)
73! q1_wake       Real           Input        specific hum(unsat draught envt)
74! qs1_wake      Real           Input        sat specific hum(unsat draughts envt)
75! s1_wake       Real           Input        fractionnal area covered by wakes
76! u1            Real           Input        u-wind
77! v1            Real           Input        v-wind
78! tra1          Real           Input        tracors
79! p1            Real           Input        full level pressure
80! ph1           Real           Input        half level pressure
81! ALE1          Real           Input        Available lifting Energy
82! ALP1          Real           Input        Available lifting Power
83! sig1feed1     Real           Input        sigma coord at lower bound of feeding layer
84! sig2feed1     Real           Input        sigma coord at upper bound of feeding layer
85! wght1         Real           Input        weight density determining the feeding mixture
86! iflag1        Integer        Output       flag for Emanuel conditions
87! ft1           Real           Output       temp tend
88! fq1           Real           Output       spec hum tend
89! fu1           Real           Output       u-wind tend
90! fv1           Real           Output       v-wind tend
91! ftra1         Real           Output       tracor tend
92! precip1       Real           Output       precipitation
93! kbas1         Integer        Output       cloud base level
94! ktop1         Integer        Output       cloud top level
95! cbmf1         Real           Output       cloud base mass flux
96! sig1          Real           In/Out       section adiabatic updraft
97! w01           Real           In/Out       vertical velocity within adiab updraft
98! ptop21        Real           In/Out       top of entraining zone
99! Ma1           Real           Output       mass flux adiabatic updraft
100! mip1          Real           Output       mass flux shed by the adiabatic updraft
101! Vprecip1      Real           Output       vertical profile of total precipitation
102! Vprecipi1     Real           Output       vertical profile of ice precipitation
103! upwd1         Real           Output       total upward mass flux (adiab+mixed)
104! dnwd1         Real           Output       saturated downward mass flux (mixed)
105! dnwd01        Real           Output       unsaturated downward mass flux
106! qcondc1       Real           Output       in-cld mixing ratio of condensed water
107! wd1           Real           Output       downdraft velocity scale for sfc fluxes
108! cape1         Real           Output       CAPE
109! cin1          Real           Output       CIN
110! tvp1          Real           Output       adiab lifted parcell virt temp
111! ftd1          Real           Output       precip temp tend
112! fqt1          Real           Output       precip spec hum tend
113! Plim11        Real           Output
114! Plim21        Real           Output
115! asupmax1      Real           Output
116! supmax01      Real           Output
117! asupmaxmin1   Real           Output
118
119! ftd1          Real           Output  Array of temperature tendency due to precipitations (K/s) of dimension ND,
120!                                      defined at same grid levels as T, Q, QS and P.
121
122! fqd1          Real           Output  Array of specific humidity tendencies due to precipitations ((gm/gm)/s)
123!                                      of dimension ND, defined at same grid levels as T, Q, QS and P.
124
125! wdtrainA1     Real           Output   precipitation detrained from adiabatic draught;
126!                                         used in tracer transport (cvltr)
127! wdtrainM1     Real           Output   precipitation detrained from mixed draughts;
128!                                         used in tracer transport (cvltr)
129! da1           Real           Output     used in tracer transport (cvltr)
130! phi1          Real           Output     used in tracer transport (cvltr)
131! mp1           Real           Output     used in tracer transport (cvltr)
132! qtc1          Real           Output     specific humidity in convection
133! sigt1         Real           Output     surface fraction in adiabatic updrafts                                         
134! phi21         Real           Output     used in tracer transport (cvltr)
135                                         
136! d1a1          Real           Output     used in tracer transport (cvltr)
137! dam1          Real           Output     used in tracer transport (cvltr)
138                                         
139! epmlmMm1      Real           Output     used in tracer transport (cvltr)
140! eplaMm1       Real           Output     used in tracer transport (cvltr)
141                                         
142! evap1         Real           Output   
143! ep1           Real           Output   
144! sigij1        Real           Output     used in tracer transport (cvltr)
145! clw1          Real           Output   condensed water content of the adiabatic updraught
146! elij1         Real           Output
147! wghti1        Real           Output   final weight of the feeding layers,
148!                                         used in tracer transport (cvltr)
149
150
151! S. Bony, Mar 2002:
152! * Several modules corresponding to different physical processes
153! * Several versions of convect may be used:
154!         - iflag_con=3: version lmd  (previously named convect3)
155!         - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
156! + tard: - iflag_con=5: version lmd with ice (previously named convectg)
157! S. Bony, Oct 2002:
158! * Vectorization of convect3 (ie version lmd)
159
160! ..............................END PROLOGUE.............................
161
162
163
164! Input
165  INTEGER, INTENT (IN)                               :: len
166  INTEGER, INTENT (IN)                               :: nd
167  INTEGER, INTENT (IN)                               :: ndp1
168  INTEGER, INTENT (IN)                               :: ntra
169  INTEGER, INTENT(IN)                                :: nloc ! (nloc=len)  pour l'instant
170  INTEGER, INTENT (IN)                               :: k_upper
171  INTEGER, INTENT (IN)                               :: iflag_con
172  INTEGER, INTENT (IN)                               :: iflag_mix
173  INTEGER, INTENT (IN)                               :: iflag_ice_thermo
174  INTEGER, INTENT (IN)                               :: iflag_clos
175  LOGICAL, INTENT (IN)                               :: ok_conserv_q
176  REAL, INTENT (IN)                                  :: tau_cld_cv
177  REAL, INTENT (IN)                                  :: coefw_cld_cv
178  REAL, INTENT (IN)                                  :: delt
179  REAL, INTENT (IN)                                  :: comp_threshold
180  REAL, DIMENSION (len, nd), INTENT (IN)             :: t1
181  REAL, DIMENSION (len, nd), INTENT (IN)             :: q1
182  REAL, DIMENSION (len, nd), INTENT (IN)             :: qs1
183  REAL, DIMENSION (len, nd), INTENT (IN)             :: t1_wake
184  REAL, DIMENSION (len, nd), INTENT (IN)             :: q1_wake
185  REAL, DIMENSION (len, nd), INTENT (IN)             :: qs1_wake
186  REAL, DIMENSION (len), INTENT (IN)                 :: s1_wake
187  REAL, DIMENSION (len, nd), INTENT (IN)             :: u1
188  REAL, DIMENSION (len, nd), INTENT (IN)             :: v1
189  REAL, DIMENSION (len, nd, ntra), INTENT (IN)       :: tra1
190  REAL, DIMENSION (len, nd), INTENT (IN)             :: p1
191  REAL, DIMENSION (len, ndp1), INTENT (IN)           :: ph1
192  REAL, DIMENSION (len), INTENT (IN)                 :: Ale1
193  REAL, DIMENSION (len), INTENT (IN)                 :: Alp1
194  REAL, DIMENSION (len, nd), INTENT (IN)             :: omega1
195  REAL, INTENT (IN)                                  :: sig1feed1 ! pressure at lower bound of feeding layer
196  REAL, INTENT (IN)                                  :: sig2feed1 ! pressure at upper bound of feeding layer
197  REAL, DIMENSION (nd), INTENT (IN)                  :: wght1     ! weight density determining the feeding mixture
198  INTEGER, DIMENSION (len), INTENT (IN)              :: lalim_conv1
199
200! Input/Output
201  REAL, DIMENSION (len, nd), INTENT (INOUT)          :: sig1
202  REAL, DIMENSION (len, nd), INTENT (INOUT)          :: w01
203
204! Output
205  INTEGER, DIMENSION (len), INTENT (OUT)             :: iflag1
206  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ft1
207  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fq1
208  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fu1
209  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fv1
210  REAL, DIMENSION (len, nd, ntra), INTENT (OUT)      :: ftra1
211  REAL, DIMENSION (len), INTENT (OUT)                :: precip1
212  INTEGER, DIMENSION (len), INTENT (OUT)             :: kbas1
213  INTEGER, DIMENSION (len), INTENT (OUT)             :: ktop1
214  REAL, DIMENSION (len), INTENT (OUT)                :: cbmf1
215  REAL, DIMENSION (len), INTENT (OUT)                :: plcl1
216  REAL, DIMENSION (len), INTENT (OUT)                :: plfc1
217  REAL, DIMENSION (len), INTENT (OUT)                :: wbeff1
218  REAL, DIMENSION (len), INTENT (OUT)                :: ptop21
219  REAL, DIMENSION (len), INTENT (OUT)                :: sigd1
220  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ma1        ! adiab. asc. mass flux (staggered grid)
221  REAL, DIMENSION (len, nd), INTENT (OUT)            :: mip1       ! mass flux shed from adiab. ascent (extensive)
222! real Vprecip1(len,nd)
223  REAL, DIMENSION (len, ndp1), INTENT (OUT)          :: vprecip1   ! tot precipitation flux (staggered grid)
224  REAL, DIMENSION (len, ndp1), INTENT (OUT)          :: vprecipi1  ! ice precipitation flux (staggered grid)
225  REAL, DIMENSION (len, nd), INTENT (OUT)            :: upwd1      ! upwd sat. mass flux (staggered grid)
226  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dnwd1      ! dnwd sat. mass flux (staggered grid)
227  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dnwd01     ! unsat. mass flux (staggered grid)
228  REAL, DIMENSION (len, nd), INTENT (OUT)            :: qcondc1    ! max cloud condensate (intensive)  ! cld
229  REAL, DIMENSION (len), INTENT (OUT)                :: wd1             ! gust
230  REAL, DIMENSION (len), INTENT (OUT)                :: cape1
231  REAL, DIMENSION (len), INTENT (OUT)                :: cin1
232  REAL, DIMENSION (len, nd), INTENT (OUT)            :: tvp1       ! Virt. temp. in the adiab. ascent
233
234!AC!
235!!      real da1(len,nd),phi1(len,nd,nd)
236!!      real da(len,nd),phi(len,nd,nd)
237!AC!
238  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ftd1       ! Temp. tendency due to the sole unsat. drafts
239  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fqd1       ! Moist. tendency due to the sole unsat. drafts
240  REAL, DIMENSION (len), INTENT (OUT)                :: Plim11
241  REAL, DIMENSION (len), INTENT (OUT)                :: Plim21
242  REAL, DIMENSION (len, nd), INTENT (OUT)            :: asupmax1   ! Highest mixing fraction of mixed updraughts
243  REAL, DIMENSION (len), INTENT (OUT)                :: supmax01
244  REAL, DIMENSION (len), INTENT (OUT)                :: asupmaxmin1
245  REAL, DIMENSION (len, nd), INTENT (OUT)            :: qtc1    ! in cloud water content (intensive)   ! cld
246  REAL, DIMENSION (len, nd), INTENT (OUT)            :: sigt1   ! fract. cloud area (intensive)        ! cld
247
248! RomP >>>
249  REAL, DIMENSION (len, nd), INTENT (OUT)            :: wdtrainA1, wdtrainM1 ! precipitation sources (extensive)
250  REAL, DIMENSION (len, nd), INTENT (OUT)            :: mp1  ! unsat. mass flux (staggered grid)
251  REAL, DIMENSION (len, nd), INTENT (OUT)            :: da1  ! detrained mass flux of adiab. asc. air (extensive)
252  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: phi1 ! mass flux of envt. air in mixed draughts (extensive)
253  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: epmlmMm1  ! (extensive)
254  REAL, DIMENSION (len, nd), INTENT (OUT)            :: eplaMm1   ! (extensive)
255  REAL, DIMENSION (len, nd), INTENT (OUT)            :: evap1 ! evaporation rate in precip. downdraft. (intensive)
256  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ep1
257  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: sigij1 ! mass fraction of env. air in mixed draughts (intensive)
258  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: elij1! cond. water per unit mass of mixed draughts (intensive)
259  REAL, DIMENSION (len, nd), INTENT (OUT)            :: clw1 ! cond. water per unit mass of the adiab. asc. (intensive)
260!JYG,RL
261  REAL, DIMENSION (len, nd), INTENT (OUT)            :: wghti1   ! final weight of the feeding layers (extensive)
262!JYG,RL
263  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: phi21    ! (extensive)
264  REAL, DIMENSION (len, nd), INTENT (OUT)            :: d1a1     ! (extensive)
265  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dam1     ! (extensive)
266! RomP <<<
267  REAL, DIMENSION (len ), INTENT (OUT)               :: epmax_diag1     
268
269! -------------------------------------------------------------------
270! Prolog by Kerry Emanuel.
271! -------------------------------------------------------------------
272! --- ARGUMENTS
273! -------------------------------------------------------------------
274! --- On input:
275
276! t:   Array of absolute temperature (K) of dimension ND, with first
277! index corresponding to lowest model level. Note that this array
278! will be altered by the subroutine if dry convective adjustment
279! occurs and if IPBL is not equal to 0.
280
281! q:   Array of specific humidity (gm/gm) of dimension ND, with first
282! index corresponding to lowest model level. Must be defined
283! at same grid levels as T. Note that this array will be altered
284! if dry convective adjustment occurs and if IPBL is not equal to 0.
285
286! qs:  Array of saturation specific humidity of dimension ND, with first
287! index corresponding to lowest model level. Must be defined
288! at same grid levels as T. Note that this array will be altered
289! if dry convective adjustment occurs and if IPBL is not equal to 0.
290
291! t_wake: Array of absolute temperature (K), seen by unsaturated draughts,
292! of dimension ND, with first index corresponding to lowest model level.
293
294! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts,
295! of dimension ND, with first index corresponding to lowest model level.
296! Must be defined at same grid levels as T.
297
298! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts,
299! of dimension ND, with first index corresponding to lowest model level.
300! Must be defined at same grid levels as T.
301
302! s_wake: Array of fractionnal area occupied by the wakes.
303
304! u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
305! index corresponding with the lowest model level. Defined at
306! same levels as T. Note that this array will be altered if
307! dry convective adjustment occurs and if IPBL is not equal to 0.
308
309! v:   Same as u but for meridional velocity.
310
311! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
312! where NTRA is the number of different tracers. If no
313! convective tracer transport is needed, define a dummy
314! input array of dimension (ND,1). Tracers are defined at
315! same vertical levels as T. Note that this array will be altered
316! if dry convective adjustment occurs and if IPBL is not equal to 0.
317
318! p:   Array of pressure (mb) of dimension ND, with first
319! index corresponding to lowest model level. Must be defined
320! at same grid levels as T.
321
322! ph:  Array of pressure (mb) of dimension ND+1, with first index
323! corresponding to lowest level. These pressures are defined at
324! levels intermediate between those of P, T, Q and QS. The first
325! value of PH should be greater than (i.e. at a lower level than)
326! the first value of the array P.
327
328! ALE:  Available lifting Energy
329
330! ALP:  Available lifting Power
331
332! nl:  The maximum number of levels to which convection can penetrate, plus 1.
333!       NL MUST be less than or equal to ND-1.
334
335! delt: The model time step (sec) between calls to CONVECT
336
337! ----------------------------------------------------------------------------
338! ---   On Output:
339
340! iflag: An output integer whose value denotes the following:
341!       VALUE   INTERPRETATION
342!       -----   --------------
343!         0     Moist convection occurs.
344!         1     Moist convection occurs, but a CFL condition
345!               on the subsidence warming is violated. This
346!               does not cause the scheme to terminate.
347!         2     Moist convection, but no precip because ep(inb) lt 0.0001
348!         3     No moist convection because new cbmf is 0 and old cbmf is 0.
349!         4     No moist convection; atmosphere is not
350!               unstable
351!         6     No moist convection because ihmin le minorig.
352!         7     No moist convection because unreasonable
353!               parcel level temperature or specific humidity.
354!         8     No moist convection: lifted condensation
355!               level is above the 200 mb level.
356!         9     No moist convection: cloud base is higher
357!               then the level NL-1.
358!        10     No moist convection: cloud top is too warm.
359!
360
361! ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
362!       grid levels as T, Q, QS and P.
363
364! fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
365!       defined at same grid levels as T, Q, QS and P.
366
367! fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
368!      defined at same grid levels as T.
369
370! fv:   Same as FU, but for forcing of meridional velocity.
371
372! ftra: Array of forcing of tracer content, in tracer mixing ratio per
373!       second, defined at same levels as T. Dimensioned (ND,NTRA).
374
375! precip: Scalar convective precipitation rate (mm/day).
376
377! wd:   A convective downdraft velocity scale. For use in surface
378!       flux parameterizations. See convect.ps file for details.
379
380! tprime: A convective downdraft temperature perturbation scale (K).
381!         For use in surface flux parameterizations. See convect.ps
382!         file for details.
383
384! qprime: A convective downdraft specific humidity
385!         perturbation scale (gm/gm).
386!         For use in surface flux parameterizations. See convect.ps
387!         file for details.
388
389! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
390!       BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
391!       ITS NEXT CALL. That is, the value of CBMF must be "remembered"
392!       by the calling program between calls to CONVECT.
393
394! det:   Array of detrainment mass flux of dimension ND.
395! -------------------------------------------------------------------
396
397! Local (non compressed) arrays
398
399
400  INTEGER i, k, il
401  INTEGER nword1, nword2, nword3, nword4
402  INTEGER icbmax
403  INTEGER nk1(len)
404  INTEGER icb1(len)
405  INTEGER icbs1(len)
406
407  LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
408  LOGICAL, SAVE :: debut = .TRUE.
409!$OMP THREADPRIVATE(debut)
410
411  REAL coef_convective(len)   ! = 1 for convective points, = 0 otherwise
412  REAL tnk1(len)
413  REAL thnk1(len)
414  REAL qnk1(len)
415  REAL gznk1(len)
416  REAL qsnk1(len)
417  REAL unk1(len)
418  REAL vnk1(len)
419  REAL cpnk1(len)
420  REAL hnk1(len)
421  REAL pbase1(len)
422  REAL buoybase1(len)
423
424  REAL lf1(len, nd), lf1_wake(len, nd)
425  REAL lv1(len, nd), lv1_wake(len, nd)
426  REAL cpn1(len, nd), cpn1_wake(len, nd)
427  REAL tv1(len, nd), tv1_wake(len, nd)
428  REAL gz1(len, nd), gz1_wake(len, nd)
429  REAL hm1(len, nd)
430  REAL h1(len, nd), h1_wake(len, nd)
431  REAL tp1(len, nd)
432  REAL th1(len, nd), th1_wake(len, nd)
433
434  REAL bid(len, nd) ! dummy array
435
436  INTEGER ncum
437
438  REAL p1feed1(len) ! pressure at lower bound of feeding layer
439  REAL p2feed1(len) ! pressure at upper bound of feeding layer
440!JYG,RL
441!!      real wghti1(len,nd) ! weights of the feeding layers
442!JYG,RL
443
444! (local) compressed fields:
445
446
447  INTEGER idcum(nloc)
448!jyg<
449  LOGICAL compress    ! True if compression occurs
450!>jyg
451  INTEGER iflag(nloc), nk(nloc), icb(nloc)
452  INTEGER nent(nloc, nd)
453  INTEGER icbs(nloc)
454  INTEGER inb(nloc), inbis(nloc)
455
456  REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc)
457  REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd)
458  REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd)
459  REAL s_wake(nloc)
460  REAL u(nloc, nd), v(nloc, nd)
461  REAL gz(nloc, nd), h(nloc, nd)
462  REAL h_wake(nloc, nd)
463  REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd)
464  REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd)
465  REAL p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd)
466  REAL tv_wake(nloc, nd)
467  REAL clw(nloc, nd)
468  REAL dph(nloc, nd)
469  REAL pbase(nloc), buoybase(nloc), th(nloc, nd)
470  REAL th_wake(nloc, nd)
471  REAL tvp(nloc, nd)
472  REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc)
473  REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd)
474  REAL buoy(nloc, nd)
475  REAL cape(nloc)
476  REAL cin(nloc)
477  REAL m(nloc, nd)
478  REAL ment(nloc, nd, nd), sigij(nloc, nd, nd)
479  REAL qent(nloc, nd, nd)
480  REAL hent(nloc, nd, nd)
481  REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
482  REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
483  REAL elij(nloc, nd, nd)
484  REAL supmax(nloc, nd)
485  REAL Ale(nloc), Alp(nloc), coef_clos(nloc)
486  REAL omega(nloc,nd)
487  REAL sigd(nloc)
488! real mp(nloc,nd), qp(nloc,nd), up(nloc,nd), vp(nloc,nd)
489! real wt(nloc,nd), water(nloc,nd), evap(nloc,nd), ice(nloc,nd)
490! real b(nloc,nd), sigd(nloc)
491! save mp,qp,up,vp,wt,water,evap,b
492  REAL, DIMENSION(len,nd)     :: mp, qp, up, vp
493  REAL, DIMENSION(len,nd)     :: wt, water, evap
494  REAL, DIMENSION(len,nd)     :: ice, fondue, b
495  REAL, DIMENSION(len,nd)     :: frac, faci
496  REAL ft(nloc, nd), fq(nloc, nd)
497  REAL ftd(nloc, nd), fqd(nloc, nd)
498  REAL fu(nloc, nd), fv(nloc, nd)
499  REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd)
500  REAL ma(nloc, nd), mip(nloc, nd)
501!!  REAL tls(nloc, nd), tps(nloc, nd)                 ! unused . jyg
502  REAL qprime(nloc), tprime(nloc)
503  REAL precip(nloc)
504! real Vprecip(nloc,nd)
505  REAL vprecip(nloc, nd+1)
506  REAL vprecipi(nloc, nd+1)
507  REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
508  REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
509  REAL qcondc(nloc, nd)      ! cld
510  REAL wd(nloc)                ! gust
511  REAL Plim1(nloc), plim2(nloc)
512  REAL asupmax(nloc, nd)
513  REAL supmax0(nloc)
514  REAL asupmaxmin(nloc)
515
516  REAL tnk(nloc), qnk(nloc), gznk(nloc)
517  REAL wghti(nloc, nd)
518  REAL hnk(nloc), unk(nloc), vnk(nloc)
519
520  REAL qtc(nloc, nd)         ! cld
521  REAL sigt(nloc, nd)        ! cld
522 
523! RomP >>>
524  REAL wdtrainA(nloc, nd), wdtrainM(nloc, nd)
525  REAL da(len, nd), phi(len, nd, nd)
526  REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd)
527  REAL phi2(len, nd, nd)
528  REAL d1a(len, nd), dam(len, nd)
529! RomP <<<
530  REAL epmax_diag(nloc) ! epmax_cape
531
532  CHARACTER (LEN=20) :: modname = 'cva_driver'
533  CHARACTER (LEN=80) :: abort_message
534
535  INTEGER,SAVE                                       :: igout=1
536!$OMP THREADPRIVATE(igout)
537
538
539! print *, 't1, t1_wake ',(k,t1(1,k),t1_wake(1,k),k=1,nd)
540! print *, 'q1, q1_wake ',(k,q1(1,k),q1_wake(1,k),k=1,nd)
541
542! -------------------------------------------------------------------
543! --- SET CONSTANTS AND PARAMETERS
544! -------------------------------------------------------------------
545
546! -- set simulation flags:
547! (common cvflag)
548
549  CALL cv_flag(iflag_ice_thermo)
550
551! -- set thermodynamical constants:
552! (common cvthermo)
553
554  CALL cv_thermo(iflag_con)
555
556! -- set convect parameters
557
558! includes microphysical parameters and parameters that
559! control the rate of approach to quasi-equilibrium)
560! (common cvparam)
561
562  IF (iflag_con==3) THEN
563    CALL cv3_param(nd, k_upper, delt)
564
565  END IF
566
567  IF (iflag_con==4) THEN
568    CALL cv_param(nd)
569  END IF
570
571! ---------------------------------------------------------------------
572! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
573! ---------------------------------------------------------------------
574  nword1 = len
575  nword2 = len*nd
576  nword3 = len*nd*ntra
577  nword4 = len*nd*nd
578
579  iflag1(:) = 0
580  ktop1(:) = 0
581  kbas1(:) = 0
582  ft1(:, :) = 0.0
583  fq1(:, :) = 0.0
584  fu1(:, :) = 0.0
585  fv1(:, :) = 0.0
586  ftra1(:, :, :) = 0.
587  precip1(:) = 0.
588  cbmf1(:) = 0.
589  plcl1(:) = 0.
590  plfc1(:) = 0.
591  wbeff1(:) = 0.
592  ptop21(:) = 0.
593  sigd1(:) = 0.
594  ma1(:, :) = 0.
595  mip1(:, :) = 0.
596  vprecip1(:, :) = 0.
597  vprecipi1(:, :) = 0.
598  upwd1(:, :) = 0.
599  dnwd1(:, :) = 0.
600  dnwd01(:, :) = 0.
601  qcondc1(:, :) = 0.
602  wd1(:) = 0.
603  cape1(:) = 0.
604  cin1(:) = 0.
605  tvp1(:, :) = 0.
606  ftd1(:, :) = 0.
607  fqd1(:, :) = 0.
608  Plim11(:) = 0.
609  Plim21(:) = 0.
610  asupmax1(:, :) = 0.
611  supmax01(:) = 0.
612  asupmaxmin1(:) = 0.
613
614  DO il = 1, len
615    cin1(il) = -100000.
616    cape1(il) = -1.
617  END DO
618
619!!  IF (iflag_con==3) THEN
620!!    DO il = 1, len
621!!      sig1(il, nd) = sig1(il, nd) + 1.
622!!      sig1(il, nd) = amin1(sig1(il,nd), 12.1)
623!!    END DO
624!!  END IF
625
626  IF (iflag_con==3) THEN
627      CALL cv3_incrcount(len,nd,delt,sig1)
628  END IF  ! (iflag_con==3)
629
630! RomP >>>
631  sigt1(:, :) = 0.
632  qtc1(:, :) = 0.
633  wdtrainA1(:, :) = 0.
634  wdtrainM1(:, :) = 0.
635  da1(:, :) = 0.
636  phi1(:, :, :) = 0.
637  epmlmMm1(:, :, :) = 0.
638  eplaMm1(:, :) = 0.
639  mp1(:, :) = 0.
640  evap1(:, :) = 0.
641  ep1(:, :) = 0.
642  sigij1(:, :, :) = 0.
643  elij1(:, :, :) = 0.
644  clw1(:,:) = 0.
645  wghti1(:,:) = 0.
646  phi21(:, :, :) = 0.
647  d1a1(:, :) = 0.
648  dam1(:, :) = 0.
649! RomP <<<
650! ---------------------------------------------------------------------
651! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS
652! ---------------------------------------------------------------------
653
654  DO il = 1, nloc
655    coef_clos(il) = 1.
656  END DO
657
658! --------------------------------------------------------------------
659! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
660! --------------------------------------------------------------------
661
662  IF (iflag_con==3) THEN
663
664    IF (debut) THEN
665      PRINT *, 'Emanuel version 3 nouvelle'
666    END IF
667! print*,'t1, q1 ',t1,q1
668        if (prt_level >= 9) &
669             PRINT *, 'cva_driver -> cv3_prelim'
670    CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, &           ! nd->na
671                    lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1)
672
673
674        if (prt_level >= 9) &
675             PRINT *, 'cva_driver -> cv3_prelim'
676    CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, & ! nd->na
677                    lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, &
678                    h1_wake, bid, th1_wake)
679
680  END IF
681
682  IF (iflag_con==4) THEN
683    PRINT *, 'Emanuel version 4 '
684        if (prt_level >= 9) &
685             PRINT *, 'cva_driver -> cv_prelim'
686    CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1, &
687                   lv1, cpn1, tv1, gz1, h1, hm1)
688  END IF
689
690! --------------------------------------------------------------------
691! --- CONVECTIVE FEED
692! --------------------------------------------------------------------
693
694! compute feeding layer potential temperature and mixing ratio :
695
696! get bounds of feeding layer
697
698! test niveaux couche alimentation KE
699  IF (sig1feed1==sig2feed1) THEN
700    WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed'
701    WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def'
702    abort_message = ''
703    CALL abort_physic(modname, abort_message, 1)
704  END IF
705
706  DO i = 1, len
707    p1feed1(i) = sig1feed1*ph1(i, 1)
708    p2feed1(i) = sig2feed1*ph1(i, 1)
709!test maf
710!   p1feed1(i)=ph1(i,1)
711!   p2feed1(i)=ph1(i,2)
712!   p2feed1(i)=ph1(i,3)
713!testCR: on prend la couche alim des thermiques
714!   p2feed1(i)=ph1(i,lalim_conv1(i)+1)
715!   print*,'lentr=',lentr(i),ph1(i,lentr(i)+1),ph1(i,2)
716  END DO
717
718  IF (iflag_con==3) THEN
719  END IF
720  DO i = 1, len
721! print*,'avant cv3_feed Plim',p1feed1(i),p2feed1(i)
722  END DO
723  IF (iflag_con==3) THEN
724
725! print*, 'IFLAG1 avant cv3_feed'
726! print*,'len,nd',len,nd
727! write(*,'(64i1)') iflag1(2:len-1)
728
729        if (prt_level >= 9) &
730             PRINT *, 'cva_driver -> cv3_feed'
731    CALL cv3_feed(len, nd, ok_conserv_q, &                 ! nd->na
732                  t1, q1, u1, v1, p1, ph1, hm1, gz1, &
733                  p1feed1, p2feed1, wght1, &
734                  wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, &
735                  cpnk1, hnk1, nk1, icb1, icbmax, iflag1, gznk1, plcl1)
736  END IF
737
738! print*, 'IFLAG1 apres cv3_feed'
739! print*,'len,nd',len,nd
740! write(*,'(64i1)') iflag1(2:len-1)
741
742  IF (iflag_con==4) THEN
743        if (prt_level >= 9) &
744             PRINT *, 'cva_driver -> cv_feed'
745    CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1, &
746                 nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
747  END IF
748
749! print *, 'cv3_feed-> iflag1, plcl1 ',iflag1(1),plcl1(1)
750
751! --------------------------------------------------------------------
752! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
753! (up through ICB for convect4, up through ICB+1 for convect3)
754! Calculates the lifted parcel virtual temperature at nk, the
755! actual temperature, and the adiabatic liquid water content.
756! --------------------------------------------------------------------
757
758  IF (iflag_con==3) THEN
759
760        if (prt_level >= 9) &
761             PRINT *, 'cva_driver -> cv3_undilute1'
762    CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, icb1, tnk1, qnk1, & ! nd->na
763                       gznk1, tp1, tvp1, clw1, icbs1)
764  END IF
765
766
767  IF (iflag_con==4) THEN
768        if (prt_level >= 9) &
769             PRINT *, 'cva_driver -> cv_undilute1'
770    CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, &
771                      tp1, tvp1, clw1)
772  END IF
773
774! -------------------------------------------------------------------
775! --- TRIGGERING
776! -------------------------------------------------------------------
777
778! print *,' avant triggering, iflag_con ',iflag_con
779
780  IF (iflag_con==3) THEN
781
782        if (prt_level >= 9) &
783             PRINT *, 'cva_driver -> cv3_trigger'
784    CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na
785                      pbase1, buoybase1, iflag1, sig1, w01)
786
787
788! print*, 'IFLAG1 apres cv3_triger'
789! print*,'len,nd',len,nd
790! write(*,'(64i1)') iflag1(2:len-1)
791
792! call dump2d(iim,jjm-1,sig1(2)
793  END IF
794
795  IF (iflag_con==4) THEN
796        if (prt_level >= 9) &
797             PRINT *, 'cva_driver -> cv_trigger'
798    CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
799  END IF
800
801
802! =====================================================================
803! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
804! =====================================================================
805
806!  Determine the number "ncum" of convective gridpoints, the list "idcum" of convective
807!  gridpoints and the weights "coef_convective" (= 1. for convective gridpoints and = 0.
808!  elsewhere).
809  ncum = 0
810  coef_convective(:) = 0.
811  DO i = 1, len
812    IF (iflag1(i)==0) THEN
813      coef_convective(i) = 1.
814      ncum = ncum + 1
815      idcum(ncum) = i
816    END IF
817  END DO
818
819! print*,'len, ncum = ',len,ncum
820
821  IF (ncum>0) THEN
822
823! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
824! --- COMPRESS THE FIELDS
825!       (-> vectorization over convective gridpoints)
826! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
827
828    IF (iflag_con==3) THEN
829! print*,'ncum tv1 ',ncum,tv1
830! print*,'tvp1 ',tvp1
831!jyg<
832!   If the fraction of convective points is larger than comp_threshold, then compression
833!   is assumed useless.
834!
835  compress = ncum .lt. len*comp_threshold
836!
837  IF (.not. compress) THEN
838    DO i = 1,len
839      idcum(i) = i
840    ENDDO
841  ENDIF
842!
843!>jyg
844        if (prt_level >= 9) &
845             PRINT *, 'cva_driver -> cv3a_compress'
846      CALL cv3a_compress(len, nloc, ncum, nd, ntra, compress, &
847                         iflag1, nk1, icb1, icbs1, &
848                         plcl1, tnk1, qnk1, gznk1, hnk1, unk1, vnk1, &
849                         wghti1, pbase1, buoybase1, &
850                         t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, &
851                         u1, v1, gz1, th1, th1_wake, &
852                         tra1, &
853                         h1, lv1, lf1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
854                         h1_wake, lv1_wake, lf1_wake, cpn1_wake, tv1_wake, &
855                         sig1, w01, ptop21, &
856                         Ale1, Alp1, omega1, &
857                         iflag, nk, icb, icbs, &
858                         plcl, tnk, qnk, gznk, hnk, unk, vnk, &
859                         wghti, pbase, buoybase, &
860                         t, q, qs, t_wake, q_wake, qs_wake, s_wake, &
861                         u, v, gz, th, th_wake, &
862                         tra, &
863                         h, lv, lf, cpn, p, ph, tv, tp, tvp, clw, &
864                         h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, &
865                         sig, w0, ptop2, &
866                         Ale, Alp, omega)
867
868! print*,'tv ',tv
869! print*,'tvp ',tvp
870
871    END IF
872
873    IF (iflag_con==4) THEN
874        if (prt_level >= 9) &
875             PRINT *, 'cva_driver -> cv_compress'
876      CALL cv_compress(len, nloc, ncum, nd, &
877                       iflag1, nk1, icb1, &
878                       cbmf1, plcl1, tnk1, qnk1, gznk1, &
879                       t1, q1, qs1, u1, v1, gz1, &
880                       h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
881                       iflag, nk, icb, &
882                       cbmf, plcl, tnk, qnk, gznk, &
883                       t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw, &
884                       dph)
885    END IF
886
887! -------------------------------------------------------------------
888! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
889! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
890! ---   &
891! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
892! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
893! ---   &
894! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
895! -------------------------------------------------------------------
896
897    IF (iflag_con==3) THEN
898        if (prt_level >= 9) &
899             PRINT *, 'cva_driver -> cv3_undilute2'
900      CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, &        !na->nd
901                         tnk, qnk, gznk, hnk, t, q, qs, gz, &
902                         p, ph, h, tv, lv, lf, pbase, buoybase, plcl, &
903                         inb, tp, tvp, clw, hp, ep, sigp, buoy, &
904                         frac)
905    END IF
906
907    IF (iflag_con==4) THEN
908        if (prt_level >= 9) &
909             PRINT *, 'cva_driver -> cv_undilute2'
910      CALL cv_undilute2(nloc, ncum, nd, icb, nk, &
911                        tnk, qnk, gznk, t, q, qs, gz, &
912                        p, dph, h, tv, lv, &
913                        inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
914    END IF
915
916    ! epmax_cape
917    ! on recalcule ep et hp   
918        if (prt_level >= 9) &
919             PRINT *, 'cva_driver -> cv3_epmax_cape'
920    call cv3_epmax_fn_cape(nloc,ncum,nd &
921                , ep,hp,icb,inb,clw,nk,t,h,hnk,lv,lf,frac &
922                , pbase, p, ph, tv, buoy, sig, w0,iflag &
923                , epmax_diag)
924
925! -------------------------------------------------------------------
926! --- MIXING(1)   (if iflag_mix .ge. 1)
927! -------------------------------------------------------------------
928    IF (iflag_con==3) THEN
929!      IF ((iflag_ice_thermo==1) .AND. (iflag_mix/=0)) THEN
930!        WRITE (*, *) ' iflag_ice_thermo==1 requires iflag_mix==0', ' but iflag_mix=', iflag_mix, &
931!          '. Might as well stop here.'
932!        STOP
933!      END IF
934      IF (iflag_mix>=1) THEN
935        CALL zilch(supmax, nloc*nd)
936        if (prt_level >= 9) &
937             PRINT *, 'cva_driver -> cv3p_mixing'
938        CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &           ! na->nd
939                         ph, t, q, qs, u, v, tra, h, lv, lf, frac, qnk, &
940                         unk, vnk, hp, tv, tvp, ep, clw, sig, &
941                         ment, qent, hent, uent, vent, nent, &
942                         sigij, elij, supmax, ments, qents, traent)
943! print*, 'cv3p_mixing-> supmax ', (supmax(1,k), k=1,nd)
944
945      ELSE
946        CALL zilch(supmax, nloc*nd)
947      END IF
948    END IF
949! -------------------------------------------------------------------
950! --- CLOSURE
951! -------------------------------------------------------------------
952
953
954    IF (iflag_con==3) THEN
955      IF (iflag_clos==0) THEN
956        if (prt_level >= 9) &
957             PRINT *, 'cva_driver -> cv3_closure'
958        CALL cv3_closure(nloc, ncum, nd, icb, inb, &           ! na->nd
959                         pbase, p, ph, tv, buoy, &
960                         sig, w0, cape, m, iflag)
961      END IF   ! iflag_clos==0
962
963      ok_inhib = iflag_mix == 2
964
965      IF (iflag_clos==1) THEN
966        PRINT *, ' pas d appel cv3p_closure'
967! c        CALL cv3p_closure(nloc,ncum,nd,icb,inb              ! na->nd
968! c    :                       ,pbase,plcl,p,ph,tv,tvp,buoy
969! c    :                       ,supmax
970! c    o                       ,sig,w0,ptop2,cape,cin,m)
971      END IF   ! iflag_clos==1
972
973      IF (iflag_clos==2) THEN
974        if (prt_level >= 9) &
975             PRINT *, 'cva_driver -> cv3p1_closure'
976        CALL cv3p1_closure(nloc, ncum, nd, icb, inb, &         ! na->nd
977                           pbase, plcl, p, ph, tv, tvp, buoy, &
978                           supmax, ok_inhib, Ale, Alp, omega, &
979                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
980                           Plim1, plim2, asupmax, supmax0, &
981                           asupmaxmin, cbmf, plfc, wbeff)
982        if (prt_level >= 10) &
983             PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1)
984      END IF   ! iflag_clos==2
985
986      IF (iflag_clos==3) THEN
987        if (prt_level >= 9) &
988             PRINT *, 'cva_driver -> cv3p2_closure'
989        CALL cv3p2_closure(nloc, ncum, nd, icb, inb, &         ! na->nd
990                           pbase, plcl, p, ph, tv, tvp, buoy, &
991                           supmax, ok_inhib, Ale, Alp, omega, &
992                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
993                           Plim1, plim2, asupmax, supmax0, &
994                           asupmaxmin, cbmf, plfc, wbeff)
995        if (prt_level >= 10) &
996             PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1)
997      END IF   ! iflag_clos==3
998    END IF ! iflag_con==3
999
1000    IF (iflag_con==4) THEN
1001        if (prt_level >= 9) &
1002             PRINT *, 'cva_driver -> cv_closure'
1003      CALL cv_closure(nloc, ncum, nd, nk, icb, &
1004                         tv, tvp, p, ph, dph, plcl, cpn, &
1005                         iflag, cbmf)
1006    END IF
1007
1008! print *,'cv_closure-> cape ',cape(1)
1009
1010! -------------------------------------------------------------------
1011! --- MIXING(2)
1012! -------------------------------------------------------------------
1013
1014    IF (iflag_con==3) THEN
1015      IF (iflag_mix==0) THEN
1016        if (prt_level >= 9) &
1017             PRINT *, 'cva_driver -> cv3_mixing'
1018        CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &             ! na->nd
1019                        ph, t, q, qs, u, v, tra, h, lv, lf, frac, qnk, &
1020                        unk, vnk, hp, tv, tvp, ep, clw, m, sig, &
1021                        ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent)
1022        CALL zilch(hent, nloc*nd*nd)
1023      ELSE
1024        CALL cv3_mixscale(nloc, ncum, nd, ment, m)
1025        IF (debut) THEN
1026          PRINT *, ' cv3_mixscale-> '
1027        END IF !(debut) THEN
1028      END IF
1029    END IF
1030
1031    IF (iflag_con==4) THEN
1032        if (prt_level >= 9) &
1033             PRINT *, 'cva_driver -> cv_mixing'
1034      CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis, &
1035                     ph, t, q, qs, u, v, h, lv, qnk, &
1036                     hp, tv, tvp, ep, clw, cbmf, &
1037                     m, ment, qent, uent, vent, nent, sigij, elij)
1038    END IF                                                                                         
1039
1040    IF (debut) THEN
1041      PRINT *, ' cv_mixing ->'
1042    END IF !(debut) THEN
1043! do i = 1,nd
1044! print*,'cv_mixing-> i,ment ',i,(ment(1,i,j),j=1,nd)
1045! enddo
1046
1047! -------------------------------------------------------------------
1048! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
1049! -------------------------------------------------------------------
1050    IF (iflag_con==3) THEN
1051      IF (debut) THEN
1052        PRINT *, ' cva_driver -> cv3_unsat '
1053      END IF !(debut) THEN
1054
1055        if (prt_level >= 9) &
1056             PRINT *, 'cva_driver -> cv3_unsat'
1057      CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, &              ! na->nd
1058                     t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, &
1059                     th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, &
1060                     ep, sigp, clw, &
1061                     m, ment, elij, delt, plcl, coef_clos, &
1062                     mp, qp, up, vp, trap, wt, water, evap, fondue, ice, &
1063                     faci, b, sigd, &
1064                     wdtrainA, wdtrainM)                                       ! RomP
1065!
1066      IF (prt_level >= 10) THEN
1067        Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue '
1068        DO k = 1,nd
1069        write (6, '(i4,5(1x,e13.6))'), &
1070          k, mp(igout,k), water(igout,k), ice(igout,k), &
1071           evap(igout,k), fondue(igout,k)
1072        ENDDO
1073        Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainM '
1074        DO k = 1,nd
1075        write (6, '(i4,2(1x,e13.6))'), &
1076           k, wdtrainA(igout,k), wdtrainM(igout,k)
1077        ENDDO
1078      ENDIF
1079!
1080    END IF  !(iflag_con==3)
1081
1082    IF (iflag_con==4) THEN
1083        if (prt_level >= 9) &
1084             PRINT *, 'cva_driver -> cv_unsat'
1085      CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph, &
1086                     h, lv, ep, sigp, clw, m, ment, elij, &
1087                     iflag, mp, qp, up, vp, wt, water, evap)
1088    END IF
1089
1090    IF (debut) THEN
1091      PRINT *, 'cv_unsat-> '
1092    END IF !(debut) THEN
1093
1094! print *,'cv_unsat-> mp ',mp
1095! print *,'cv_unsat-> water ',water
1096! -------------------------------------------------------------------
1097! --- YIELD
1098! (tendencies, precipitation, variables of interface with other
1099! processes, etc)
1100! -------------------------------------------------------------------
1101
1102    IF (iflag_con==3) THEN
1103
1104        if (prt_level >= 9) &
1105             PRINT *, 'cva_driver -> cv3_yield'
1106      CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, &                      ! na->nd
1107                     icb, inb, delt, &
1108                     t, q, t_wake, q_wake, s_wake, u, v, tra, &
1109                     gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, &
1110                     ep, clw, m, tp, mp, qp, up, vp, trap, &
1111                     wt, water, ice, evap, fondue, faci, b, sigd, &
1112                     ment, qent, hent, iflag_mix, uent, vent, &
1113                     nent, elij, traent, sig, &
1114                     tv, tvp, wghti, &
1115                     iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, &      ! jyg
1116                     cbmf, upwd, dnwd, dnwd0, ma, mip, &
1117!!                     tls, tps, &                            ! useless . jyg
1118                     qcondc, wd, &
1119                     ftd, fqd, qnk, qtc, sigt, tau_cld_cv, coefw_cld_cv)
1120!
1121      IF (debut) THEN
1122        PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1)
1123      END IF !(debut) THEN
1124!   
1125      IF (prt_level >= 10) THEN
1126        Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', &
1127                    ft(igout,1), ftd(igout,1)
1128        Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', &
1129                    fq(igout,1), fqd(igout,1)
1130      ENDIF
1131!   
1132    END IF
1133
1134    IF (iflag_con==4) THEN
1135        if (prt_level >= 9) &
1136             PRINT *, 'cva_driver -> cv_yield'
1137      CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt, &
1138                     t, q, u, v, &
1139                     gz, p, ph, h, hp, lv, cpn, &
1140                     ep, clw, frac, m, mp, qp, up, vp, &
1141                     wt, water, evap, &
1142                     ment, qent, uent, vent, nent, elij, &
1143                     tv, tvp, &
1144                     iflag, wd, qprime, tprime, &
1145                     precip, cbmf, ft, fq, fu, fv, ma, qcondc)
1146    END IF
1147
1148!AC!
1149!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1150!--- passive tracers
1151!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1152
1153    IF (iflag_con==3) THEN
1154!RomP >>>
1155        if (prt_level >= 9) &
1156             PRINT *, 'cva_driver -> cv3_tracer'
1157      CALL cv3_tracer(nloc, len, ncum, nd, nd, &
1158                     ment, sigij, da, phi, phi2, d1a, dam, &
1159                     ep, vprecip, elij, clw, epmlmMm, eplaMm, &
1160                     icb, inb)
1161!RomP <<<
1162    END IF
1163
1164!AC!
1165
1166! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1167! --- UNCOMPRESS THE FIELDS
1168! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1169
1170
1171    IF (iflag_con==3) THEN
1172        if (prt_level >= 9) &
1173             PRINT *, 'cva_driver -> cv3a_uncompress'
1174      CALL cv3a_uncompress(nloc, len, ncum, nd, ntra, idcum, compress, &
1175                           iflag, icb, inb, &
1176                           precip, cbmf, plcl, plfc, wbeff, sig, w0, ptop2, &
1177                           ft, fq, fu, fv, ftra, &
1178                           sigd, ma, mip, vprecip, vprecipi, upwd, dnwd, dnwd0, &
1179                           qcondc, wd, cape, cin, &
1180                           tvp, &
1181                           ftd, fqd, &
1182                           Plim1, plim2, asupmax, supmax0, &
1183                           asupmaxmin, &
1184                           da, phi, mp, phi2, d1a, dam, sigij, &         ! RomP
1185                           clw, elij, evap, ep, epmlmMm, eplaMm, &       ! RomP
1186                           wdtrainA, wdtrainM, &                         ! RomP
1187                           qtc, sigt, epmax_diag, & ! epmax_cape
1188                           iflag1, kbas1, ktop1, &
1189                           precip1, cbmf1, plcl1, plfc1, wbeff1, sig1, w01, ptop21, &
1190                           ft1, fq1, fu1, fv1, ftra1, &
1191                           sigd1, ma1, mip1, vprecip1, vprecipi1, upwd1, dnwd1, dnwd01, &
1192                           qcondc1, wd1, cape1, cin1, &
1193                           tvp1, &
1194                           ftd1, fqd1, &
1195                           Plim11, plim21, asupmax1, supmax01, &
1196                           asupmaxmin1, &
1197                           da1, phi1, mp1, phi21, d1a1, dam1, sigij1,  & ! RomP
1198                           clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & ! RomP
1199                           wdtrainA1, wdtrainM1,                       & ! RomP
1200                           qtc1, sigt1, epmax_diag1) ! epmax_cape
1201    END IF
1202
1203    IF (iflag_con==4) THEN
1204        if (prt_level >= 9) &
1205             PRINT *, 'cva_driver -> cv_uncompress'
1206      CALL cv_uncompress(nloc, len, ncum, nd, idcum, &
1207                           iflag, &
1208                           precip, cbmf, &
1209                           ft, fq, fu, fv, &
1210                           ma, qcondc, &
1211                           iflag1, &
1212                           precip1,cbmf1, &
1213                           ft1, fq1, fu1, fv1, &
1214                           ma1, qcondc1)
1215    END IF
1216
1217  END IF ! ncum>0
1218
1219!
1220! In order take into account the possibility of changing the compression,
1221! reset m, sig and w0 to zero for non-convective points.
1222  DO k = 1,nd-1
1223        sig1(:, k) = sig1(:, k)*coef_convective(:)
1224        w01(:, k)  = w01(:, k)*coef_convective(:)
1225  ENDDO
1226
1227  IF (debut) THEN
1228    PRINT *, ' cv_uncompress -> '
1229    debut = .FALSE.
1230  END IF  !(debut) THEN
1231
1232
1233  RETURN
1234END SUBROUTINE cva_driver
Note: See TracBrowser for help on using the repository browser.