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

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

Use cvl_comp_threshold to control whether to use compression or not

cvl_comp_threshold now defaults to 0 (never compress)

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