source: LMDZ6/branches/Amaury_dev/libf/phylmdiso/cva_driver.F90 @ 5116

Last change on this file since 5116 was 5116, checked in by abarral, 4 months ago

rename modules properly lmdz_*
move ismin, ismax, minmax into new lmdz_libmath.f90
(lint) uppercase fortran keywords

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