source: LMDZ6/trunk/libf/phylmd/cva_driver.f90 @ 5491

Last change on this file since 5491 was 5491, checked in by jyg, 5 hours ago

New outputs :

+ coef_clos = [conv mass flux given by Alp closure]/[conv mass flux given by Emanuel scheme closure]
+ coef_clos_eff = effective coefficient used in the convective scheme.

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