source: LMDZ6/branches/LMDZ_cdrag_LSCE/libf/phylmdiso/cva_driver.F90 @ 5456

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