source: LMDZ6/branches/Optimisation_LMDZ/libf/phylmd/cva_driver.F90 @ 3757

Last change on this file since 3757 was 3757, checked in by adurocher, 4 years ago

Create separate subroutines for cv3a_driver and cv4a_driver

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