source: LMDZ6/trunk/libf/phylmdiso/cva_driver.F90 @ 5322

Last change on this file since 5322 was 5279, checked in by abarral, 3 weeks ago

fix: comma before i/o item

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