source: LMDZ6/branches/cirrus/libf/phylmd/cva_driver.F90 @ 5329

Last change on this file since 5329 was 4613, checked in by fhourdin, 17 months ago

Integration/replay_isation travail Louis (ratqs)

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