source: LMDZ6/branches/Amaury_dev/libf/phylmd/cva_driver.F90 @ 5151

Last change on this file since 5151 was 5142, checked in by abarral, 8 weeks ago

Put cvparam.h, fcg_gcssold.h, planete.h, tsoilnudge.h, YOECUMF.h into modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 52.3 KB
Line 
1
2! $Id: cva_driver.F90 5142 2024-07-29 13:07:34Z 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! **************************************************************
32! *
33! CV_DRIVER                                                   *
34! *
35! *
36! written by   : Sandrine Bony-Lena , 17/05/2003, 11.19.41    *
37! modified by :                                               *
38! **************************************************************
39! **************************************************************
40
41  USE lmdz_print_control, ONLY: prt_level, lunout
42  USE add_phys_tend_mod, ONLY: fl_cor_ebil
43  USE lmdz_abort_physic, ONLY: abort_physic
44  USE lmdz_cv, ONLY: cv_param, cv_prelim, cv_feed, cv_undilute1, cv_trigger, cv_compress, &
45          cv_undilute2, cv_closure, cv_mixing, cv_unsat, cv_yield, cv_uncompress
46
47  IMPLICIT NONE
48
49! .............................START PROLOGUE............................
50
51
52! All argument names (except len,nd,ntra,nloc,delt and the flags) have a "1" appended.
53! The "1" is removed for the corresponding compressed variables.
54! PARAMETERS:
55! Name            Type         Usage            Description
56! ----------      ----------     -------  ----------------------------
57
58! len           Integer        Input        first (i) dimension
59! nd            Integer        Input        vertical (k) dimension
60! ndp1          Integer        Input        nd + 1
61! ntra          Integer        Input        number of tracors
62! nloc          Integer        Input        dimension of arrays for compressed fields
63! k_upper       Integer        Input        upmost level for vertical loops
64! iflag_con     Integer        Input        version of convect (3/4)
65! iflag_mix     Integer        Input        version of mixing  (0/1/2)
66! iflag_ice_thermo Integer        Input        accounting for ice thermodynamics (0/1)
67! iflag_clos    Integer        Input        version of closure (0/1)
68! tau_cld_cv    Real           Input        characteristic time of dissipation of mixing fluxes
69! coefw_cld_cv  Real           Input        coefficient for updraft velocity in convection
70! ok_conserv_q  Logical        Input        when true corrections for water conservation are swtiched on
71! delt          Real           Input        time step
72! comp_threshold Real           Input       threshold on the fraction of convective points below which
73!                                            fields  are compressed
74! t1            Real           Input        temperature (sat draught envt)
75! q1            Real           Input        specific hum (sat draught envt)
76! qs1           Real           Input        sat specific hum (sat draught envt)
77! t1_wake       Real           Input        temperature (unsat draught envt)
78! q1_wake       Real           Input        specific hum(unsat draught envt)
79! qs1_wake      Real           Input        sat specific hum(unsat draughts envt)
80! s1_wake       Real           Input        fractionnal area covered by wakes
81! u1            Real           Input        u-wind
82! v1            Real           Input        v-wind
83! tra1          Real           Input        tracors
84! p1            Real           Input        full level pressure
85! ph1           Real           Input        half level pressure
86! ALE1          Real           Input        Available lifting Energy
87! ALP1          Real           Input        Available lifting Power
88! sig1feed1     Real           Input        sigma coord at lower bound of feeding layer
89! sig2feed1     Real           Input        sigma coord at upper bound of feeding layer
90! wght1         Real           Input        weight density determining the feeding mixture
91! iflag1        Integer        Output       flag for Emanuel conditions
92! ft1           Real           Output       temp tend
93! fq1           Real           Output       spec hum tend
94! fqcomp1       Real           Output       spec hum tend (only mixed draughts)
95! fu1           Real           Output       u-wind tend
96! fv1           Real           Output       v-wind tend
97! ftra1         Real           Output       tracor tend
98! precip1       Real           Output       precipitation
99! kbas1         Integer        Output       cloud base level
100! ktop1         Integer        Output       cloud top level
101! cbmf1         Real           Output       cloud base mass flux
102! sig1          Real           In/Out       section adiabatic updraft
103! w01           Real           In/Out       vertical velocity within adiab updraft
104! ptop21        Real           In/Out       top of entraining zone
105! Ma1           Real           Output       mass flux adiabatic updraft
106! mip1          Real           Output       mass flux shed by the adiabatic updraft
107! Vprecip1      Real           Output       vertical profile of total precipitation
108! Vprecipi1     Real           Output       vertical profile of ice precipitation
109! upwd1         Real           Output       total upward mass flux (adiab+mixed)
110! dnwd1         Real           Output       saturated downward mass flux (mixed)
111! dnwd01        Real           Output       unsaturated downward mass flux
112! qcondc1       Real           Output       in-cld mixing ratio of condensed water
113! wd1           Real           Output       downdraft velocity scale for sfc fluxes
114! cape1         Real           Output       CAPE
115! cin1          Real           Output       CIN
116! tvp1          Real           Output       adiab lifted parcell virt temp
117! ftd1          Real           Output       precip temp tend
118! fqt1          Real           Output       precip spec hum tend
119! Plim11        Real           Output
120! Plim21        Real           Output
121! asupmax1      Real           Output
122! supmax01      Real           Output
123! asupmaxmin1   Real           Output
124
125! ftd1          Real           Output  Array of temperature tendency due to precipitations (K/s) of dimension ND,
126!                                      defined at same grid levels as T, Q, QS and P.
127
128! fqd1          Real           Output  Array of specific humidity tendencies due to precipitations ((gm/gm)/s)
129!                                      of dimension ND, defined at same grid levels as T, Q, QS and P.
130
131! wdtrainA1     Real           Output   precipitation ejected from adiabatic draught;
132!                                         should be used in tracer transport (cvltr)
133! wdtrainS1     Real           Output   precipitation detrained from shedding of adiabatic draught;
134!                                         used in tracer transport (cvltr)
135! wdtrainM1     Real           Output   precipitation detrained from mixed draughts;
136!                                         used in tracer transport (cvltr)
137! da1           Real           Output     used in tracer transport (cvltr)
138! phi1          Real           Output     used in tracer transport (cvltr)
139! mp1           Real           Output     used in tracer transport (cvltr)
140! qtc1          Real           Output     specific humidity in convection
141! sigt1         Real           Output     surface fraction in adiabatic updrafts                                         
142! detrain1      Real           Output     detrainment terme klein
143! phi21         Real           Output     used in tracer transport (cvltr)
144                                         
145! d1a1          Real           Output     used in tracer transport (cvltr)
146! dam1          Real           Output     used in tracer transport (cvltr)
147                                         
148! epmlmMm1      Real           Output     used in tracer transport (cvltr)
149! eplaMm1       Real           Output     used in tracer transport (cvltr)
150                                         
151! evap1         Real           Output   
152! ep1           Real           Output   
153! sigij1        Real           Output     used in tracer transport (cvltr)
154! clw1          Real           Output   condensed water content of the adiabatic updraught
155! elij1         Real           Output
156! wghti1        Real           Output   final weight of the feeding layers,
157!                                         used in tracer transport (cvltr)
158
159
160! S. Bony, Mar 2002:
161! * Several modules corresponding to different physical processes
162! * Several versions of convect may be used:
163!         - iflag_con=3: version lmd  (previously named convect3)
164!         - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
165! + tard: - iflag_con=5: version lmd with ice (previously named convectg)
166! S. Bony, Oct 2002:
167! * Vectorization of convect3 (ie version lmd)
168
169! ..............................END PROLOGUE.............................
170
171
172
173! Input
174  INTEGER, INTENT (IN)                               :: len
175  INTEGER, INTENT (IN)                               :: nd
176  INTEGER, INTENT (IN)                               :: ndp1
177  INTEGER, INTENT (IN)                               :: ntra
178  INTEGER, INTENT(IN)                                :: nloc ! (nloc=len)  pour l'instant
179  INTEGER, INTENT (IN)                               :: k_upper
180  INTEGER, INTENT (IN)                               :: iflag_con
181  INTEGER, INTENT (IN)                               :: iflag_mix
182  INTEGER, INTENT (IN)                               :: iflag_ice_thermo
183  INTEGER, INTENT (IN)                               :: iflag_clos
184  LOGICAL, INTENT (IN)                               :: ok_conserv_q
185  REAL, INTENT (IN)                                  :: tau_cld_cv
186  REAL, INTENT (IN)                                  :: coefw_cld_cv
187  REAL, INTENT (IN)                                  :: delt
188  REAL, INTENT (IN)                                  :: comp_threshold
189  REAL, DIMENSION (len, nd), INTENT (IN)             :: t1
190  REAL, DIMENSION (len, nd), INTENT (IN)             :: q1
191  REAL, DIMENSION (len, nd), INTENT (IN)             :: qs1
192  REAL, DIMENSION (len, nd), INTENT (IN)             :: t1_wake
193  REAL, DIMENSION (len, nd), INTENT (IN)             :: q1_wake
194  REAL, DIMENSION (len, nd), INTENT (IN)             :: qs1_wake
195  REAL, DIMENSION (len), INTENT (IN)                 :: s1_wake
196  REAL, DIMENSION (len, nd), INTENT (IN)             :: u1
197  REAL, DIMENSION (len, nd), INTENT (IN)             :: v1
198  REAL, DIMENSION (len, nd, ntra), INTENT (IN)       :: tra1
199  REAL, DIMENSION (len, nd), INTENT (IN)             :: p1
200  REAL, DIMENSION (len, ndp1), INTENT (IN)           :: ph1
201  REAL, DIMENSION (len), INTENT (IN)                 :: Ale1
202  REAL, DIMENSION (len), INTENT (IN)                 :: Alp1
203  REAL, DIMENSION (len, nd), INTENT (IN)             :: omega1
204  REAL, INTENT (IN)                                  :: sig1feed1 ! pressure at lower bound of feeding layer
205  REAL, INTENT (IN)                                  :: sig2feed1 ! pressure at upper bound of feeding layer
206  REAL, DIMENSION (nd), INTENT (IN)                  :: wght1     ! weight density determining the feeding mixture
207  INTEGER, DIMENSION (len), INTENT (IN)              :: lalim_conv1
208
209! Input/Output
210  REAL, DIMENSION (len, nd), INTENT (INOUT)          :: sig1
211  REAL, DIMENSION (len, nd), INTENT (INOUT)          :: w01
212
213! Output
214  INTEGER, DIMENSION (len), INTENT (OUT)             :: iflag1
215  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ft1
216  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fq1
217  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fqcomp1
218  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fu1
219  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fv1
220  REAL, DIMENSION (len, nd, ntra), INTENT (OUT)      :: ftra1
221  REAL, DIMENSION (len), INTENT (OUT)                :: precip1
222  INTEGER, DIMENSION (len), INTENT (OUT)             :: kbas1
223  INTEGER, DIMENSION (len), INTENT (OUT)             :: ktop1
224  REAL, DIMENSION (len), INTENT (OUT)                :: cbmf1
225  REAL, DIMENSION (len), INTENT (OUT)                :: plcl1
226  REAL, DIMENSION (len), INTENT (OUT)                :: plfc1
227  REAL, DIMENSION (len), INTENT (OUT)                :: wbeff1
228  REAL, DIMENSION (len), INTENT (OUT)                :: ptop21
229  REAL, DIMENSION (len), INTENT (OUT)                :: sigd1
230  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ma1        ! adiab. asc. mass flux (staggered grid)
231  REAL, DIMENSION (len, nd), INTENT (OUT)            :: mip1       ! mass flux shed from adiab. ascent (extensive)
232! real Vprecip1(len,nd)
233  REAL, DIMENSION (len, ndp1), INTENT (OUT)          :: vprecip1   ! tot precipitation flux (staggered grid)
234  REAL, DIMENSION (len, ndp1), INTENT (OUT)          :: vprecipi1  ! ice precipitation flux (staggered grid)
235  REAL, DIMENSION (len, nd), INTENT (OUT)            :: upwd1      ! upwd sat. mass flux (staggered grid)
236  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dnwd1      ! dnwd sat. mass flux (staggered grid)
237  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dnwd01     ! unsat. mass flux (staggered grid)
238  REAL, DIMENSION (len, nd), INTENT (OUT)            :: qcondc1    ! max cloud condensate (intensive)  ! cld
239  REAL, DIMENSION (len), INTENT (OUT)                :: wd1             ! gust
240  REAL, DIMENSION (len), INTENT (OUT)                :: cape1
241  REAL, DIMENSION (len), INTENT (OUT)                :: cin1
242  REAL, DIMENSION (len, nd), INTENT (OUT)            :: tvp1       ! Virt. temp. in the adiab. ascent
243
244!AC!
245!!      real da1(len,nd),phi1(len,nd,nd)
246!!      real da(len,nd),phi(len,nd,nd)
247!AC!
248  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ftd1       ! Temp. tendency due to the sole unsat. drafts
249  REAL, DIMENSION (len, nd), INTENT (OUT)            :: fqd1       ! Moist. tendency due to the sole unsat. drafts
250  REAL, DIMENSION (len), INTENT (OUT)                :: Plim11
251  REAL, DIMENSION (len), INTENT (OUT)                :: Plim21
252  REAL, DIMENSION (len, nd), INTENT (OUT)            :: asupmax1   ! Highest mixing fraction of mixed updraughts
253  REAL, DIMENSION (len), INTENT (OUT)                :: supmax01
254  REAL, DIMENSION (len), INTENT (OUT)                :: asupmaxmin1
255  REAL, DIMENSION (len, nd), INTENT (OUT)            :: qtc1    ! in cloud water content (intensive)   ! cld
256  REAL, DIMENSION (len, nd), INTENT (OUT)            :: sigt1   ! fract. cloud area (intensive)        ! cld
257  REAL, DIMENSION (len, nd), INTENT (OUT)            :: detrain1   ! detrainement term of mixed draughts in environment
258
259! RomP >>>
260  REAL, DIMENSION (len, nd), INTENT (OUT)            :: wdtrainA1, wdtrainS1, wdtrainM1 ! precipitation sources (extensive)
261  REAL, DIMENSION (len, nd), INTENT (OUT)            :: mp1  ! unsat. mass flux (staggered grid)
262  REAL, DIMENSION (len, nd), INTENT (OUT)            :: da1  ! detrained mass flux of adiab. asc. air (extensive)
263  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: phi1 ! mass flux of envt. air in mixed draughts (extensive)
264  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: epmlmMm1  ! (extensive)
265  REAL, DIMENSION (len, nd), INTENT (OUT)            :: eplaMm1   ! (extensive)
266  REAL, DIMENSION (len, nd), INTENT (OUT)            :: evap1 ! evaporation rate in precip. downdraft. (intensive)
267  REAL, DIMENSION (len, nd), INTENT (OUT)            :: ep1
268  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: sigij1 ! mass fraction of env. air in mixed draughts (intensive)
269  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: elij1! cond. water per unit mass of mixed draughts (intensive)
270  REAL, DIMENSION (len, nd), INTENT (OUT)            :: qta1 ! total water per unit mass of the adiab. asc. (intensive)
271  REAL, DIMENSION (len, nd), INTENT (OUT)            :: clw1 ! cond. water per unit mass of the adiab. asc. (intensive)
272!JYG,RL
273  REAL, DIMENSION (len, nd), INTENT (OUT)            :: wghti1   ! final weight of the feeding layers (extensive)
274!JYG,RL
275  REAL, DIMENSION (len, nd, nd), INTENT (OUT)        :: phi21    ! (extensive)
276  REAL, DIMENSION (len, nd), INTENT (OUT)            :: d1a1     ! (extensive)
277  REAL, DIMENSION (len, nd), INTENT (OUT)            :: dam1     ! (extensive)
278! RomP <<<
279  REAL, DIMENSION (len ), INTENT (OUT)               :: epmax_diag1     
280
281! -------------------------------------------------------------------
282! Prolog by Kerry Emanuel.
283! -------------------------------------------------------------------
284! --- ARGUMENTS
285! -------------------------------------------------------------------
286! --- On input:
287
288! t:   Array of absolute temperature (K) of dimension ND, with first
289! index corresponding to lowest model level. Note that this array
290! will be altered by the SUBROUTINE if dry convective adjustment
291! occurs and if IPBL is not equal to 0.
292
293! q:   Array of specific humidity (gm/gm) of dimension ND, with first
294! index corresponding to lowest model level. Must be defined
295! at same grid levels as T. Note that this array will be altered
296! if dry convective adjustment occurs and if IPBL is not equal to 0.
297
298! qs:  Array of saturation specific humidity of dimension ND, with first
299! index corresponding to lowest model level. Must be defined
300! at same grid levels as T. Note that this array will be altered
301! if dry convective adjustment occurs and if IPBL is not equal to 0.
302
303! t_wake: Array of absolute temperature (K), seen by unsaturated draughts,
304! of dimension ND, with first index corresponding to lowest model level.
305
306! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts,
307! of dimension ND, with first index corresponding to lowest model level.
308! Must be defined at same grid levels as T.
309
310! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts,
311! of dimension ND, with first index corresponding to lowest model level.
312! Must be defined at same grid levels as T.
313
314! s_wake: Array of fractionnal area occupied by the wakes.
315
316! u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
317! index corresponding with the lowest model level. Defined at
318! same levels as T. Note that this array will be altered if
319! dry convective adjustment occurs and if IPBL is not equal to 0.
320
321! v:   Same as u but for meridional velocity.
322
323! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
324! where NTRA is the number of different tracers. If no
325! convective tracer transport is needed, define a dummy
326! input array of dimension (ND,1). Tracers are defined at
327! same vertical levels as T. Note that this array will be altered
328! if dry convective adjustment occurs and if IPBL is not equal to 0.
329
330! p:   Array of pressure (mb) of dimension ND, with first
331! index corresponding to lowest model level. Must be defined
332! at same grid levels as T.
333
334! ph:  Array of pressure (mb) of dimension ND+1, with first index
335! corresponding to lowest level. These pressures are defined at
336! levels intermediate between those of P, T, Q and QS. The first
337! value of PH should be greater than (i.e. at a lower level than)
338! the first value of the array P.
339
340! ALE:  Available lifting Energy
341
342! ALP:  Available lifting Power
343
344! nl:  The maximum number of levels to which convection can penetrate, plus 1.
345!       NL MUST be less than or equal to ND-1.
346
347! delt: The model time step (sec) between calls to CONVECT
348
349! ----------------------------------------------------------------------------
350! ---   On Output:
351
352! iflag: An output integer whose value denotes the following:
353!       VALUE   INTERPRETATION
354!       -----   --------------
355!         0     Moist convection occurs.
356!         1     Moist convection occurs, but a CFL condition
357!               on the subsidence warming is violated. This
358!               does not cause the scheme to terminate.
359!         2     Moist convection, but no precip because ep(inb) lt 0.0001
360!         3     No moist convection because new cbmf is 0 and old cbmf is 0.
361!         4     No moist convection; atmosphere is not
362!               unstable
363!         6     No moist convection because ihmin le minorig.
364!         7     No moist convection because unreasonable
365!               parcel level temperature or specific humidity.
366!         8     No moist convection: lifted condensation
367!               level is above the 200 mb level.
368!         9     No moist convection: cloud base is higher
369!               then the level NL-1.
370!        10     No moist convection: cloud top is too warm.
371!        14     No moist convection; atmosphere is very
372!               stable (=> no computation)
373
374! ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
375!       grid levels as T, Q, QS and P.
376
377! fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
378!       defined at same grid levels as T, Q, QS and P.
379
380! fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
381!      defined at same grid levels as T.
382
383! fv:   Same as FU, but for forcing of meridional velocity.
384
385! ftra: Array of forcing of tracer content, in tracer mixing ratio per
386!       second, defined at same levels as T. Dimensioned (ND,NTRA).
387
388! precip: Scalar convective precipitation rate (mm/day).
389
390! wd:   A convective downdraft velocity scale. For use in surface
391!       flux parameterizations. See convect.ps file for details.
392
393! tprime: A convective downdraft temperature perturbation scale (K).
394!         For use in surface flux parameterizations. See convect.ps
395!         file for details.
396
397! qprime: A convective downdraft specific humidity
398!         perturbation scale (gm/gm).
399!         For use in surface flux parameterizations. See convect.ps
400!         file for details.
401
402! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
403!       BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
404!       ITS NEXT CALL. That is, the value of CBMF must be "remembered"
405!       by the calling program between calls to CONVECT.
406
407! det:   Array of detrainment mass flux of dimension ND.
408! -------------------------------------------------------------------
409
410! Local (non compressed) arrays
411
412
413  INTEGER i, k, il
414  INTEGER nword1, nword2, nword3, nword4
415  INTEGER icbmax
416  INTEGER nk1(len)
417  INTEGER icb1(len)
418  INTEGER icbs1(len)
419
420  LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
421  LOGICAL, SAVE :: debut = .TRUE.
422!$OMP THREADPRIVATE(debut)
423
424  REAL coef_convective(len)   ! = 1 for convective points, = 0 otherwise
425  REAL tnk1(len)
426  REAL thnk1(len)
427  REAL qnk1(len)
428  REAL gznk1(len)
429  REAL qsnk1(len)
430  REAL unk1(len)
431  REAL vnk1(len)
432  REAL cpnk1(len)
433  REAL hnk1(len)
434  REAL pbase1(len)
435  REAL buoybase1(len)
436
437  REAL lf1(len, nd), lf1_wake(len, nd)
438  REAL lv1(len, nd), lv1_wake(len, nd)
439  REAL cpn1(len, nd), cpn1_wake(len, nd)
440  REAL tv1(len, nd), tv1_wake(len, nd)
441  REAL gz1(len, nd), gz1_wake(len, nd)
442  REAL hm1(len, nd)
443  REAL h1(len, nd), h1_wake(len, nd)
444  REAL tp1(len, nd)
445  REAL th1(len, nd), th1_wake(len, nd)
446
447  REAL bid(len, nd) ! dummy array
448
449  INTEGER ncum
450
451  REAL p1feed1(len) ! pressure at lower bound of feeding layer
452  REAL p2feed1(len) ! pressure at upper bound of feeding layer
453!JYG,RL
454!!      real wghti1(len,nd) ! weights of the feeding layers
455!JYG,RL
456
457! (local) compressed fields:
458
459
460  INTEGER idcum(nloc)
461!jyg<
462  LOGICAL compress    ! True if compression occurs
463!>jyg
464  INTEGER iflag(nloc), nk(nloc), icb(nloc)
465  INTEGER nent(nloc, nd)
466  INTEGER icbs(nloc)
467  INTEGER inb(nloc), inbis(nloc)
468
469  REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc)
470  REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd)
471  REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd)
472  REAL s_wake(nloc)
473  REAL u(nloc, nd), v(nloc, nd)
474  REAL gz(nloc, nd), h(nloc, nd)
475  REAL h_wake(nloc, nd)
476  REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd)
477  REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd)
478  REAL p(nloc, nd), ph(nloc, nd+1), tv(nloc, nd), tp(nloc, nd)
479  REAL tv_wake(nloc, nd)
480  REAL clw(nloc, nd)
481  REAL, DIMENSION(nloc, nd)    :: qta, qpreca                       !!jygprl
482  REAL dph(nloc, nd)
483  REAL pbase(nloc), buoybase(nloc), th(nloc, nd)
484  REAL th_wake(nloc, nd)
485  REAL tvp(nloc, nd)
486  REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc)
487  REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd)
488  REAL buoy(nloc, nd)
489  REAL cape(nloc)
490  REAL cin(nloc)
491  REAL m(nloc, nd)
492  REAL mm(nloc, nd)
493  REAL ment(nloc, nd, nd), sigij(nloc, nd, nd)
494  REAL qent(nloc, nd, nd)
495  REAL hent(nloc, nd, nd)
496  REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
497  REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
498  REAL elij(nloc, nd, nd)
499  REAL supmax(nloc, nd)
500  REAL Ale(nloc), Alp(nloc), coef_clos(nloc)
501  REAL omega(nloc,nd)
502  REAL sigd(nloc)
503! real mp(nloc,nd), qp(nloc,nd), up(nloc,nd), vp(nloc,nd)
504! real wt(nloc,nd), water(nloc,nd), evap(nloc,nd), ice(nloc,nd)
505! real b(nloc,nd), sigd(nloc)
506! save mp,qp,up,vp,wt,water,evap,b
507  REAL, DIMENSION(len,nd)     :: mp, qp, up, vp
508  REAL, DIMENSION(len,nd)     :: wt, water, evap
509  REAL, DIMENSION(len,nd)     :: ice, fondue, b
510  REAL, DIMENSION(len,nd)     :: frac_a, frac_s, faci               !!jygprl
511  REAL ft(nloc, nd), fq(nloc, nd), fqcomp(nloc, nd)
512  REAL ftd(nloc, nd), fqd(nloc, nd)
513  REAL fu(nloc, nd), fv(nloc, nd)
514  REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd)
515  REAL ma(nloc, nd), mip(nloc, nd)
516!!  REAL tls(nloc, nd), tps(nloc, nd)                 ! unused . jyg
517  REAL qprime(nloc), tprime(nloc)
518  REAL precip(nloc)
519! real Vprecip(nloc,nd)
520  REAL vprecip(nloc, nd+1)
521  REAL vprecipi(nloc, nd+1)
522  REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
523  REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
524  REAL qcondc(nloc, nd)      ! cld
525  REAL wd(nloc)                ! gust
526  REAL Plim1(nloc), plim2(nloc)
527  REAL asupmax(nloc, nd)
528  REAL supmax0(nloc)
529  REAL asupmaxmin(nloc)
530
531  REAL tnk(nloc), qnk(nloc), gznk(nloc)
532  REAL wghti(nloc, nd)
533  REAL hnk(nloc), unk(nloc), vnk(nloc)
534
535  REAL qtc(nloc, nd)         ! cld
536  REAL sigt(nloc, nd)        ! cld
537  REAL detrain(nloc, nd)     ! cld
538 
539! RomP >>>
540  REAL wdtrainA(nloc, nd), wdtrainS(nloc, nd), wdtrainM(nloc, nd)   !!jygprl
541  REAL da(len, nd), phi(len, nd, nd)
542  REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd)
543  REAL phi2(len, nd, nd)
544  REAL d1a(len, nd), dam(len, nd)
545! RomP <<<
546  REAL epmax_diag(nloc) ! epmax_cape
547
548  CHARACTER (LEN=20) :: modname = 'cva_driver'
549  CHARACTER (LEN=80) :: abort_message
550
551  REAL, PARAMETER    :: Cin_noconv = -100000.
552  REAL, PARAMETER    :: Cape_noconv = -1.
553
554  INTEGER,SAVE                                       :: igout=1
555!$OMP THREADPRIVATE(igout)
556
557
558! print *, 't1, t1_wake ',(k,t1(1,k),t1_wake(1,k),k=1,nd)
559! print *, 'q1, q1_wake ',(k,q1(1,k),q1_wake(1,k),k=1,nd)
560
561! -------------------------------------------------------------------
562! --- SET CONSTANTS AND PARAMETERS
563! -------------------------------------------------------------------
564
565! -- set simulation flags:
566! (common cvflag)
567
568  CALL cv_flag(iflag_ice_thermo)
569
570! -- set thermodynamical constants:
571! (common cvthermo)
572
573  CALL cv_thermo(iflag_con)
574
575! -- set convect parameters
576
577! includes microphysical parameters and parameters that
578! control the rate of approach to quasi-equilibrium)
579! (common cvparam)
580
581  IF (iflag_con==3) THEN
582    CALL cv3_param(nd, k_upper, delt)
583
584  END IF
585
586  IF (iflag_con==4) THEN
587    CALL cv_param(nd)
588  END IF
589
590! ---------------------------------------------------------------------
591! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
592! ---------------------------------------------------------------------
593  nword1 = len
594  nword2 = len*nd
595  nword3 = len*nd*ntra
596  nword4 = len*nd*nd
597
598  iflag1(:) = 0
599  ktop1(:) = 0
600  kbas1(:) = 0
601  ft1(:, :) = 0.0
602  fq1(:, :) = 0.0
603  fqcomp1(:, :) = 0.0
604  fu1(:, :) = 0.0
605  fv1(:, :) = 0.0
606  ftra1(:, :, :) = 0.
607  precip1(:) = 0.
608  cbmf1(:) = 0.
609  plcl1(:) = 0.
610  plfc1(:) = 0.
611  wbeff1(:) = 0.
612  ptop21(:) = 0.
613  sigd1(:) = 0.
614  ma1(:, :) = 0.
615  mip1(:, :) = 0.
616  vprecip1(:, :) = 0.
617  vprecipi1(:, :) = 0.
618  upwd1(:, :) = 0.
619  dnwd1(:, :) = 0.
620  dnwd01(:, :) = 0.
621  qcondc1(:, :) = 0.
622  wd1(:) = 0.
623  cape1(:) = 0.
624  cin1(:) = 0.
625  tvp1(:, :) = 0.
626  ftd1(:, :) = 0.
627  fqd1(:, :) = 0.
628  Plim11(:) = 0.
629  Plim21(:) = 0.
630  asupmax1(:, :) = 0.
631  supmax01(:) = 0.
632  asupmaxmin1(:) = 0.
633
634  tvp(:, :) = 0. !ym missing init, need to have a look by developpers
635  tv(:, :) = 0. !ym missing init, need to have a look by developpers
636
637  DO il = 1, len
638!!    cin1(il) = -100000.
639!!    cape1(il) = -1.
640    cin1(il) = Cin_noconv
641    cape1(il) = Cape_noconv
642  END DO
643
644!!  IF (iflag_con==3) THEN
645!!    DO il = 1, len
646!!      sig1(il, nd) = sig1(il, nd) + 1.
647!!      sig1(il, nd) = amin1(sig1(il,nd), 12.1)
648!!    END DO
649!!  END IF
650
651  IF (iflag_con==3) THEN
652      CALL cv3_incrcount(len,nd,delt,sig1)
653  END IF  ! (iflag_con==3)
654
655! RomP >>>
656  sigt1(:, :) = 0.
657  detrain1(:, :) = 0.
658  qtc1(:, :) = 0.
659  wdtrainA1(:, :) = 0.
660  wdtrainS1(:, :) = 0.
661  wdtrainM1(:, :) = 0.
662  da1(:, :) = 0.
663  phi1(:, :, :) = 0.
664  epmlmMm1(:, :, :) = 0.
665  eplaMm1(:, :) = 0.
666  mp1(:, :) = 0.
667  evap1(:, :) = 0.
668  ep1(:, :) = 0.
669  sigij1(:, :, :) = 0.
670  elij1(:, :, :) = 0.
671  qta1(:,:) = 0.
672  clw1(:,:) = 0.
673  wghti1(:,:) = 0.
674  phi21(:, :, :) = 0.
675  d1a1(:, :) = 0.
676  dam1(:, :) = 0.
677! RomP <<<
678! ---------------------------------------------------------------------
679! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS
680! ---------------------------------------------------------------------
681
682  DO il = 1, nloc
683    coef_clos(il) = 1.
684  END DO
685
686! --------------------------------------------------------------------
687! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
688! --------------------------------------------------------------------
689
690  IF (iflag_con==3) THEN
691
692    IF (debut) THEN
693      PRINT *, 'Emanuel version 3 nouvelle'
694    END IF
695! PRINT*,'t1, q1 ',t1,q1
696        IF (prt_level >= 9) &
697             PRINT *, 'cva_driver -> cv3_prelim'
698    CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, &           ! nd->na
699                    lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1)
700
701
702        IF (prt_level >= 9) &
703             PRINT *, 'cva_driver -> cv3_prelim'
704    CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, & ! nd->na
705                    lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, &
706                    h1_wake, bid, th1_wake)
707
708  END IF
709
710  IF (iflag_con==4) THEN
711    PRINT *, 'Emanuel version 4 '
712        IF (prt_level >= 9) &
713             PRINT *, 'cva_driver -> cv_prelim'
714    CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1, &
715                   lv1, cpn1, tv1, gz1, h1, hm1)
716  END IF
717
718! --------------------------------------------------------------------
719! --- CONVECTIVE FEED
720! --------------------------------------------------------------------
721
722! compute feeding layer potential temperature and mixing ratio :
723
724! get bounds of feeding layer
725
726! test niveaux couche alimentation KE
727  IF (sig1feed1==sig2feed1) THEN
728    WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed'
729    WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def'
730    abort_message = ''
731    CALL abort_physic(modname, abort_message, 1)
732  END IF
733
734  DO i = 1, len
735    p1feed1(i) = sig1feed1*ph1(i, 1)
736    p2feed1(i) = sig2feed1*ph1(i, 1)
737!test maf
738!   p1feed1(i)=ph1(i,1)
739!   p2feed1(i)=ph1(i,2)
740!   p2feed1(i)=ph1(i,3)
741!testCR: on prend la couche alim des thermiques
742!   p2feed1(i)=ph1(i,lalim_conv1(i)+1)
743!   PRINT*,'lentr=',lentr(i),ph1(i,lentr(i)+1),ph1(i,2)
744  END DO
745
746  IF (iflag_con==3) THEN
747  END IF
748  DO i = 1, len
749! PRINT*,'avant cv3_feed Plim',p1feed1(i),p2feed1(i)
750  END DO
751  IF (iflag_con==3) THEN
752
753! PRINT*, 'IFLAG1 avant cv3_feed'
754! PRINT*,'len,nd',len,nd
755! WRITE(*,'(64i1)') iflag1(2:len-1)
756
757        IF (prt_level >= 9) &
758             PRINT *, 'cva_driver -> cv3_feed'
759    CALL cv3_feed(len, nd, ok_conserv_q, &                 ! nd->na
760                  t1, q1, u1, v1, p1, ph1, h1, gz1, &
761                  p1feed1, p2feed1, wght1, &
762                  wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, &
763                  cpnk1, hnk1, nk1, icb1, icbmax, iflag1, gznk1, plcl1)
764  END IF
765
766! PRINT*, 'IFLAG1 apres cv3_feed'
767! PRINT*,'len,nd',len,nd
768! WRITE(*,'(64i1)') iflag1(2:len-1)
769
770  IF (iflag_con==4) THEN
771        IF (prt_level >= 9) &
772             PRINT *, 'cva_driver -> cv_feed'
773    CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1, &
774                 nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
775  END IF
776
777! print *, 'cv3_feed-> iflag1, plcl1 ',iflag1(1),plcl1(1)
778
779! --------------------------------------------------------------------
780! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
781! (up through ICB for convect4, up through ICB+1 for convect3)
782! Calculates the lifted parcel virtual temperature at nk, the
783! actual temperature, and the adiabatic liquid water content.
784! --------------------------------------------------------------------
785
786  IF (iflag_con==3) THEN
787
788        IF (prt_level >= 9) &
789             PRINT *, 'cva_driver -> cv3_undilute1'
790    CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, icb1, tnk1, qnk1, & ! nd->na
791                       gznk1, tp1, tvp1, clw1, icbs1)
792  END IF
793
794
795  IF (iflag_con==4) THEN
796        IF (prt_level >= 9) &
797             PRINT *, 'cva_driver -> cv_undilute1'
798    CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, &
799                      tp1, tvp1, clw1)
800  END IF
801
802! -------------------------------------------------------------------
803! --- TRIGGERING
804! -------------------------------------------------------------------
805
806! print *,' avant triggering, iflag_con ',iflag_con
807
808  IF (iflag_con==3) THEN
809
810        IF (prt_level >= 9) &
811             PRINT *, 'cva_driver -> cv3_trigger'
812    CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na
813                      pbase1, buoybase1, iflag1, sig1, w01)
814
815
816! PRINT*, 'IFLAG1 apres cv3_triger'
817! PRINT*,'len,nd',len,nd
818! WRITE(*,'(64i1)') iflag1(2:len-1)
819
820! CALL dump2d(iim,jjm-1,sig1(2)
821  END IF
822
823  IF (iflag_con==4) THEN
824        IF (prt_level >= 9) &
825             PRINT *, 'cva_driver -> cv_trigger'
826    CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
827  END IF
828
829
830! =====================================================================
831! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
832! =====================================================================
833
834!  Determine the number "ncum" of convective gridpoints, the list "idcum" of convective
835!  gridpoints and the weights "coef_convective" (= 1. for convective gridpoints and = 0.
836!  elsewhere).
837  ncum = 0
838  coef_convective(:) = 0.
839  DO i = 1, len
840    IF (iflag1(i)==0) THEN
841      coef_convective(i) = 1.
842      ncum = ncum + 1
843      idcum(ncum) = i
844    END IF
845  END DO
846
847! PRINT*,'len, ncum = ',len,ncum
848
849  IF (ncum>0) THEN
850
851! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
852! --- COMPRESS THE FIELDS
853!       (-> vectorization over convective gridpoints)
854! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
855
856    IF (iflag_con==3) THEN
857! PRINT*,'ncum tv1 ',ncum,tv1
858! PRINT*,'tvp1 ',tvp1
859!jyg<
860!   If the fraction of convective points is larger than comp_threshold, then compression
861!   is assumed useless.
862
863  compress = ncum < len*comp_threshold
864
865  IF (.NOT. compress) THEN
866    DO i = 1,len
867      idcum(i) = i
868    ENDDO
869  ENDIF
870
871!>jyg
872        IF (prt_level >= 9) &
873             PRINT *, 'cva_driver -> cv3a_compress'
874      CALL cv3a_compress(len, nloc, ncum, nd, ntra, compress, &
875                         iflag1, nk1, icb1, icbs1, &
876                         plcl1, tnk1, qnk1, gznk1, hnk1, unk1, vnk1, &
877                         wghti1, pbase1, buoybase1, &
878                         t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, &
879                         u1, v1, gz1, th1, th1_wake, &
880                         tra1, &
881                         h1, lv1, lf1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
882                         h1_wake, lv1_wake, lf1_wake, cpn1_wake, tv1_wake, &
883                         sig1, w01, ptop21, &
884                         Ale1, Alp1, omega1, &
885                         iflag, nk, icb, icbs, &
886                         plcl, tnk, qnk, gznk, hnk, unk, vnk, &
887                         wghti, pbase, buoybase, &
888                         t, q, qs, t_wake, q_wake, qs_wake, s_wake, &
889                         u, v, gz, th, th_wake, &
890                         tra, &
891                         h, lv, lf, cpn, p, ph, tv, tp, tvp, clw, &
892                         h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, &
893                         sig, w0, ptop2, &
894                         Ale, Alp, omega)
895
896! PRINT*,'tv ',tv
897! PRINT*,'tvp ',tvp
898
899    END IF
900
901    IF (iflag_con==4) THEN
902        IF (prt_level >= 9) &
903             PRINT *, 'cva_driver -> cv_compress'
904      CALL cv_compress(len, nloc, ncum, nd, &
905                       iflag1, nk1, icb1, &
906                       cbmf1, plcl1, tnk1, qnk1, gznk1, &
907                       t1, q1, qs1, u1, v1, gz1, &
908                       h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
909                       iflag, nk, icb, &
910                       cbmf, plcl, tnk, qnk, gznk, &
911                       t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw, &
912                       dph)
913    END IF
914
915! -------------------------------------------------------------------
916! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
917! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
918! ---   &
919! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
920! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
921! ---   &
922! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
923! -------------------------------------------------------------------
924
925    IF (iflag_con==3) THEN
926        IF (prt_level >= 9) &
927             PRINT *, 'cva_driver -> cv3_undilute2'
928      CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, &        !na->nd
929                         tnk, qnk, gznk, hnk, t, q, qs, gz, &
930                         p, ph, h, tv, lv, lf, pbase, buoybase, plcl, &
931                         inb, tp, tvp, clw, hp, ep, sigp, buoy, &
932                         frac_a, frac_s, qpreca, qta)                        !!jygprl
933    END IF
934
935    IF (iflag_con==4) THEN
936        IF (prt_level >= 9) &
937             PRINT *, 'cva_driver -> cv_undilute2'
938      CALL cv_undilute2(nloc, ncum, nd, icb, nk, &
939                        tnk, qnk, gznk, t, q, qs, gz, &
940                        p, dph, h, tv, lv, &
941                        inb, inbis, tp, tvp, clw, hp, ep, sigp, frac_s)
942    END IF
943
944    ! epmax_cape
945    ! on recalcule ep et hp   
946        IF (prt_level >= 9) &
947             PRINT *, 'cva_driver -> cv3_epmax_cape'
948    CALL cv3_epmax_fn_cape(nloc,ncum,nd &
949                , ep,hp,icb,inb,clw,nk,t,h,hnk,lv,lf,frac_s &
950                , pbase, p, ph, tv, buoy, sig, w0,iflag &
951                , epmax_diag)
952
953! -------------------------------------------------------------------
954! --- MIXING(1)   (if iflag_mix .ge. 1)
955! -------------------------------------------------------------------
956    IF (iflag_con==3) THEN
957!      IF ((iflag_ice_thermo==1) .AND. (iflag_mix/=0)) THEN
958!        WRITE (*, *) ' iflag_ice_thermo==1 requires iflag_mix==0', ' but iflag_mix=', iflag_mix, &
959!          '. Might as well stop here.'
960!        STOP
961!      END IF
962      IF (iflag_mix>=1) THEN
963        CALL zilch(supmax, nloc*nd)
964        IF (prt_level >= 9) &
965             PRINT *, 'cva_driver -> cv3p_mixing'
966        CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &           ! na->nd
967!!                         ph, t, q, qs, u, v, tra, h, lv, lf, frac, qnk, &
968                         ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qta, &      !!jygprl
969                         unk, vnk, hp, tv, tvp, ep, clw, sig, &
970                         ment, qent, hent, uent, vent, nent, &
971                         sigij, elij, supmax, ments, qents, traent)
972! PRINT*, 'cv3p_mixing-> supmax ', (supmax(1,k), k=1,nd)
973
974      ELSE
975        CALL zilch(supmax, nloc*nd)
976      END IF
977    END IF
978! -------------------------------------------------------------------
979! --- CLOSURE
980! -------------------------------------------------------------------
981
982
983    IF (iflag_con==3) THEN
984      IF (iflag_clos==0) THEN
985        IF (prt_level >= 9) &
986             PRINT *, 'cva_driver -> cv3_closure'
987        CALL cv3_closure(nloc, ncum, nd, icb, inb, &           ! na->nd
988                         pbase, p, ph, tv, buoy, &
989                         sig, w0, cape, m, iflag)
990      END IF   ! iflag_clos==0
991
992      ok_inhib = iflag_mix == 2
993
994      IF (iflag_clos==1) THEN
995        PRINT *, ' pas d appel cv3p_closure'
996! c        CALL cv3p_closure(nloc,ncum,nd,icb,inb              ! na->nd
997! c    :                       ,pbase,plcl,p,ph,tv,tvp,buoy
998! c    :                       ,supmax
999! c    o                       ,sig,w0,ptop2,cape,cin,m)
1000      END IF   ! iflag_clos==1
1001
1002      IF (iflag_clos==2) THEN
1003        IF (prt_level >= 9) &
1004             PRINT *, 'cva_driver -> cv3p1_closure'
1005        CALL cv3p1_closure(nloc, ncum, nd, icb, inb, &         ! na->nd
1006                           pbase, plcl, p, ph, tv, tvp, buoy, &
1007                           supmax, ok_inhib, Ale, Alp, omega, &
1008                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
1009                           Plim1, plim2, asupmax, supmax0, &
1010                           asupmaxmin, cbmf, plfc, wbeff)
1011        IF (prt_level >= 10) &
1012             PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1)
1013      END IF   ! iflag_clos==2
1014
1015      IF (iflag_clos==3) THEN
1016        IF (prt_level >= 9) &
1017             PRINT *, 'cva_driver -> cv3p2_closure'
1018        CALL cv3p2_closure(nloc, ncum, nd, icb, inb, &         ! na->nd
1019                           pbase, plcl, p, ph, tv, tvp, buoy, &
1020                           supmax, ok_inhib, Ale, Alp, omega, &
1021                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
1022                           Plim1, plim2, asupmax, supmax0, &
1023                           asupmaxmin, cbmf, plfc, wbeff)
1024        IF (prt_level >= 10) &
1025             PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1)
1026      END IF   ! iflag_clos==3
1027    END IF ! iflag_con==3
1028
1029    IF (iflag_con==4) THEN
1030        IF (prt_level >= 9) &
1031             PRINT *, 'cva_driver -> cv_closure'
1032      CALL cv_closure(nloc, ncum, nd, nk, icb, &
1033                         tv, tvp, p, ph, dph, plcl, cpn, &
1034                         iflag, cbmf)
1035    END IF
1036
1037! print *,'cv_closure-> cape ',cape(1)
1038
1039! -------------------------------------------------------------------
1040! --- MIXING(2)
1041! -------------------------------------------------------------------
1042
1043    IF (iflag_con==3) THEN
1044      IF (iflag_mix==0) THEN
1045        IF (prt_level >= 9) &
1046             PRINT *, 'cva_driver -> cv3_mixing'
1047        CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &             ! na->nd
1048                        ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qnk, &
1049                        unk, vnk, hp, tv, tvp, ep, clw, m, sig, &
1050                        ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent)
1051        CALL zilch(hent, nloc*nd*nd)
1052      ELSE
1053!!jyg:  Essais absurde pour voir
1054!!        mm(:,1) = 0.
1055!!        DO  i = 2,nd
1056!!          mm(:,i) = m(:,i)*(1.-qta(:,i-1))
1057!!        ENDDO
1058        mm(:,:) = m(:,:)
1059        CALL cv3_mixscale(nloc, ncum, nd, ment, mm)
1060        IF (debut) THEN
1061          PRINT *, ' cv3_mixscale-> '
1062        END IF !(debut) THEN
1063      END IF
1064    END IF
1065
1066    IF (iflag_con==4) THEN
1067        IF (prt_level >= 9) &
1068             PRINT *, 'cva_driver -> cv_mixing'
1069      CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis, &
1070                     ph, t, q, qs, u, v, h, lv, qnk, &
1071                     hp, tv, tvp, ep, clw, cbmf, &
1072                     m, ment, qent, uent, vent, nent, sigij, elij)
1073    END IF                                                                                         
1074
1075    IF (debut) THEN
1076      PRINT *, ' cv_mixing ->'
1077    END IF !(debut) THEN
1078! do i = 1,nd
1079! PRINT*,'cv_mixing-> i,ment ',i,(ment(1,i,j),j=1,nd)
1080! enddo
1081
1082! -------------------------------------------------------------------
1083! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
1084! -------------------------------------------------------------------
1085    IF (iflag_con==3) THEN
1086      IF (debut) THEN
1087        PRINT *, ' cva_driver -> cv3_unsat '
1088      END IF !(debut) THEN
1089
1090        IF (prt_level >= 9) &
1091             PRINT *, 'cva_driver -> cv3_unsat'
1092      CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, &              ! na->nd
1093                     t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, &
1094                     th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, &
1095                     ep, sigp, clw, frac_s, qpreca, frac_a, qta, &                    !!jygprl
1096                     m, ment, elij, delt, plcl, coef_clos, &
1097                     mp, qp, up, vp, trap, wt, water, evap, fondue, ice, &
1098                     faci, b, sigd, &
1099!!                     wdtrainA, wdtrainM)                                       ! RomP
1100                     wdtrainA, wdtrainS, wdtrainM)                               !!jygprl
1101
1102      IF (prt_level >= 10) THEN
1103        Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue '
1104        DO k = 1,nd
1105        write (6, '(i4,5(1x,e13.6))') &
1106          k, mp(igout,k), water(igout,k), ice(igout,k), &
1107           evap(igout,k), fondue(igout,k)
1108        ENDDO
1109        Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainS, wdtrainM '     !!jygprl
1110        DO k = 1,nd
1111        write (6, '(i4,3(1x,e13.6))') &
1112           k, wdtrainA(igout,k), wdtrainS(igout,k), wdtrainM(igout,k)            !!jygprl
1113        ENDDO
1114      ENDIF
1115
1116    END IF  !(iflag_con==3)
1117
1118    IF (iflag_con==4) THEN
1119        IF (prt_level >= 9) &
1120             PRINT *, 'cva_driver -> cv_unsat'
1121      CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph, &
1122                     h, lv, ep, sigp, clw, m, ment, elij, &
1123                     iflag, mp, qp, up, vp, wt, water, evap)
1124    END IF
1125
1126    IF (debut) THEN
1127      PRINT *, 'cv_unsat-> '
1128    END IF !(debut) THEN
1129
1130! print *,'cv_unsat-> mp ',mp
1131! print *,'cv_unsat-> water ',water
1132! -------------------------------------------------------------------
1133! --- YIELD
1134! (tendencies, precipitation, variables of interface with other
1135! processes, etc)
1136! -------------------------------------------------------------------
1137
1138    IF (iflag_con==3) THEN
1139
1140        IF (prt_level >= 9) &
1141             PRINT *, 'cva_driver -> cv3_yield'
1142      CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, &                      ! na->nd
1143                     icb, inb, delt, &
1144                     t, q, t_wake, q_wake, s_wake, u, v, tra, &
1145                     gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, &
1146                     ep, clw, qpreca, m, tp, mp, qp, up, vp, trap, &
1147                     wt, water, ice, evap, fondue, faci, b, sigd, &
1148                     ment, qent, hent, iflag_mix, uent, vent, &
1149                     nent, elij, traent, sig, &
1150                     tv, tvp, wghti, &
1151                     iflag, precip, Vprecip, Vprecipi, ft, fq, fqcomp, fu, fv, ftra, &      ! jyg
1152                     cbmf, upwd, dnwd, dnwd0, ma, mip, &
1153!!                     tls, tps, &                            ! useless . jyg
1154                     qcondc, wd, &
1155!!                     ftd, fqd, qnk, qtc, sigt, tau_cld_cv, coefw_cld_cv)
1156                     ftd, fqd, qta, qtc, sigt, detrain, tau_cld_cv, coefw_cld_cv)         !!jygprl
1157
1158!         Test conseravtion de l'eau
1159
1160      IF (debut) THEN
1161        PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1)
1162      END IF !(debut) THEN
1163
1164      IF (prt_level >= 10) THEN
1165        Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', &
1166                    ft(igout,1), ftd(igout,1)
1167        Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', &
1168                    fq(igout,1), fqd(igout,1)
1169      ENDIF
1170
1171    END IF
1172
1173    IF (iflag_con==4) THEN
1174        IF (prt_level >= 9) &
1175             PRINT *, 'cva_driver -> cv_yield'
1176      CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt, &
1177                     t, q, u, v, &
1178                     gz, p, ph, h, hp, lv, cpn, &
1179                     ep, clw, frac_s, m, mp, qp, up, vp, &
1180                     wt, water, evap, &
1181                     ment, qent, uent, vent, nent, elij, &
1182                     tv, tvp, &
1183                     iflag, wd, qprime, tprime, &
1184                     precip, cbmf, ft, fq, fu, fv, ma, qcondc)
1185    END IF
1186
1187!AC!
1188!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1189!--- passive tracers
1190!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1191
1192    IF (iflag_con==3) THEN
1193!RomP >>>
1194        IF (prt_level >= 9) &
1195             PRINT *, 'cva_driver -> cv3_tracer'
1196      CALL cv3_tracer(nloc, len, ncum, nd, nd, &
1197                     ment, sigij, da, phi, phi2, d1a, dam, &
1198                     ep, vprecip, elij, clw, epmlmMm, eplaMm, &
1199                     icb, inb)
1200!RomP <<<
1201    END IF
1202
1203!AC!
1204
1205! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1206! --- UNCOMPRESS THE FIELDS
1207! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
1208
1209
1210    IF (iflag_con==3) THEN
1211        IF (prt_level >= 9) &
1212             PRINT *, 'cva_driver -> cv3a_uncompress'
1213      CALL cv3a_uncompress(nloc, len, ncum, nd, ntra, idcum, compress, &
1214                           iflag, icb, inb, &
1215                           precip, cbmf, plcl, plfc, wbeff, sig, w0, ptop2, &
1216                           ft, fq, fqcomp, fu, fv, ftra, &
1217                           sigd, ma, mip, vprecip, vprecipi, upwd, dnwd, dnwd0, &
1218                           qcondc, wd, cape, cin, &
1219                           tvp, &
1220                           ftd, fqd, &
1221                           Plim1, plim2, asupmax, supmax0, &
1222                           asupmaxmin, &
1223                           da, phi, mp, phi2, d1a, dam, sigij, &         ! RomP
1224                           qta, clw, elij, evap, ep, epmlmMm, eplaMm, &  ! RomP
1225                           wdtrainA, wdtrainS, wdtrainM, &                         ! RomP
1226                           qtc, sigt, detrain, epmax_diag, & ! epmax_cape
1227                           iflag1, kbas1, ktop1, &
1228                           precip1, cbmf1, plcl1, plfc1, wbeff1, sig1, w01, ptop21, &
1229                           ft1, fq1, fqcomp1, fu1, fv1, ftra1, &
1230                           sigd1, ma1, mip1, vprecip1, vprecipi1, upwd1, dnwd1, dnwd01, &
1231                           qcondc1, wd1, cape1, cin1, &
1232                           tvp1, &
1233                           ftd1, fqd1, &
1234                           Plim11, plim21, asupmax1, supmax01, &
1235                           asupmaxmin1, &
1236                           da1, phi1, mp1, phi21, d1a1, dam1, sigij1,  &       ! RomP
1237                           qta1, clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, & ! RomP
1238                           wdtrainA1, wdtrainS1, wdtrainM1,                  & ! RomP
1239                           qtc1, sigt1, detrain1, epmax_diag1) ! epmax_cape
1240
1241      IF (prt_level >= 10) THEN
1242        Print *, 'cva_driver after cv3_uncompress:ft1(1) , ftd1(1) ', &
1243                    ft1(igout,1), ftd1(igout,1)
1244        Print *, 'cva_driver after cv3_uncompress:fq1(1) , fqd1(1) ', &
1245                    fq1(igout,1), fqd1(igout,1)
1246      ENDIF
1247
1248    END IF
1249
1250    IF (iflag_con==4) THEN
1251        IF (prt_level >= 9) &
1252             PRINT *, 'cva_driver -> cv_uncompress'
1253      CALL cv_uncompress(nloc, len, ncum, nd, idcum, &
1254                           iflag, &
1255                           precip, cbmf, &
1256                           ft, fq, fu, fv, &
1257                           ma, qcondc, &
1258                           iflag1, &
1259                           precip1,cbmf1, &
1260                           ft1, fq1, fu1, fv1, &
1261                           ma1, qcondc1)
1262    END IF
1263
1264  END IF ! ncum>0
1265
1266
1267  DO i = 1,len
1268    IF (iflag1(i) == 14) THEN
1269      Cin1(i) = Cin_noconv
1270      Cape1(i) = Cape_noconv
1271    ENDIF
1272  ENDDO
1273
1274! In order take into account the possibility of changing the compression,
1275! reset m, sig and w0 to zero for non-convective points.
1276  DO k = 1,nd-1
1277        sig1(:, k) = sig1(:, k)*coef_convective(:)
1278        w01(:, k)  = w01(:, k)*coef_convective(:)
1279  ENDDO
1280
1281  IF (debut) THEN
1282    PRINT *, ' cv_uncompress -> '
1283    debut = .FALSE.
1284  END IF  !(debut) THEN
1285
1286
1287
1288END SUBROUTINE cva_driver
Note: See TracBrowser for help on using the repository browser.