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

Last change on this file since 5746 was 5712, checked in by yann meurdesoif, 5 weeks ago

Convection GPU porting : Compression of active convection point is now optional (default remain to true). For GPU runs, convection is not compressed and is computed on each column. The update is done only for column where convection is active

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