source: LMDZ6/branches/Optimisation_LMDZ/libf/phylmd/cv3a_driver.f90 @ 5040

Last change on this file since 5040 was 3765, checked in by adurocher, 4 years ago

Use cvl_comp_threshold to control whether to use compression or not

cvl_comp_threshold now defaults to 0 (never compress)

  • 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: 46.1 KB
Line 
1module cv3a_driver_mod
2  implicit none
3contains
4
5! -------------------------------------------------------------------
6! Prolog by Kerry Emanuel.
7! -------------------------------------------------------------------
8! --- ARGUMENTS
9! -------------------------------------------------------------------
10! --- On input:
11!
12! t:   Array of absolute temperature (K) of dimension ND, with first
13! index corresponding to lowest model level. Note that this array
14! will be altered by the subroutine if dry convective adjustment
15! occurs and if IPBL is not equal to 0.
16!
17! q:   Array of specific humidity (gm/gm) of dimension ND, with first
18! index corresponding to lowest model level. Must be defined
19! at same grid levels as T. Note that this array will be altered
20! if dry convective adjustment occurs and if IPBL is not equal to 0.
21!
22! qs:  Array of saturation specific humidity of dimension ND, with first
23! index corresponding to lowest model level. Must be defined
24! at same grid levels as T. Note that this array will be altered
25! if dry convective adjustment occurs and if IPBL is not equal to 0.
26!
27! t_wake: Array of absolute temperature (K), seen by unsaturated draughts,
28! of dimension ND, with first index corresponding to lowest model level.
29!
30! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts,
31! of dimension ND, with first index corresponding to lowest model level.
32! Must be defined at same grid levels as T.
33!
34! qs_wake: Array of saturation specific humidity, seen by unsaturated draughts,
35! of dimension ND, with first index corresponding to lowest model level.
36! Must be defined at same grid levels as T.
37!
38! s_wake: Array of fractionnal area occupied by the wakes.
39!
40! u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
41! index corresponding with the lowest model level. Defined at
42! same levels as T. Note that this array will be altered if
43! dry convective adjustment occurs and if IPBL is not equal to 0.
44!
45! v:   Same as u but for meridional velocity.
46!
47! tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
48! where NTRA is the number of different tracers. If no
49! convective tracer transport is needed, define a dummy
50! input array of dimension (ND,1). Tracers are defined at
51! same vertical levels as T. Note that this array will be altered
52! if dry convective adjustment occurs and if IPBL is not equal to 0.
53!
54! p:   Array of pressure (mb) of dimension ND, with first
55! index corresponding to lowest model level. Must be defined
56! at same grid levels as T.
57!
58! ph:  Array of pressure (mb) of dimension ND+1, with first index
59! corresponding to lowest level. These pressures are defined at
60! levels intermediate between those of P, T, Q and QS. The first
61! value of PH should be greater than (i.e. at a lower level than)
62! the first value of the array P.
63!
64! ALE:  Available lifting Energy
65!
66! ALP:  Available lifting Power
67!
68! nl:  The maximum number of levels to which convection can penetrate, plus 1.
69!       NL MUST be less than or equal to ND-1.
70!
71! delt: The model time step (sec) between calls to CONVECT
72!
73! ----------------------------------------------------------------------------
74! ---   On Output:
75!
76! iflag: An output integer whose value denotes the following:
77!       VALUE   INTERPRETATION
78!       -----   --------------
79!         0     Moist convection occurs.
80!         1     Moist convection occurs, but a CFL condition
81!               on the subsidence warming is violated. This
82!               does not cause the scheme to terminate.
83!         2     Moist convection, but no precip because ep(inb) lt 0.0001
84!         3     No moist convection because new cbmf is 0 and old cbmf is 0.
85!         4     No moist convection; atmosphere is not
86!               unstable
87!         6     No moist convection because ihmin le minorig.
88!         7     No moist convection because unreasonable
89!               parcel level temperature or specific humidity.
90!         8     No moist convection: lifted condensation
91!               level is above the 200 mb level.
92!         9     No moist convection: cloud base is higher
93!               then the level NL-1.
94!        10     No moist convection: cloud top is too warm.
95!
96!
97! ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
98!       grid levels as T, Q, QS and P.
99!
100! fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
101!       defined at same grid levels as T, Q, QS and P.
102!
103! fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
104!      defined at same grid levels as T.
105!
106! fv:   Same as FU, but for forcing of meridional velocity.
107!
108! ftra: Array of forcing of tracer content, in tracer mixing ratio per
109!       second, defined at same levels as T. Dimensioned (ND,NTRA).
110!
111! precip: Scalar convective precipitation rate (mm/day).
112!
113! wd:   A convective downdraft velocity scale. For use in surface
114!       flux parameterizations. See convect.ps file for details.
115!
116! tprime: A convective downdraft temperature perturbation scale (K).
117!         For use in surface flux parameterizations. See convect.ps
118!         file for details.
119!
120! qprime: A convective downdraft specific humidity
121!         perturbation scale (gm/gm).
122!         For use in surface flux parameterizations. See convect.ps
123!         file for details.
124
125! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
126!       BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
127!       ITS NEXT CALL. That is, the value of CBMF must be "remembered"
128!       by the calling program between calls to CONVECT.
129!
130! det:   Array of detrainment mass flux of dimension ND.
131! -------------------------------------------------------------------
132! written by   : Sandrine Bony-Lena , 17/05/2003, 11.19.41
133! S. Bony, Mar 2002:
134! * Several modules corresponding to different physical processes
135! * Several versions of convect may be used:
136!         - iflag_con=3: version lmd  (previously named convect3)
137!         - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
138! + tard: - iflag_con=5: version lmd with ice (previously named convectg)
139! S. Bony, Oct 2002:
140! * Vectorization of convect3 (ie version lmd)
141  SUBROUTINE cv3a_driver(len, nd, ndp1, ntra, nloc, k_upper, &
142                         iflag_con, iflag_mix, iflag_ice_thermo, iflag_clos, ok_conserv_q, &
143                         delt, comp_threshold, &
144                         t1, q1, qs1, t1_wake, q1_wake, qs1_wake, s1_wake, &
145                         u1, v1, &
146                         p1, ph1, &
147                         Ale1, Alp1, omega1, &
148                         sig1feed1, sig2feed1, wght1, &
149                         iflag1, ft1, fq1, fu1, fv1, ftra1, &
150                         precip1, kbas1, ktop1, &
151                         cbmf1, plcl1, plfc1, wbeff1, &
152                         sig1, w01, & !input/output
153                         ptop21, sigd1, &
154                         ma1, mip1, Vprecip1, Vprecipi1, upwd1, dnwd1, dnwd01, &
155                         qcondc1, wd1, &
156                         cape1, cin1, tvp1, &
157                         ftd1, fqd1, &
158                         Plim11, Plim21, asupmax1, supmax01, asupmaxmin1, &
159                         da1, phi1, mp1, phi21, d1a1, dam1, sigij1, wghti1, &
160                         qta1, clw1, elij1, evap1, ep1, epmlmMm1, eplaMm1, &
161                         wdtrainA1, wdtrainS1, wdtrainM1, qtc1, sigt1, tau_cld_cv, &
162                         coefw_cld_cv, &
163                         epmax_diag1)
164
165    USE print_control_mod, ONLY: prt_level, lunout
166    USE cv3a_compress_mod
167    USE profiling_physic_mod
168    IMPLICIT NONE
169! Input
170    INTEGER, INTENT(IN)                               :: len ! first (i) dimension
171    INTEGER, INTENT(IN)                               :: nd ! vertical (k) dimension
172    INTEGER, INTENT(IN)                               :: ndp1 ! nd + 1
173    INTEGER, INTENT(IN)                               :: ntra ! number of tracors
174    INTEGER, INTENT(IN)                               :: k_upper ! upmost level for vertical loops
175    INTEGER, INTENT(IN)                               :: iflag_con ! version of convect (3)
176    INTEGER, INTENT(IN)                               :: iflag_mix ! version of mixing  (0/1/2)
177    INTEGER, INTENT(IN)                               :: iflag_ice_thermo ! accounting for ice thermodynamics (0/1)
178    INTEGER, INTENT(IN)                               :: iflag_clos ! version of closure (0/1)
179    LOGICAL, INTENT(IN)                               :: ok_conserv_q ! when true corrections for water conservation are swtiched on
180    REAL, INTENT(IN)                                  :: tau_cld_cv ! characteristic time of dissipation of mixing fluxes
181    REAL, INTENT(IN)                                  :: coefw_cld_cv ! coefficient for updraft velocity in convection
182    REAL, INTENT(IN)                                  :: delt ! time step
183    REAL, INTENT (IN)                                 :: comp_threshold
184    REAL, DIMENSION(len, nd), INTENT(IN)             :: t1 ! temperature (sat draught envt)
185    REAL, DIMENSION(len, nd), INTENT(IN)             :: q1 ! specific hum (sat draught envt)
186    REAL, DIMENSION(len, nd), INTENT(IN)             :: qs1 ! sat specific hum (sat draught envt)
187    REAL, DIMENSION(len, nd), INTENT(IN)             :: t1_wake ! temperature (unsat draught envt)
188    REAL, DIMENSION(len, nd), INTENT(IN)             :: q1_wake ! specific hum(unsat draught envt)
189    REAL, DIMENSION(len, nd), INTENT(IN)             :: qs1_wake ! sat specific hum(unsat draughts envt)
190    REAL, DIMENSION(len), INTENT(IN)                 :: s1_wake ! fractionnal area covered by wakes
191    REAL, DIMENSION(len, nd), INTENT(IN)             :: u1 ! u-wind
192    REAL, DIMENSION(len, nd), INTENT(IN)             :: v1 ! v-wind
193    REAL, DIMENSION(len, nd), INTENT(IN)             :: p1 ! full level pressure
194    REAL, DIMENSION(len, ndp1), INTENT(IN)           :: ph1 ! half level pressure
195    REAL, DIMENSION(len), INTENT(IN)                 :: Ale1 ! Available lifting Energy
196    REAL, DIMENSION(len), INTENT(IN)                 :: Alp1 ! Available lifting Power
197    REAL, DIMENSION(len, nd), INTENT(IN)             :: omega1
198    REAL, INTENT(IN)                                  :: sig1feed1 ! sigma coord/pressure at lower bound of feeding layer
199    REAL, INTENT(IN)                                  :: sig2feed1 ! sigma coord/pressure at upper bound of feeding layer
200    REAL, DIMENSION(nd), INTENT(IN)                  :: wght1     ! weight density determining the feeding mixture
201
202! Input/Output
203    REAL, DIMENSION(len, nd), INTENT(INOUT)          :: sig1 ! section adiabatic updraft
204    REAL, DIMENSION(len, nd), INTENT(INOUT)          :: w01 ! vertical velocity within adiab updraft
205
206! Output
207    INTEGER, DIMENSION(len), INTENT(OUT)             :: iflag1 ! flag for Emanuel conditions
208    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ft1 ! temp tend
209    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fq1 ! spec hum tend
210    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fu1 ! u-wind tend
211    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fv1 ! v-wind tend
212    REAL, DIMENSION(len, nd, ntra), INTENT(OUT)      :: ftra1 ! tracor tend
213    REAL, DIMENSION(len), INTENT(OUT)                :: precip1 ! precipitation
214    INTEGER, DIMENSION(len), INTENT(OUT)             :: kbas1 ! cloud base level
215    INTEGER, DIMENSION(len), INTENT(OUT)             :: ktop1 ! cloud top level
216    REAL, DIMENSION(len), INTENT(OUT)                :: cbmf1 ! cloud base mass flux
217    REAL, DIMENSION(len), INTENT(OUT)                :: plcl1
218    REAL, DIMENSION(len), INTENT(OUT)                :: plfc1
219    REAL, DIMENSION(len), INTENT(OUT)                :: wbeff1
220    REAL, DIMENSION(len), INTENT(OUT)                :: ptop21 ! top of entraining zone
221    REAL, DIMENSION(len), INTENT(OUT)                :: sigd1
222    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ma1        ! mass flux adiabatic updraft (staggered grid)
223    REAL, DIMENSION(len, nd), INTENT(OUT)            :: mip1       ! mass flux shed by the adiabatic updraft (extensive)
224    REAL, DIMENSION(len, ndp1), INTENT(OUT)          :: vprecip1   ! vertical profile of total precipitation (staggered grid)
225    REAL, DIMENSION(len, ndp1), INTENT(OUT)          :: vprecipi1  ! vertical profile of ice precipitation (staggered grid)
226    REAL, DIMENSION(len, nd), INTENT(OUT)            :: upwd1      ! total upward mass flux (adiab+mixed) (staggered grid)
227    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dnwd1      ! saturated downward mass flux (mixed) (staggered grid)
228    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dnwd01     ! unsaturated downward mass flux (staggered grid)
229    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qcondc1    ! in-cld mixing ratio of condensed water (intensive)
230    REAL, DIMENSION(len), INTENT(OUT)                :: wd1        ! downdraft velocity scale for sfc fluxes
231    REAL, DIMENSION(len), INTENT(OUT)                :: cape1
232    REAL, DIMENSION(len), INTENT(OUT)                :: cin1
233    REAL, DIMENSION(len, nd), INTENT(OUT)            :: tvp1       ! adiab lifted parcell virt temp
234    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ftd1       ! temperature tendency due to precipitations (K/s)
235    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fqd1       ! specific humidity tendencies due to precipitations ((gm/gm)/s)
236    REAL, DIMENSION(len), INTENT(OUT)                :: Plim11
237    REAL, DIMENSION(len), INTENT(OUT)                :: Plim21
238    REAL, DIMENSION(len, nd), INTENT(OUT)            :: asupmax1   ! Highest mixing fraction of mixed updraughts
239    REAL, DIMENSION(len), INTENT(OUT)                :: supmax01
240    REAL, DIMENSION(len), INTENT(OUT)                :: asupmaxmin1
241    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qtc1    ! in cloud water content / specific humidity in convection (intensive)
242    REAL, DIMENSION(len, nd), INTENT(OUT)            :: sigt1   ! surface fraction in adiabatic updrafts / fract. cloud area (intensive)
243    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainA1 ! precipitation ejected from adiabatic draught
244    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainS1 ! precipitation detrained from shedding of adiabatic draught
245    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainM1 ! precipitation detrained from mixed draughts
246    REAL, DIMENSION(len, nd), INTENT(OUT)            :: mp1  ! unsat. mass flux (staggered grid)
247    REAL, DIMENSION(len, nd), INTENT(OUT)            :: da1  ! detrained mass flux of adiab. asc. air (extensive)
248    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: phi1 ! mass flux of envt. air in mixed draughts (extensive)
249    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: epmlmMm1  ! (extensive)
250    REAL, DIMENSION(len, nd), INTENT(OUT)            :: eplaMm1   ! (extensive)
251    REAL, DIMENSION(len, nd), INTENT(OUT)            :: evap1 ! evaporation rate in precip. downdraft. (intensive)
252    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ep1
253    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: sigij1 ! mass fraction of env. air in mixed draughts (intensive)
254    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: elij1! cond. water per unit mass of mixed draughts (intensive)
255    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qta1 ! total water per unit mass of the adiab. asc. (intensive)
256    REAL, DIMENSION(len, nd), INTENT(OUT)            :: clw1 ! condensed water content of the adiabatic updraught / cond. water per unit mass of the adiab. asc. (intensive)
257    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wghti1   ! final weight of the feeding layers (extensive)
258    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: phi21    ! (extensive)
259    REAL, DIMENSION(len, nd), INTENT(OUT)            :: d1a1     ! (extensive)
260    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dam1     ! (extensive)
261    REAL, DIMENSION(len), INTENT(OUT)               :: epmax_diag1
262
263! Local (non compressed) arrays
264    INTEGER nk1(len)
265    INTEGER icbs1(len)
266
267    LOGICAL, SAVE :: debut = .TRUE.
268!$omp THREADPRIVATE(debut)
269
270    REAL tnk1(len)
271    REAL qnk1(len)
272    REAL gznk1(len)
273    REAL unk1(len)
274    REAL vnk1(len)
275    REAL hnk1(len)
276    REAL pbase1(len)
277    REAL buoybase1(len)
278
279    REAL lf1(len, nd), lf1_wake(len, nd)
280    REAL lv1(len, nd), lv1_wake(len, nd)
281    REAL cpn1(len, nd), cpn1_wake(len, nd)
282    REAL tv1(len, nd), tv1_wake(len, nd)
283    REAL gz1(len, nd), gz1_wake(len, nd)
284    REAL h1(len, nd), h1_wake(len, nd)
285    REAL tp1(len, nd)
286    REAL th1(len, nd), th1_wake(len, nd)
287
288    CHARACTER(LEN=20) :: modname = 'cva_driver'
289    CHARACTER(LEN=80) :: abort_message
290
291    INTEGER :: coef_convective(len)
292    INTEGER :: k
293    REAL, parameter :: Cin_noconv = -100000., Cape_noconv = -1
294
295    INTEGER, SAVE :: igout = 1
296!$omp THREADPRIVATE(igout)
297
298    INTEGER nloc ! Allocation size for compressed arrays
299    logical :: compress
300
301    call enter_profile("cv3a_driver")
302
303    if (iflag_con /= 3) call abort_physic("cv3a_driver", "iflag_con must be 3", 1)
304
305! -------------------------------------------------------------------
306! --- SET CONSTANTS AND PARAMETERS
307! -------------------------------------------------------------------
308! -- set simulation flags:
309! (common cvflag
310    CALL cv_flag(iflag_ice_thermo)
311! -- set thermodynamical constants:
312! (common cvthermo)
313    CALL cv_thermo(iflag_con)
314! -- set convect parameters
315! includes microphysical parameters and parameters that
316! control the rate of approach to quasi-equilibrium)
317! (common cvparam)
318    CALL cv3_param(nd, k_upper, delt)
319
320    CALL cv3_incrcount(len, nd, delt, sig1)
321
322! --------------------------------------------------------------------
323! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
324! --------------------------------------------------------------------
325
326    IF (debut) PRINT *, 'Emanuel version 3 nouvelle'
327
328    block
329      REAL thnk1(len)
330      REAL qsnk1(len)
331      REAL cpnk1(len)
332      REAL hm1(len, nd)
333      REAL bid(len, nd) ! dummy array
334      REAL p1feed1(len) ! pressure at lower bound of feeding layer
335      REAL p2feed1(len) ! pressure at upper bound of feeding layer
336
337      integer :: i, icbmax
338
339      call driver_log('cv3_prelim')
340      CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, &
341                      lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1)
342
343      call driver_log('cv3_prelim')
344      CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, &
345                      lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, &
346                      h1_wake, bid, th1_wake)
347
348! --------------------------------------------------------------------
349! --- CONVECTIVE FEED
350! --------------------------------------------------------------------
351! compute feeding layer potential temperature and mixing ratio :
352! get bounds of feeding layer
353! test niveaux couche alimentation KE
354      IF (sig1feed1 == sig2feed1) THEN
355        WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed'
356        WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def'
357        abort_message = ''
358        CALL abort_physic(modname, abort_message, 1)
359      END IF
360
361      DO i = 1, len
362        p1feed1(i) = sig1feed1*ph1(i, 1)
363        p2feed1(i) = sig2feed1*ph1(i, 1)
364      END DO
365
366      call driver_log('cv3_feed')
367
368      ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed
369      iflag1(:) = 0
370      plcl1(:) = 0.
371      wghti1(:, :) = 0.
372      CALL cv3_feed(len, nd, ok_conserv_q, &
373                    t1, q1, u1, v1, p1, ph1, h1, gz1, &
374                    p1feed1, p2feed1, wght1, &
375                    wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, &
376                    cpnk1, hnk1, nk1, kbas1, icbmax, iflag1, gznk1, plcl1)
377
378! --------------------------------------------------------------------
379! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
380! (up through ICB for convect4, up through ICB+1 for convect3)
381! Calculates the lifted parcel virtual temperature at nk, the
382! actual temperature, and the adiabatic liquid water content.
383! --------------------------------------------------------------------
384      call driver_log('cv3_undilute1')
385      ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed
386      tvp1(:, :) = 0.
387      clw1(:, :) = 0.
388      CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, kbas1, tnk1, qnk1, & ! nd->na
389                         gznk1, tp1, tvp1, clw1, icbs1)
390
391! -------------------------------------------------------------------
392! --- TRIGGERING
393! -------------------------------------------------------------------
394
395      call driver_log('cv3_trigger')
396      CALL cv3_trigger(len, nd, kbas1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na
397                       pbase1, buoybase1, iflag1, sig1, w01)
398    end block
399
400! =====================================================================
401! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
402! =====================================================================
403
404! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
405! --- COMPRESS THE FIELDS
406!       (-> vectorization over convective gridpoints)
407! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
408
409    ! Compression has 3 possible modes:
410    ! 1) Compress = true :
411    compress_mode = COMPRESS_MODE_COMPRESS ! 1.1) Copy convective cells in contiguous array
412    !compress_mode = COMPRESS_MODE_COPY ! 1.2) Copy all cells and also compute on non-convective cells
413    ! 2) Compress = false : Don't copy and use original arrays, compute on non-convective cells
414    ! Never compress when comp_threshold = 0, always compress when comp_threshold = 1
415    ! comp_threshold = 0 (default) offers best performance when using many CPUs
416    compress = (comp_threshold == 1) .or. ( get_compress_size(len, (iflag1(:) == 0)) < len*comp_threshold )
417
418    if (compress) then
419      nloc = get_compress_size(len, (iflag1(:) == 0))
420
421      block
422        ! (local) compressed fields:
423        INTEGER iflag(nloc), nk(nloc), icb(nloc)
424        INTEGER nent(nloc, nd)
425        INTEGER icbs(nloc)
426        INTEGER inb(nloc)
427
428        REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc)
429        REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd)
430        REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd)
431        REAL s_wake(nloc)
432        REAL u(nloc, nd), v(nloc, nd)
433        REAL gz(nloc, nd), h(nloc, nd)
434        REAL h_wake(nloc, nd)
435        REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd)
436        REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd)
437        REAL p(nloc, nd), ph(nloc, nd + 1), tv(nloc, nd), tp(nloc, nd)
438        REAL tv_wake(nloc, nd)
439        REAL clw(nloc, nd)
440        REAL qta(nloc, nd), qpreca(nloc, nd)
441        REAL pbase(nloc), buoybase(nloc), th(nloc, nd)
442        REAL th_wake(nloc, nd)
443        REAL tvp(nloc, nd)
444        REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc)
445        REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd)
446        REAL buoy(nloc, nd)
447        REAL cape(nloc)
448        REAL cin(nloc)
449        REAL m(nloc, nd)
450        REAL mm(nloc, nd)
451        REAL ment(nloc, nd, nd), sigij(nloc, nd, nd)
452        REAL qent(nloc, nd, nd)
453        REAL hent(nloc, nd, nd)
454        REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
455        REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
456        REAL elij(nloc, nd, nd)
457        REAL supmax(nloc, nd)
458        REAL Ale(nloc), Alp(nloc), coef_clos(nloc)
459        REAL omega(nloc, nd)
460        REAL sigd(nloc)
461        REAL, DIMENSION(nloc, nd)     :: mp, qp, up, vp
462        REAL, DIMENSION(nloc, nd)     :: wt, water, evap
463        REAL, DIMENSION(nloc, nd)     :: ice, fondue, b
464        REAL, DIMENSION(nloc, nd)     :: frac_a, frac_s, faci
465        REAL ft(nloc, nd), fq(nloc, nd)
466        REAL ftd(nloc, nd), fqd(nloc, nd)
467        REAL fu(nloc, nd), fv(nloc, nd)
468        REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd)
469        REAL ma(nloc, nd), mip(nloc, nd)
470        REAL precip(nloc)
471        REAL vprecip(nloc, nd + 1)
472        REAL vprecipi(nloc, nd + 1)
473        REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
474        REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
475        REAL qcondc(nloc, nd)
476        REAL wd(nloc)
477        REAL Plim1(nloc), plim2(nloc)
478        REAL asupmax(nloc, nd)
479        REAL supmax0(nloc)
480        REAL asupmaxmin(nloc)
481
482        REAL tnk(nloc), qnk(nloc), gznk(nloc)
483        REAL wghti(nloc, nd)
484        REAL hnk(nloc), unk(nloc), vnk(nloc)
485
486        REAL qtc(nloc, nd)
487        REAL sigt(nloc, nd)
488
489        REAL wdtrainA(nloc, nd), wdtrainS(nloc, nd), wdtrainM(nloc, nd)
490        REAL da(nloc, nd), phi(nloc, nd, nd)
491        REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd)
492        REAL phi2(nloc, nd, nd)
493        REAL d1a(nloc, nd), dam(nloc, nd)
494        REAL epmax_diag(nloc)
495
496        LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
497
498        integer :: ncum ! Number of convective cells
499        type(compress_data_t) :: compress_data
500        type(array_list) :: cv3a_compress_list, cv3a_uncompress_list
501        logical, save :: timers_first = .true.
502        !$omp threadprivate(timers_first)
503
504        call enter_profile("cv3a_compress")
505        call driver_log('cv3a_compress')
506
507        call add_array_i1(cv3a_compress_list, iflag1, iflag)
508        call add_array_i1(cv3a_compress_list, nk1, nk)
509        call add_array_i1(cv3a_compress_list, kbas1, icb)
510        call add_array_i1(cv3a_compress_list, icbs1, icbs)
511        call add_array_r1(cv3a_compress_list, plcl1, plcl)
512        call add_array_r1(cv3a_compress_list, tnk1, tnk)
513        call add_array_r1(cv3a_compress_list, qnk1, qnk)
514        call add_array_r1(cv3a_compress_list, gznk1, gznk)
515        call add_array_r1(cv3a_compress_list, hnk1, hnk)
516        call add_array_r1(cv3a_compress_list, unk1, unk)
517        call add_array_r1(cv3a_compress_list, vnk1, vnk)
518        call add_array_r2(cv3a_compress_list, wghti1, wghti)
519        call add_array_r1(cv3a_compress_list, pbase1, pbase)
520        call add_array_r1(cv3a_compress_list, buoybase1, buoybase)
521        call add_array_r2(cv3a_compress_list, th1, th)
522        call add_array_r2(cv3a_compress_list, t1, t)
523        call add_array_r2(cv3a_compress_list, q1, q)
524        call add_array_r2(cv3a_compress_list, qs1, qs)
525        call add_array_r2(cv3a_compress_list, t1_wake, t_wake)
526        call add_array_r2(cv3a_compress_list, q1_wake, q_wake)
527        call add_array_r2(cv3a_compress_list, qs1_wake, qs_wake)
528        call add_array_r1(cv3a_compress_list, s1_wake, s_wake)
529        call add_array_r2(cv3a_compress_list, u1, u)
530        call add_array_r2(cv3a_compress_list, v1, v)
531        call add_array_r2(cv3a_compress_list, gz1, gz)
532        call add_array_r2(cv3a_compress_list, h1, h)
533        call add_array_r2(cv3a_compress_list, th1_wake, th_wake)
534        call add_array_r2(cv3a_compress_list, lv1, lv)
535        call add_array_r2(cv3a_compress_list, lf1, lf)
536        call add_array_r2(cv3a_compress_list, cpn1, cpn)
537        call add_array_r2(cv3a_compress_list, p1, p)
538        call add_array_r2(cv3a_compress_list, ph1, ph)
539        call add_array_r2(cv3a_compress_list, tv1, tv)
540        call add_array_r2(cv3a_compress_list, tp1, tp)
541        call add_array_r2(cv3a_compress_list, tvp1, tvp)
542        call add_array_r2(cv3a_compress_list, clw1, clw)
543        call add_array_r2(cv3a_compress_list, h1_wake, h_wake)
544        call add_array_r2(cv3a_compress_list, lv1_wake, lv_wake)
545        call add_array_r2(cv3a_compress_list, lf1_wake, lf_wake)
546        call add_array_r2(cv3a_compress_list, cpn1_wake, cpn_wake)
547        call add_array_r2(cv3a_compress_list, tv1_wake, tv_wake)
548        call add_array_r2(cv3a_compress_list, sig1, sig)
549        call add_array_r1(cv3a_compress_list, sig1(:, nd), sig(:, nd))
550        call add_array_r2(cv3a_compress_list, w01, w0)
551        call add_array_r1(cv3a_compress_list, Ale1, Ale)
552        call add_array_r1(cv3a_compress_list, Alp1, Alp)
553        call add_array_r2(cv3a_compress_list, omega1, omega)
554
555        call cv3a_compress(len, (iflag1 == 0), cv3a_compress_list, compress_data)
556        ncum = compress_data%ncum
557        call exit_profile("cv3a_compress")
558
559        call cv3a_driver_compressed(nloc, nd, ntra, &
560                                    debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
561                                    ! Input fields
562                                    delt, tau_cld_cv, coefw_cld_cv, nk, icbs, tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp, wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega, ph, &
563                                    ! Input/output
564                                    iflag, icb, plcl, sig, w0, tvp, clw, &
565                                    ! Output fields
566                                    inb, precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag, ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta,  evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt, vprecip, vprecipi, phi, phi2, sigij, elij, epmlmMm)
567
568        ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
569        ! --- UNCOMPRESS THE FIELDS
570        ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
571        call enter_profile("cv3a_uncompress")
572        call driver_log('cv3a_uncompress')
573        call add_array_i1(cv3a_uncompress_list, iflag, iflag1, init=.false.)
574        call add_array_i1(cv3a_uncompress_list, icb, kbas1)
575        call add_array_i1(cv3a_uncompress_list, inb, ktop1)
576        call add_array_r1(cv3a_uncompress_list, precip, precip1)
577        call add_array_r1(cv3a_uncompress_list, cbmf, cbmf1)
578        call add_array_r1(cv3a_uncompress_list, plcl, plcl1, init=.false.)
579        call add_array_r1(cv3a_uncompress_list, plfc, plfc1)
580        call add_array_r1(cv3a_uncompress_list, wbeff, wbeff1)
581        call add_array_r2(cv3a_uncompress_list, sig, sig1)
582        call add_array_r2(cv3a_uncompress_list, w0, w01)
583        call add_array_r1(cv3a_uncompress_list, ptop2, ptop21)
584        call add_array_r2(cv3a_uncompress_list, ft, ft1)
585        call add_array_r2(cv3a_uncompress_list, fq, fq1)
586        call add_array_r2(cv3a_uncompress_list, fu, fu1)
587        call add_array_r2(cv3a_uncompress_list, fv, fv1)
588        call add_array_r1(cv3a_uncompress_list, sigd, sigd1)
589        call add_array_r2(cv3a_uncompress_list, ma, ma1)
590        call add_array_r2(cv3a_uncompress_list, mip, mip1)
591        call add_array_r2(cv3a_uncompress_list, vprecip, vprecip1)
592        call add_array_r2(cv3a_uncompress_list, vprecipi, vprecipi1)
593        call add_array_r2(cv3a_uncompress_list, upwd, upwd1)
594        call add_array_r2(cv3a_uncompress_list, dnwd, dnwd1)
595        call add_array_r2(cv3a_uncompress_list, dnwd0, dnwd01)
596        call add_array_r2(cv3a_uncompress_list, qcondc, qcondc1)
597        call add_array_r1(cv3a_uncompress_list, wd, wd1)
598        cape1(:) = Cape_noconv
599        call add_array_r1(cv3a_uncompress_list, cape, cape1, init=.false.)
600        cin1(:) = Cin_noconv
601        call add_array_r1(cv3a_uncompress_list, cin, cin1, init=.false.)
602        call add_array_r2(cv3a_uncompress_list, tvp, tvp1, init=.false.)
603        call add_array_r2(cv3a_uncompress_list, ftd, ftd1)
604        call add_array_r2(cv3a_uncompress_list, fqd, fqd1)
605        call add_array_r1(cv3a_uncompress_list, Plim1, Plim11)
606        call add_array_r1(cv3a_uncompress_list, plim2, plim21)
607        call add_array_r2(cv3a_uncompress_list, asupmax, asupmax1)
608        call add_array_r1(cv3a_uncompress_list, supmax0, supmax01)
609        call add_array_r1(cv3a_uncompress_list, asupmaxmin, asupmaxmin1)
610        call add_array_r2(cv3a_uncompress_list, da, da1)
611        call add_array_r3(cv3a_uncompress_list, phi, phi1)
612        call add_array_r2(cv3a_uncompress_list, mp, mp1)
613        call add_array_r3(cv3a_uncompress_list, phi2, phi21)
614        call add_array_r2(cv3a_uncompress_list, d1a, d1a1)
615        call add_array_r2(cv3a_uncompress_list, dam, dam1)
616        call add_array_r3(cv3a_uncompress_list, sigij, sigij1)
617        call add_array_r2(cv3a_uncompress_list, qta, qta1)
618        call add_array_r2(cv3a_uncompress_list, clw, clw1, init=.false.)
619        call add_array_r3(cv3a_uncompress_list, elij, elij1)
620        call add_array_r2(cv3a_uncompress_list, evap, evap1)
621        call add_array_r2(cv3a_uncompress_list, ep, ep1)
622        call add_array_r3(cv3a_uncompress_list, epmlmMm, epmlmMm1)
623        call add_array_r2(cv3a_uncompress_list, eplaMm, eplaMm1)
624        call add_array_r2(cv3a_uncompress_list, wdtrainA, wdtrainA1)
625        call add_array_r2(cv3a_uncompress_list, wdtrainS, wdtrainS1)
626        call add_array_r2(cv3a_uncompress_list, wdtrainM, wdtrainM1)
627        call add_array_r2(cv3a_uncompress_list, qtc, qtc1)
628        call add_array_r2(cv3a_uncompress_list, sigt, sigt1)
629        call add_array_r1(cv3a_uncompress_list, epmax_diag, epmax_diag1)
630        call add_array_r1(cv3a_uncompress_list, sig(:, nd), sig1(:, nd))
631        call cv3a_uncompress(len, compress_data, cv3a_uncompress_list)
632        call exit_profile("cv3a_uncompress")
633
634      end block
635    else !compress
636      nloc = len
637      coef_convective(:) = merge(1, 0, iflag1(:) == 0)
638      call cv3a_driver_compressed(nloc, nd, ntra, &
639                                  debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
640                                  ! Input fields
641                                  delt, tau_cld_cv, coefw_cld_cv, nk1, icbs1, tnk1, qnk1, gznk1, hnk1, unk1, vnk1, pbase1, buoybase1, s1_wake, Ale1, Alp1, wghti1, th1, t1, q1, qs1, t1_wake, q1_wake, qs1_wake, u1, v1, gz1, h1, th1_wake, lv1, lf1, cpn1, p1, tv1, tp1, h1_wake, lv1_wake, lf1_wake, cpn1_wake, tv1_wake, omega1, ph1, &
642                                  ! Input/output
643                                  iflag1, kbas1, plcl1, sig1, w01, tvp1, clw1, &
644                                  ! Output fields
645                                  ktop1, precip1, cbmf1, plfc1, wbeff1, ptop21, sigd1, wd1, cape1, cin1, Plim11, plim21, supmax01, asupmaxmin1, epmax_diag1, ft1, fq1, fu1, fv1, ma1, mip1, upwd1, dnwd1, dnwd01, qcondc1, ftd1, fqd1, asupmax1, da1, mp1, d1a1, dam1, qta1,  evap1, ep1, eplaMm1, wdtrainA1, wdtrainS1, wdtrainM1, qtc1, sigt1, vprecip1, vprecipi1, phi1, phi21, sigij1, elij1, epmlmMm1)
646      !  Reset Cin1 and Cape1 to their default values for non-convective grid-points.
647      Cin1(:) = Cin1(:)*coef_convective(:) + Cin_noconv*(1.-coef_convective(:))
648      Cape1(:) = Cape1(:)*coef_convective(:) + Cape_noconv*(1.-coef_convective(:))
649      precip1(:) = precip1(:)*coef_convective(:)
650      kbas1(:) = kbas1(:)*coef_convective(:)
651      ktop1(:) = ktop1(:)*coef_convective(:)
652      sigd1(:) = sigd1(:)*coef_convective(:)
653      wbeff1(:) = wbeff1(:)*coef_convective(:)
654      DO k = 1, nd
655        qcondc1(:, k) = qcondc1(:, k)*coef_convective(:)
656        vprecip1(:, k) = vprecip1(:, k)*coef_convective(:)
657        vprecipi1(:, k) = vprecipi1(:, k)*coef_convective(:)
658      ENDDO
659    endif
660
661    IF (prt_level >= 10) THEN
662      Print *, 'cva_driver after cv3_uncompress:ft1(1) , ftd1(1) ', &
663        ft1(igout, 1), ftd1(igout, 1)
664      Print *, 'cva_driver after cv3_uncompress:fq1(1) , fqd1(1) ', &
665        fq1(igout, 1), fqd1(igout, 1)
666    ENDIF
667    IF (debut) THEN
668      PRINT *, ' cv_uncompress -> '
669      debut = .FALSE.
670    END IF
671
672    ftra1(:, :, :) = 0
673
674    call exit_profile("cv3a_driver")
675  END SUBROUTINE
676
677  subroutine cv3a_driver_compressed(nloc, nd, ntra, &
678                                    debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
679                                    ! Input fields
680                                    delt, tau_cld_cv, coefw_cld_cv, nk, icbs, tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp, wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega, ph, &
681                                    ! Input/output
682                                    iflag, icb, plcl, sig, w0, tvp, clw, &
683                                    ! Output fields
684                                    inb, precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag, ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta,  evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt, vprecip, vprecipi, phi, phi2, sigij, elij, epmlmMm)
685    use profiling_physic_mod
686    USE cv3p_mixing_mod, ONLY : cv3p_mixing
687
688    integer, intent(in) :: nloc, nd, ntra ! Arrays sizes
689    logical, intent(in) :: debut, ok_conserv_q
690    integer, intent(in) :: iflag_mix, iflag_clos, prt_level
691    real, intent(in)    :: delt, tau_cld_cv, coefw_cld_cv
692
693    ! Input arrays
694    integer, dimension(nloc), intent(in) :: nk, icbs
695    real, dimension(nloc), intent(in) :: tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp
696    real, dimension(nloc, nd), intent(in) :: wghti, th, t, q, qs, t_wake, q_wake, qs_wake, u, v, gz, h, th_wake, lv, lf, cpn, p, tv, tp, h_wake, lv_wake, lf_wake, cpn_wake, tv_wake, omega
697    real, dimension(nloc, nd + 1), intent(in) :: ph
698
699    ! Input/output arrays
700    integer, dimension(nloc), intent(inout) :: iflag, icb
701    real, dimension(nloc), intent(inout) :: plcl
702    real, dimension(nloc, nd), intent(inout) :: sig, w0, tvp, clw
703
704    ! Output arrays
705    integer, dimension(nloc), intent(out) :: inb
706    real, dimension(nloc), intent(out) :: precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag
707    real, dimension(nloc, nd), intent(out) :: ft, fq, fu, fv, ma, mip, upwd, dnwd, dnwd0, qcondc, ftd, fqd, asupmax, da, mp, d1a, dam, qta,  evap, ep, eplaMm, wdtrainA, wdtrainS, wdtrainM, qtc, sigt
708    real, dimension(nloc, nd + 1), intent(out) :: vprecip, vprecipi
709    real, dimension(nloc, nd, nd), intent(out) :: phi, phi2, sigij, elij, epmlmMm
710
711    ! Local fields
712    INTEGER nent(nloc, nd)
713    REAL hp(nloc, nd), sigp(nloc, nd), qpreca(nloc, nd)
714    REAL buoy(nloc, nd)
715    REAL m(nloc, nd)
716    REAL mm(nloc, nd)
717    REAL ment(nloc, nd, nd)
718    REAL qent(nloc, nd, nd)
719    REAL hent(nloc, nd, nd)
720    REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
721    REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
722    REAL supmax(nloc, nd)
723    REAL coef_clos(nloc)
724    REAL, DIMENSION(nloc, nd)     :: qp, up, vp
725    REAL, DIMENSION(nloc, nd)     :: wt, water
726    REAL, DIMENSION(nloc, nd)     :: ice, fondue, b
727    REAL, DIMENSION(nloc, nd)     :: frac_a, frac_s, faci
728    REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
729    REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
730
731    LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
732    integer :: k
733
734    integer :: ncum ! Number of convective cells
735    logical, save :: timers_first = .true.
736    !$omp threadprivate(timers_first)
737    INTEGER, PARAMETER :: igout = 1
738
739    ncum = nloc
740    call enter_profile("cv3a_compresd")
741    IF (ncum > 0) THEN
742
743! -------------------------------------------------------------------
744! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
745! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
746! ---   &
747! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
748! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
749! ---   &
750! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
751! -------------------------------------------------------------------
752
753      call driver_log('cv3_undilute2')
754      qta(:,:) = 0
755      ep(:,:) = 0
756      CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, &
757                         tnk, qnk, gznk, hnk, t, q, qs, gz, &
758                         p, ph, h, tv, lv, lf, pbase, buoybase, plcl, &
759                         inb, tp, tvp, clw, hp, ep, sigp, buoy, &
760                         frac_a, frac_s, qpreca, qta)     
761     
762! epmax_cape
763! on recalcule ep et hp
764      call driver_log('cv3_epmax_cape')
765      epmax_diag(:) = 0
766      call cv3_epmax_fn_cape(nloc, ncum, nd &
767                             , ep, hp, icb, inb, clw, nk, t, h, hnk, lv, lf, frac_s &
768                             , pbase, p, ph, tv, buoy, sig, w0, iflag &
769                             , epmax_diag)
770
771! -------------------------------------------------------------------
772! --- MIXING(1)   (if iflag_mix .ge. 1)
773! -------------------------------------------------------------------
774      elij(:,:,:) = 0
775      supmax(:,:) = 0
776      IF (iflag_mix >= 1) THEN
777        call driver_log('cv3p_mixing')
778        CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &
779                         ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qta, &
780                         unk, vnk, hp, tv, tvp, ep, clw, sig, &
781                         ment, qent, hent, uent, vent, nent, &
782                         sigij, elij, supmax, ments, qents, traent)
783      END IF
784
785! -------------------------------------------------------------------
786! --- CLOSURE
787! -------------------------------------------------------------------
788
789      cape(:) = -1
790      cin(:) = -100000.
791      ptop2(:) = 0
792      coef_clos(:) = 1.
793      Plim1(:) = 0
794      plim2(:) = 0
795      supmax0(:) = 0
796      asupmaxmin(:) = 0
797      cbmf(:) = 0
798      plfc(:) = 0
799      wbeff(:) = 0
800      asupmax(:,:) = 0
801
802      ok_inhib = (iflag_mix == 2)
803     
804      IF (iflag_clos == 0) THEN
805        call driver_log('cv3_closure')
806        CALL cv3_closure(nloc, ncum, nd, icb, inb, &
807                         pbase, p, ph, tv, buoy, &
808                         sig, w0, cape, m, iflag)
809      END IF   ! iflag_clos==0     
810
811      IF (iflag_clos == 1) PRINT *, ' pas d appel cv3p_closure'
812
813      IF (iflag_clos == 2) THEN
814        call driver_log('cv3p1_closure')
815        CALL cv3p1_closure(nloc, ncum, nd, icb, inb, &
816                           pbase, plcl, p, ph, tv, tvp, buoy, &
817                           supmax, ok_inhib, Ale, Alp, omega, &
818                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
819                           Plim1, plim2, asupmax, supmax0, &
820                           asupmaxmin, cbmf, plfc, wbeff)
821        if (prt_level >= 10) PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1)
822      END IF ! iflag_clos==2
823
824      IF (iflag_clos == 3) THEN
825        call driver_log('cv3p2_closure')
826        CALL cv3p2_closure(nloc, ncum, nd, icb, inb, &
827                           pbase, plcl, p, ph, tv, tvp, buoy, &
828                           supmax, ok_inhib, Ale, Alp, omega, &
829                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
830                           Plim1, plim2, asupmax, supmax0, &
831                           asupmaxmin, cbmf, plfc, wbeff)
832        if (prt_level >= 10) PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1)
833      END IF ! iflag_clos==3
834
835! -------------------------------------------------------------------
836! --- MIXING(2)
837! -------------------------------------------------------------------
838
839      IF (iflag_mix == 0) THEN
840        call driver_log('cv3_mixing')
841        CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &             ! na->nd
842                        ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qnk, &
843                        unk, vnk, hp, tv, tvp, ep, clw, m, sig, &
844                        ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent)
845        CALL zilch(hent, nloc*nd*nd)
846      ELSE
847        mm(:, :) = m(:, :)
848        CALL cv3_mixscale(nloc, ncum, nd, ment, mm)
849        IF (debut) PRINT *, ' cv3_mixscale-> '
850      END IF
851
852      IF (debut) PRINT *, ' cv_mixing ->'
853
854! -------------------------------------------------------------------
855! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
856! -------------------------------------------------------------------
857      IF (debut) PRINT *, ' cva_driver -> cv3_unsat '
858
859      call driver_log('cv3_unsat')
860      CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, &
861                     t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, &
862                     th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, &
863                     ep, sigp, clw, frac_s, qpreca, frac_a, qta, &
864                     m, ment, elij, delt, plcl, coef_clos, &
865                     mp, qp, up, vp, trap, wt, water, evap, fondue, ice, &
866                     faci, b, sigd, &
867                     wdtrainA, wdtrainS, wdtrainM)
868      IF (prt_level >= 10) THEN
869        Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue '
870        DO k = 1, nd
871          write (6, '(i4,5(1x,e13.6))'), &
872            k, mp(igout, k), water(igout, k), ice(igout, k), &
873            evap(igout, k), fondue(igout, k)
874        ENDDO
875        Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainS, wdtrainM '
876        DO k = 1, nd
877          write (6, '(i4,3(1x,e13.6))'), &
878            k, wdtrainA(igout, k), wdtrainS(igout, k), wdtrainM(igout, k)
879        ENDDO
880      ENDIF
881
882      IF (debut) PRINT *, 'cv_unsat-> '
883! -------------------------------------------------------------------
884! YIELD
885! (tendencies, precipitation, variables of interface with other processes, etc)
886! -------------------------------------------------------------------
887
888      call driver_log('cv3_yield')
889      CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, &
890                     icb, inb, delt, &
891                     t, q, t_wake, q_wake, s_wake, u, v, tra, &
892                     gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, &
893                     ep, clw, qpreca, m, tp, mp, qp, up, vp, trap, &
894                     wt, water, ice, evap, fondue, faci, b, sigd, &
895                     ment, qent, hent, iflag_mix, uent, vent, &
896                     nent, elij, traent, sig, &
897                     tv, tvp, wghti, &
898                     iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, &
899                     cbmf, upwd, dnwd, dnwd0, ma, mip, &
900                     qcondc, wd, &
901                     ftd, fqd, qta, qtc, sigt, tau_cld_cv, coefw_cld_cv)
902      ! Test conseravtion de l'eau
903      IF (debut) PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1)
904      IF (prt_level >= 10) THEN
905        Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', &
906          ft(igout, 1), ftd(igout, 1)
907        Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', &
908          fq(igout, 1), fqd(igout, 1)
909      ENDIF
910
911   
912
913!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
914!--- passive tracers
915!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
916
917      call driver_log('cv3_tracer')
918      sigij(:,:,:) = 0   
919      vprecip(:,:) = 0
920      ! GLITCHY : vprecip is unused in cv3_tracer and sigij is intent in
921      CALL cv3_tracer(nloc, -1, ncum, nd, nd, &
922                      ment, sigij, da, phi, phi2, d1a, dam, &
923                      ep, vprecip, elij, clw, epmlmMm, eplaMm, &
924                      icb, inb)
925
926      vprecipi = 0
927    END IF ! ncum>0
928    call exit_profile("cv3a_compresd")
929  end subroutine
930
931  subroutine driver_log(message)
932    use print_control_mod, only: prt_level
933    character(*) :: message
934    if (prt_level >= 9) PRINT *, 'cva_driver ->', message
935  end subroutine
936end module
Note: See TracBrowser for help on using the repository browser.