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

Last change on this file since 5756 was 5714, checked in by yann meurdesoif, 5 months ago

Convection GPU porting : Repair broken isotopes compilation : work a
round : modif introducing coef_clos_eff has not been reported in isotope version

YM

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