source: LMDZ6/trunk/libf/phylmd/cva_driver.f90 @ 5694

Last change on this file since 5694 was 5694, checked in by yann meurdesoif, 6 weeks ago

Convection GPU porting : set convection subroutines into module (continued)

Files will be renamed later to *_mod.f90

YM

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