source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phyiso/cva_driver.F90 @ 3363

Last change on this file since 3363 was 3331, checked in by acozic, 6 years ago

Add modification for isotopes

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