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

Last change on this file since 5101 was 5101, checked in by abarral, 2 months ago

Handle DEBUG_IO in lmdz_cppkeys_wrapper.F90
Transform some files .F -> .[fF]90
[ne compile pas à cause de writefield_u non défini - en attente de réponse Laurent]

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