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

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

Added a path to execute cv3a_driver without compression at all

+ fix : Compute coef_convective BEFORE cv3a_driver_compressed

  • 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.6 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, &
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, DIMENSION(len, nd), INTENT(IN)             :: t1 ! temperature (sat draught envt)
184    REAL, DIMENSION(len, nd), INTENT(IN)             :: q1 ! specific hum (sat draught envt)
185    REAL, DIMENSION(len, nd), INTENT(IN)             :: qs1 ! sat specific hum (sat draught envt)
186    REAL, DIMENSION(len, nd), INTENT(IN)             :: t1_wake ! temperature (unsat draught envt)
187    REAL, DIMENSION(len, nd), INTENT(IN)             :: q1_wake ! specific hum(unsat draught envt)
188    REAL, DIMENSION(len, nd), INTENT(IN)             :: qs1_wake ! sat specific hum(unsat draughts envt)
189    REAL, DIMENSION(len), INTENT(IN)                 :: s1_wake ! fractionnal area covered by wakes
190    REAL, DIMENSION(len, nd), INTENT(IN)             :: u1 ! u-wind
191    REAL, DIMENSION(len, nd), INTENT(IN)             :: v1 ! v-wind
192    REAL, DIMENSION(len, nd), INTENT(IN)             :: p1 ! full level pressure
193    REAL, DIMENSION(len, ndp1), INTENT(IN)           :: ph1 ! half level pressure
194    REAL, DIMENSION(len), INTENT(IN)                 :: Ale1 ! Available lifting Energy
195    REAL, DIMENSION(len), INTENT(IN)                 :: Alp1 ! Available lifting Power
196    REAL, DIMENSION(len, nd), INTENT(IN)             :: omega1
197    REAL, INTENT(IN)                                  :: sig1feed1 ! sigma coord/pressure at lower bound of feeding layer
198    REAL, INTENT(IN)                                  :: sig2feed1 ! sigma coord/pressure at upper bound of feeding layer
199    REAL, DIMENSION(nd), INTENT(IN)                  :: wght1     ! weight density determining the feeding mixture
200
201! Input/Output
202    REAL, DIMENSION(len, nd), INTENT(INOUT)          :: sig1 ! section adiabatic updraft
203    REAL, DIMENSION(len, nd), INTENT(INOUT)          :: w01 ! vertical velocity within adiab updraft
204
205! Output
206    INTEGER, DIMENSION(len), INTENT(OUT)             :: iflag1 ! flag for Emanuel conditions
207    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ft1 ! temp tend
208    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fq1 ! spec hum tend
209    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fu1 ! u-wind tend
210    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fv1 ! v-wind tend
211    REAL, DIMENSION(len, nd, ntra), INTENT(OUT)      :: ftra1 ! tracor tend
212    REAL, DIMENSION(len), INTENT(OUT)                :: precip1 ! precipitation
213    INTEGER, DIMENSION(len), INTENT(OUT)             :: kbas1 ! cloud base level
214    INTEGER, DIMENSION(len), INTENT(OUT)             :: ktop1 ! cloud top level
215    REAL, DIMENSION(len), INTENT(OUT)                :: cbmf1 ! cloud base mass flux
216    REAL, DIMENSION(len), INTENT(OUT)                :: plcl1
217    REAL, DIMENSION(len), INTENT(OUT)                :: plfc1
218    REAL, DIMENSION(len), INTENT(OUT)                :: wbeff1
219    REAL, DIMENSION(len), INTENT(OUT)                :: ptop21 ! top of entraining zone
220    REAL, DIMENSION(len), INTENT(OUT)                :: sigd1
221    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ma1        ! mass flux adiabatic updraft (staggered grid)
222    REAL, DIMENSION(len, nd), INTENT(OUT)            :: mip1       ! mass flux shed by the adiabatic updraft (extensive)
223    REAL, DIMENSION(len, ndp1), INTENT(OUT)          :: vprecip1   ! vertical profile of total precipitation (staggered grid)
224    REAL, DIMENSION(len, ndp1), INTENT(OUT)          :: vprecipi1  ! vertical profile of ice precipitation (staggered grid)
225    REAL, DIMENSION(len, nd), INTENT(OUT)            :: upwd1      ! total upward mass flux (adiab+mixed) (staggered grid)
226    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dnwd1      ! saturated downward mass flux (mixed) (staggered grid)
227    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dnwd01     ! unsaturated downward mass flux (staggered grid)
228    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qcondc1    ! in-cld mixing ratio of condensed water (intensive)
229    REAL, DIMENSION(len), INTENT(OUT)                :: wd1        ! downdraft velocity scale for sfc fluxes
230    REAL, DIMENSION(len), INTENT(OUT)                :: cape1
231    REAL, DIMENSION(len), INTENT(OUT)                :: cin1
232    REAL, DIMENSION(len, nd), INTENT(OUT)            :: tvp1       ! adiab lifted parcell virt temp
233    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ftd1       ! temperature tendency due to precipitations (K/s)
234    REAL, DIMENSION(len, nd), INTENT(OUT)            :: fqd1       ! specific humidity tendencies due to precipitations ((gm/gm)/s)
235    REAL, DIMENSION(len), INTENT(OUT)                :: Plim11
236    REAL, DIMENSION(len), INTENT(OUT)                :: Plim21
237    REAL, DIMENSION(len, nd), INTENT(OUT)            :: asupmax1   ! Highest mixing fraction of mixed updraughts
238    REAL, DIMENSION(len), INTENT(OUT)                :: supmax01
239    REAL, DIMENSION(len), INTENT(OUT)                :: asupmaxmin1
240    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qtc1    ! in cloud water content / specific humidity in convection (intensive)
241    REAL, DIMENSION(len, nd), INTENT(OUT)            :: sigt1   ! surface fraction in adiabatic updrafts / fract. cloud area (intensive)
242    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainA1 ! precipitation ejected from adiabatic draught
243    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainS1 ! precipitation detrained from shedding of adiabatic draught
244    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wdtrainM1 ! precipitation detrained from mixed draughts
245    REAL, DIMENSION(len, nd), INTENT(OUT)            :: mp1  ! unsat. mass flux (staggered grid)
246    REAL, DIMENSION(len, nd), INTENT(OUT)            :: da1  ! detrained mass flux of adiab. asc. air (extensive)
247    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: phi1 ! mass flux of envt. air in mixed draughts (extensive)
248    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: epmlmMm1  ! (extensive)
249    REAL, DIMENSION(len, nd), INTENT(OUT)            :: eplaMm1   ! (extensive)
250    REAL, DIMENSION(len, nd), INTENT(OUT)            :: evap1 ! evaporation rate in precip. downdraft. (intensive)
251    REAL, DIMENSION(len, nd), INTENT(OUT)            :: ep1
252    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: sigij1 ! mass fraction of env. air in mixed draughts (intensive)
253    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: elij1! cond. water per unit mass of mixed draughts (intensive)
254    REAL, DIMENSION(len, nd), INTENT(OUT)            :: qta1 ! total water per unit mass of the adiab. asc. (intensive)
255    REAL, DIMENSION(len, nd), INTENT(OUT)            :: clw1 ! condensed water content of the adiabatic updraught / cond. water per unit mass of the adiab. asc. (intensive)
256    REAL, DIMENSION(len, nd), INTENT(OUT)            :: wghti1   ! final weight of the feeding layers (extensive)
257    REAL, DIMENSION(len, nd, nd), INTENT(OUT)        :: phi21    ! (extensive)
258    REAL, DIMENSION(len, nd), INTENT(OUT)            :: d1a1     ! (extensive)
259    REAL, DIMENSION(len, nd), INTENT(OUT)            :: dam1     ! (extensive)
260    REAL, DIMENSION(len), INTENT(OUT)               :: epmax_diag1
261
262! Local (non compressed) arrays
263    INTEGER nk1(len)
264    INTEGER icbs1(len)
265
266    LOGICAL, SAVE :: debut = .TRUE.
267!$omp THREADPRIVATE(debut)
268
269    REAL tnk1(len)
270    REAL qnk1(len)
271    REAL gznk1(len)
272    REAL unk1(len)
273    REAL vnk1(len)
274    REAL hnk1(len)
275    REAL pbase1(len)
276    REAL buoybase1(len)
277
278    REAL lf1(len, nd), lf1_wake(len, nd)
279    REAL lv1(len, nd), lv1_wake(len, nd)
280    REAL cpn1(len, nd), cpn1_wake(len, nd)
281    REAL tv1(len, nd), tv1_wake(len, nd)
282    REAL gz1(len, nd), gz1_wake(len, nd)
283    REAL h1(len, nd), h1_wake(len, nd)
284    REAL tp1(len, nd)
285    REAL th1(len, nd), th1_wake(len, nd)
286
287    CHARACTER(LEN=20) :: modname = 'cva_driver'
288    CHARACTER(LEN=80) :: abort_message
289
290    INTEGER :: coef_convective(len)
291    INTEGER :: k
292    REAL, parameter :: Cin_noconv = -100000., Cape_noconv = -1
293
294    INTEGER, SAVE :: igout = 1
295!$omp THREADPRIVATE(igout)
296
297    INTEGER nloc ! Allocation size for compressed arrays
298    logical :: compress
299
300    call enter_profile("cv3a_driver")
301
302    if (iflag_con /= 3) call abort_physic("cv3a_driver", "iflag_con must be 3", 1)
303
304    call enter_profile("cv3a_params")
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    call exit_profile("cv3a_params")
322
323! --------------------------------------------------------------------
324! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
325! --------------------------------------------------------------------
326
327    IF (debut) PRINT *, 'Emanuel version 3 nouvelle'
328
329    block
330      REAL thnk1(len)
331      REAL qsnk1(len)
332      REAL cpnk1(len)
333      REAL hm1(len, nd)
334      REAL bid(len, nd) ! dummy array
335      REAL p1feed1(len) ! pressure at lower bound of feeding layer
336      REAL p2feed1(len) ! pressure at upper bound of feeding layer
337
338      integer :: i, icbmax
339
340      call enter_profile("cv3a_uncompressed")
341
342      call driver_log('cv3_prelim')
343      CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, &
344                      lv1, lf1, cpn1, tv1, gz1, h1, hm1, th1)
345
346      call driver_log('cv3_prelim')
347      CALL cv3_prelim(len, nd, ndp1, t1_wake, q1_wake, p1, ph1, &
348                      lv1_wake, lf1_wake, cpn1_wake, tv1_wake, gz1_wake, &
349                      h1_wake, bid, th1_wake)
350
351! --------------------------------------------------------------------
352! --- CONVECTIVE FEED
353! --------------------------------------------------------------------
354! compute feeding layer potential temperature and mixing ratio :
355! get bounds of feeding layer
356! test niveaux couche alimentation KE
357      IF (sig1feed1 == sig2feed1) THEN
358        WRITE (lunout, *) 'impossible de choisir sig1feed=sig2feed'
359        WRITE (lunout, *) 'changer la valeur de sig2feed dans physiq.def'
360        abort_message = ''
361        CALL abort_physic(modname, abort_message, 1)
362      END IF
363
364      DO i = 1, len
365        p1feed1(i) = sig1feed1*ph1(i, 1)
366        p2feed1(i) = sig2feed1*ph1(i, 1)
367      END DO
368
369      call driver_log('cv3_feed')
370
371      ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed
372      iflag1(:) = 0
373      plcl1(:) = 0.
374      wghti1(:, :) = 0.
375      CALL cv3_feed(len, nd, ok_conserv_q, &
376                    t1, q1, u1, v1, p1, ph1, h1, gz1, &
377                    p1feed1, p2feed1, wght1, &
378                    wghti1, tnk1, thnk1, qnk1, qsnk1, unk1, vnk1, &
379                    cpnk1, hnk1, nk1, kbas1, icbmax, iflag1, gznk1, plcl1)
380
381! --------------------------------------------------------------------
382! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
383! (up through ICB for convect4, up through ICB+1 for convect3)
384! Calculates the lifted parcel virtual temperature at nk, the
385! actual temperature, and the adiabatic liquid water content.
386! --------------------------------------------------------------------
387      call driver_log('cv3_undilute1')
388      ! GLITCHY : arrays are set to zero but are intent(out) in call to cv3_feed
389      tvp1(:, :) = 0.
390      clw1(:, :) = 0.
391      CALL cv3_undilute1(len, nd, t1, qs1, gz1, plcl1, p1, kbas1, tnk1, qnk1, & ! nd->na
392                         gznk1, tp1, tvp1, clw1, icbs1)
393
394! -------------------------------------------------------------------
395! --- TRIGGERING
396! -------------------------------------------------------------------
397
398      call driver_log('cv3_trigger')
399      CALL cv3_trigger(len, nd, kbas1, plcl1, p1, th1, tv1, tvp1, thnk1, & ! nd->na
400                       pbase1, buoybase1, iflag1, sig1, w01)
401      call exit_profile("cv3a_uncompressed")
402    end block
403
404! =====================================================================
405! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
406! =====================================================================
407
408! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
409! --- COMPRESS THE FIELDS
410!       (-> vectorization over convective gridpoints)
411! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
412
413    compress = .false.
414    if (compress) then
415      compress_mode = COMPRESS_MODE_COMPRESS
416      !compress_mode = COMPRESS_MODE_COPY
417
418      nloc = get_compress_size(len, (iflag1(:) == 0))
419
420      block
421        ! (local) compressed fields:
422        INTEGER iflag(nloc), nk(nloc), icb(nloc)
423        INTEGER nent(nloc, nd)
424        INTEGER icbs(nloc)
425        INTEGER inb(nloc)
426
427        REAL cbmf(nloc), plcl(nloc), plfc(nloc), wbeff(nloc)
428        REAL t(nloc, nd), q(nloc, nd), qs(nloc, nd)
429        REAL t_wake(nloc, nd), q_wake(nloc, nd), qs_wake(nloc, nd)
430        REAL s_wake(nloc)
431        REAL u(nloc, nd), v(nloc, nd)
432        REAL gz(nloc, nd), h(nloc, nd)
433        REAL h_wake(nloc, nd)
434        REAL lv(nloc, nd), lf(nloc, nd), cpn(nloc, nd)
435        REAL lv_wake(nloc, nd), lf_wake(nloc, nd), cpn_wake(nloc, nd)
436        REAL p(nloc, nd), ph(nloc, nd + 1), tv(nloc, nd), tp(nloc, nd)
437        REAL tv_wake(nloc, nd)
438        REAL clw(nloc, nd)
439        REAL qta(nloc, nd), qpreca(nloc, nd)
440        REAL pbase(nloc), buoybase(nloc), th(nloc, nd)
441        REAL th_wake(nloc, nd)
442        REAL tvp(nloc, nd)
443        REAL sig(nloc, nd), w0(nloc, nd), ptop2(nloc)
444        REAL hp(nloc, nd), ep(nloc, nd), sigp(nloc, nd)
445        REAL buoy(nloc, nd)
446        REAL cape(nloc)
447        REAL cin(nloc)
448        REAL m(nloc, nd)
449        REAL mm(nloc, nd)
450        REAL ment(nloc, nd, nd), sigij(nloc, nd, nd)
451        REAL qent(nloc, nd, nd)
452        REAL hent(nloc, nd, nd)
453        REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
454        REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
455        REAL elij(nloc, nd, nd)
456        REAL supmax(nloc, nd)
457        REAL Ale(nloc), Alp(nloc), coef_clos(nloc)
458        REAL omega(nloc, nd)
459        REAL sigd(nloc)
460        REAL, DIMENSION(nloc, nd)     :: mp, qp, up, vp
461        REAL, DIMENSION(nloc, nd)     :: wt, water, evap
462        REAL, DIMENSION(nloc, nd)     :: ice, fondue, b
463        REAL, DIMENSION(nloc, nd)     :: frac_a, frac_s, faci
464        REAL ft(nloc, nd), fq(nloc, nd)
465        REAL ftd(nloc, nd), fqd(nloc, nd)
466        REAL fu(nloc, nd), fv(nloc, nd)
467        REAL upwd(nloc, nd), dnwd(nloc, nd), dnwd0(nloc, nd)
468        REAL ma(nloc, nd), mip(nloc, nd)
469        REAL precip(nloc)
470        REAL vprecip(nloc, nd + 1)
471        REAL vprecipi(nloc, nd + 1)
472        REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
473        REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
474        REAL qcondc(nloc, nd)
475        REAL wd(nloc)
476        REAL Plim1(nloc), plim2(nloc)
477        REAL asupmax(nloc, nd)
478        REAL supmax0(nloc)
479        REAL asupmaxmin(nloc)
480
481        REAL tnk(nloc), qnk(nloc), gznk(nloc)
482        REAL wghti(nloc, nd)
483        REAL hnk(nloc), unk(nloc), vnk(nloc)
484
485        REAL qtc(nloc, nd)
486        REAL sigt(nloc, nd)
487
488        REAL wdtrainA(nloc, nd), wdtrainS(nloc, nd), wdtrainM(nloc, nd)
489        REAL da(nloc, nd), phi(nloc, nd, nd)
490        REAL epmlmMm(nloc, nd, nd), eplaMm(nloc, nd)
491        REAL phi2(nloc, nd, nd)
492        REAL d1a(nloc, nd), dam(nloc, nd)
493        REAL epmax_diag(nloc)
494
495        LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
496
497        integer :: ncum ! Number of convective cells
498        type(compress_data_t) :: compress_data
499        type(array_list) :: cv3a_compress_list, cv3a_uncompress_list
500        logical, save :: timers_first = .true.
501        !$omp threadprivate(timers_first)
502
503        call enter_profile("cv3a_compress")
504        call driver_log('cv3a_compress')
505
506        call add_array_i1(cv3a_compress_list, iflag1, iflag)
507        call add_array_i1(cv3a_compress_list, nk1, nk)
508        call add_array_i1(cv3a_compress_list, kbas1, icb)
509        call add_array_i1(cv3a_compress_list, icbs1, icbs)
510        call add_array_r1(cv3a_compress_list, plcl1, plcl)
511        call add_array_r1(cv3a_compress_list, tnk1, tnk)
512        call add_array_r1(cv3a_compress_list, qnk1, qnk)
513        call add_array_r1(cv3a_compress_list, gznk1, gznk)
514        call add_array_r1(cv3a_compress_list, hnk1, hnk)
515        call add_array_r1(cv3a_compress_list, unk1, unk)
516        call add_array_r1(cv3a_compress_list, vnk1, vnk)
517        call add_array_r2(cv3a_compress_list, wghti1, wghti)
518        call add_array_r1(cv3a_compress_list, pbase1, pbase)
519        call add_array_r1(cv3a_compress_list, buoybase1, buoybase)
520        call add_array_r2(cv3a_compress_list, th1, th)
521        call add_array_r2(cv3a_compress_list, t1, t)
522        call add_array_r2(cv3a_compress_list, q1, q)
523        call add_array_r2(cv3a_compress_list, qs1, qs)
524        call add_array_r2(cv3a_compress_list, t1_wake, t_wake)
525        call add_array_r2(cv3a_compress_list, q1_wake, q_wake)
526        call add_array_r2(cv3a_compress_list, qs1_wake, qs_wake)
527        call add_array_r1(cv3a_compress_list, s1_wake, s_wake)
528        call add_array_r2(cv3a_compress_list, u1, u)
529        call add_array_r2(cv3a_compress_list, v1, v)
530        call add_array_r2(cv3a_compress_list, gz1, gz)
531        call add_array_r2(cv3a_compress_list, h1, h)
532        call add_array_r2(cv3a_compress_list, th1_wake, th_wake)
533        call add_array_r2(cv3a_compress_list, lv1, lv)
534        call add_array_r2(cv3a_compress_list, lf1, lf)
535        call add_array_r2(cv3a_compress_list, cpn1, cpn)
536        call add_array_r2(cv3a_compress_list, p1, p)
537        call add_array_r2(cv3a_compress_list, ph1, ph)
538        call add_array_r2(cv3a_compress_list, tv1, tv)
539        call add_array_r2(cv3a_compress_list, tp1, tp)
540        call add_array_r2(cv3a_compress_list, tvp1, tvp)
541        call add_array_r2(cv3a_compress_list, clw1, clw)
542        call add_array_r2(cv3a_compress_list, h1_wake, h_wake)
543        call add_array_r2(cv3a_compress_list, lv1_wake, lv_wake)
544        call add_array_r2(cv3a_compress_list, lf1_wake, lf_wake)
545        call add_array_r2(cv3a_compress_list, cpn1_wake, cpn_wake)
546        call add_array_r2(cv3a_compress_list, tv1_wake, tv_wake)
547        call add_array_r2(cv3a_compress_list, sig1, sig)
548        call add_array_r1(cv3a_compress_list, sig1(:, nd), sig(:, nd))
549        call add_array_r2(cv3a_compress_list, w01, w0)
550        call add_array_r1(cv3a_compress_list, Ale1, Ale)
551        call add_array_r1(cv3a_compress_list, Alp1, Alp)
552        call add_array_r2(cv3a_compress_list, omega1, omega)
553
554        call cv3a_compress(len, (iflag1 == 0), cv3a_compress_list, compress_data)
555        ncum = compress_data%ncum
556        call exit_profile("cv3a_compress")
557
558        call cv3a_driver_compressed(nloc, nd, ntra, &
559                                    debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
560                                    ! Input fields
561                                    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, &
562                                    ! Input/output
563                                    iflag, icb, plcl, sig, w0, tvp, clw, &
564                                    ! Output fields
565                                    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)
566
567        ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
568        ! --- UNCOMPRESS THE FIELDS
569        ! ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
570        call enter_profile("cv3a_uncompress")
571        call driver_log('cv3a_uncompress')
572        call add_array_i1(cv3a_uncompress_list, iflag, iflag1, init=.false.)
573        call add_array_i1(cv3a_uncompress_list, icb, kbas1)
574        call add_array_i1(cv3a_uncompress_list, inb, ktop1)
575        call add_array_r1(cv3a_uncompress_list, precip, precip1)
576        call add_array_r1(cv3a_uncompress_list, cbmf, cbmf1)
577        call add_array_r1(cv3a_uncompress_list, plcl, plcl1, init=.false.)
578        call add_array_r1(cv3a_uncompress_list, plfc, plfc1)
579        call add_array_r1(cv3a_uncompress_list, wbeff, wbeff1)
580        call add_array_r2(cv3a_uncompress_list, sig, sig1)
581        call add_array_r2(cv3a_uncompress_list, w0, w01)
582        call add_array_r1(cv3a_uncompress_list, ptop2, ptop21)
583        call add_array_r2(cv3a_uncompress_list, ft, ft1)
584        call add_array_r2(cv3a_uncompress_list, fq, fq1)
585        call add_array_r2(cv3a_uncompress_list, fu, fu1)
586        call add_array_r2(cv3a_uncompress_list, fv, fv1)
587        call add_array_r1(cv3a_uncompress_list, sigd, sigd1)
588        call add_array_r2(cv3a_uncompress_list, ma, ma1)
589        call add_array_r2(cv3a_uncompress_list, mip, mip1)
590        call add_array_r2(cv3a_uncompress_list, vprecip, vprecip1)
591        call add_array_r2(cv3a_uncompress_list, vprecipi, vprecipi1)
592        call add_array_r2(cv3a_uncompress_list, upwd, upwd1)
593        call add_array_r2(cv3a_uncompress_list, dnwd, dnwd1)
594        call add_array_r2(cv3a_uncompress_list, dnwd0, dnwd01)
595        call add_array_r2(cv3a_uncompress_list, qcondc, qcondc1)
596        call add_array_r1(cv3a_uncompress_list, wd, wd1)
597        cape1(:) = Cape_noconv
598        call add_array_r1(cv3a_uncompress_list, cape, cape1, init=.false.)
599        cin1(:) = Cin_noconv
600        call add_array_r1(cv3a_uncompress_list, cin, cin1, init=.false.)
601        call add_array_r2(cv3a_uncompress_list, tvp, tvp1, init=.false.)
602        call add_array_r2(cv3a_uncompress_list, ftd, ftd1)
603        call add_array_r2(cv3a_uncompress_list, fqd, fqd1)
604        call add_array_r1(cv3a_uncompress_list, Plim1, Plim11)
605        call add_array_r1(cv3a_uncompress_list, plim2, plim21)
606        call add_array_r2(cv3a_uncompress_list, asupmax, asupmax1)
607        call add_array_r1(cv3a_uncompress_list, supmax0, supmax01)
608        call add_array_r1(cv3a_uncompress_list, asupmaxmin, asupmaxmin1)
609        call add_array_r2(cv3a_uncompress_list, da, da1)
610        call add_array_r3(cv3a_uncompress_list, phi, phi1)
611        call add_array_r2(cv3a_uncompress_list, mp, mp1)
612        call add_array_r3(cv3a_uncompress_list, phi2, phi21)
613        call add_array_r2(cv3a_uncompress_list, d1a, d1a1)
614        call add_array_r2(cv3a_uncompress_list, dam, dam1)
615        call add_array_r3(cv3a_uncompress_list, sigij, sigij1)
616        call add_array_r2(cv3a_uncompress_list, qta, qta1)
617        call add_array_r2(cv3a_uncompress_list, clw, clw1, init=.false.)
618        call add_array_r3(cv3a_uncompress_list, elij, elij1)
619        call add_array_r2(cv3a_uncompress_list, evap, evap1)
620        call add_array_r2(cv3a_uncompress_list, ep, ep1)
621        call add_array_r3(cv3a_uncompress_list, epmlmMm, epmlmMm1)
622        call add_array_r2(cv3a_uncompress_list, eplaMm, eplaMm1)
623        call add_array_r2(cv3a_uncompress_list, wdtrainA, wdtrainA1)
624        call add_array_r2(cv3a_uncompress_list, wdtrainS, wdtrainS1)
625        call add_array_r2(cv3a_uncompress_list, wdtrainM, wdtrainM1)
626        call add_array_r2(cv3a_uncompress_list, qtc, qtc1)
627        call add_array_r2(cv3a_uncompress_list, sigt, sigt1)
628        call add_array_r1(cv3a_uncompress_list, epmax_diag, epmax_diag1)
629        call add_array_r1(cv3a_uncompress_list, sig(:, nd), sig1(:, nd))
630        call cv3a_uncompress(len, compress_data, cv3a_uncompress_list)
631        call exit_profile("cv3a_uncompress")
632
633      end block
634    else !compress
635      nloc = len
636      ktop1 = 0
637      precip1 = 0
638      cbmf1 = 0
639      plfc1 = 0
640      wbeff1 = 0
641      ptop21 = 0
642      sigd1 = 0
643      wd1 = 0
644      cape1 = Cape_noconv
645      cin1 = Cin_noconv
646      Plim11 = 0
647      plim21 = 0
648      supmax01 = 0
649      asupmaxmin1 = 0
650      epmax_diag1 = 0
651      ft1 = 0
652      fq1 = 0
653      fu1 = 0
654      fv1 = 0
655      ma1 = 0
656      mip1 = 0
657      upwd1 = 0
658      dnwd1 = 0
659      dnwd01 = 0
660      qcondc1 = 0
661      ftd1 = 0
662      fqd1 = 0
663      asupmax1 = 0
664      da1 = 0
665      mp1 = 0
666      d1a1 = 0
667      dam1 = 0
668      qta1 = 0
669      evap1 = 0
670      ep1 = 0
671      eplaMm1 = 0
672      wdtrainA1 = 0
673      wdtrainS1 = 0
674      wdtrainM1 = 0
675      qtc1 = 0
676      sigt1 = 0
677      vprecip1 = 0
678      vprecipi1 = 0
679      phi1 = 0
680      phi21 = 0
681      sigij1 = 0
682      elij1 = 0
683      epmlmMm1 = 0
684      coef_convective(:) = merge(1, 0, iflag1(:) == 0)
685      call cv3a_driver_compressed(nloc, nd, ntra, &
686                                  debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
687                                  ! Input fields
688                                  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, &
689                                  ! Input/output
690                                  iflag1, kbas1, plcl1, sig1, w01, tvp1, clw1, &
691                                  ! Output fields
692                                  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)
693      !  Reset Cin1 and Cape1 to their default values for non-convective grid-points.
694      Cin1(:) = Cin1(:)*coef_convective(:) + Cin_noconv*(1.-coef_convective(:))
695      Cape1(:) = Cape1(:)*coef_convective(:) + Cape_noconv*(1.-coef_convective(:))
696      precip1(:) = precip1(:)*coef_convective(:)
697      kbas1(:) = kbas1(:)*coef_convective(:)
698      ktop1(:) = ktop1(:)*coef_convective(:)
699      sigd1(:) = sigd1(:)*coef_convective(:)
700      wbeff1(:) = wbeff1(:)*coef_convective(:)
701      DO k = 1, nd
702        qcondc1(:, k) = qcondc1(:, k)*coef_convective(:)
703        vprecip1(:, k) = vprecip1(:, k)*coef_convective(:)
704        vprecipi1(:, k) = vprecipi1(:, k)*coef_convective(:)
705      ENDDO
706    endif
707
708    IF (prt_level >= 10) THEN
709      Print *, 'cva_driver after cv3_uncompress:ft1(1) , ftd1(1) ', &
710        ft1(igout, 1), ftd1(igout, 1)
711      Print *, 'cva_driver after cv3_uncompress:fq1(1) , fqd1(1) ', &
712        fq1(igout, 1), fqd1(igout, 1)
713    ENDIF
714    IF (debut) THEN
715      PRINT *, ' cv_uncompress -> '
716      debut = .FALSE.
717    END IF
718
719    ftra1(:, :, :) = 0
720
721    call exit_profile("cv3a_driver")
722  END SUBROUTINE
723
724  subroutine cv3a_driver_compressed(nloc, nd, ntra, &
725                                    debut, ok_conserv_q, iflag_mix, iflag_clos, prt_level, &
726                                    ! Input fields
727                                    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, &
728                                    ! Input/output
729                                    iflag, icb, plcl, sig, w0, tvp, clw, &
730                                    ! Output fields
731                                    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)
732    use profiling_physic_mod
733    USE cv3p_mixing_mod, ONLY : cv3p_mixing
734
735    integer, intent(in) :: nloc, nd, ntra ! Arrays sizes
736    logical, intent(in) :: debut, ok_conserv_q
737    integer, intent(in) :: iflag_mix, iflag_clos, prt_level
738    real, intent(in)    :: delt, tau_cld_cv, coefw_cld_cv
739
740    ! Input arrays
741    integer, dimension(nloc), intent(in) :: nk, icbs
742    real, dimension(nloc), intent(in) :: tnk, qnk, gznk, hnk, unk, vnk, pbase, buoybase, s_wake, Ale, Alp
743    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
744    real, dimension(nloc, nd + 1), intent(in) :: ph
745
746    ! Input/output arrays
747    integer, dimension(nloc), intent(inout) :: iflag, icb
748    real, dimension(nloc), intent(inout) :: plcl
749    real, dimension(nloc, nd), intent(inout) :: sig, w0, tvp, clw
750
751    ! Output arrays
752    integer, dimension(nloc), intent(out) :: inb
753    real, dimension(nloc), intent(out) :: precip, cbmf, plfc, wbeff, ptop2, sigd, wd, cape, cin, Plim1, plim2, supmax0, asupmaxmin, epmax_diag
754    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
755    real, dimension(nloc, nd + 1), intent(out) :: vprecip, vprecipi
756    real, dimension(nloc, nd, nd), intent(out) :: phi, phi2, sigij, elij, epmlmMm
757
758    ! Local fields
759    INTEGER nent(nloc, nd)
760    REAL hp(nloc, nd), sigp(nloc, nd), qpreca(nloc, nd)
761    REAL buoy(nloc, nd)
762    REAL m(nloc, nd)
763    REAL mm(nloc, nd)
764    REAL ment(nloc, nd, nd)
765    REAL qent(nloc, nd, nd)
766    REAL hent(nloc, nd, nd)
767    REAL uent(nloc, nd, nd), vent(nloc, nd, nd)
768    REAL ments(nloc, nd, nd), qents(nloc, nd, nd)
769    REAL supmax(nloc, nd)
770    REAL coef_clos(nloc)
771    REAL, DIMENSION(nloc, nd)     :: qp, up, vp
772    REAL, DIMENSION(nloc, nd)     :: wt, water
773    REAL, DIMENSION(nloc, nd)     :: ice, fondue, b
774    REAL, DIMENSION(nloc, nd)     :: frac_a, frac_s, faci
775    REAL tra(nloc, nd, ntra), trap(nloc, nd, ntra)
776    REAL ftra(nloc, nd, ntra), traent(nloc, nd, nd, ntra)
777
778    LOGICAL ok_inhib ! True => possible inhibition of convection by dryness
779    integer :: k
780
781    integer :: ncum ! Number of convective cells
782    logical, save :: timers_first = .true.
783    !$omp threadprivate(timers_first)
784    INTEGER, PARAMETER :: igout = 1
785
786    if (timers_first) then
787      call enter_profile("cv3p_mixing")
788      call exit_profile("cv3p_mixing")
789      call enter_profile("cv3_yield")
790      call exit_profile("cv3_yield")
791      call enter_profile("cv3_tracer")
792      call exit_profile("cv3_tracer")
793
794      timers_first = .false.
795    endif
796
797    ncum = nloc
798    call enter_profile("cv3a_compressed")
799    IF (ncum > 0) THEN
800
801! -------------------------------------------------------------------
802! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
803! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
804! ---   &
805! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
806! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
807! ---   &
808! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
809! -------------------------------------------------------------------
810
811      call driver_log('cv3_undilute2')
812      CALL cv3_undilute2(nloc, ncum, nd, iflag, icb, icbs, nk, &
813                         tnk, qnk, gznk, hnk, t, q, qs, gz, &
814                         p, ph, h, tv, lv, lf, pbase, buoybase, plcl, &
815                         inb, tp, tvp, clw, hp, ep, sigp, buoy, &
816                         frac_a, frac_s, qpreca, qta)
817
818! epmax_cape
819! on recalcule ep et hp
820      call driver_log('cv3_epmax_cape')
821      call cv3_epmax_fn_cape(nloc, ncum, nd &
822                             , ep, hp, icb, inb, clw, nk, t, h, hnk, lv, lf, frac_s &
823                             , pbase, p, ph, tv, buoy, sig, w0, iflag &
824                             , epmax_diag)
825
826! -------------------------------------------------------------------
827! --- MIXING(1)   (if iflag_mix .ge. 1)
828! -------------------------------------------------------------------
829      call enter_profile("cv3p_mixing")
830      IF (iflag_mix >= 1) THEN
831        CALL zilch(supmax, nloc*nd)
832        call driver_log('cv3p_mixing')
833        CALL cv3p_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &
834                         ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qta, &
835                         unk, vnk, hp, tv, tvp, ep, clw, sig, &
836                         ment, qent, hent, uent, vent, nent, &
837                         sigij, elij, supmax, ments, qents, traent)
838      ELSE
839        CALL zilch(supmax, nloc*nd)
840      END IF
841      call exit_profile("cv3p_mixing")
842
843! -------------------------------------------------------------------
844! --- CLOSURE
845! -------------------------------------------------------------------
846
847      ptop2(:) = 0
848      coef_clos(:) = 1.
849      IF (iflag_clos == 0) THEN
850        call driver_log('cv3_closure')
851        CALL cv3_closure(nloc, ncum, nd, icb, inb, &
852                         pbase, p, ph, tv, buoy, &
853                         sig, w0, cape, m, iflag)
854      END IF   ! iflag_clos==0
855
856      ok_inhib = iflag_mix == 2
857
858      IF (iflag_clos == 1) PRINT *, ' pas d appel cv3p_closure'
859
860      IF (iflag_clos == 2) THEN
861        call driver_log('cv3p1_closure')
862        CALL cv3p1_closure(nloc, ncum, nd, icb, inb, &
863                           pbase, plcl, p, ph, tv, tvp, buoy, &
864                           supmax, ok_inhib, Ale, Alp, omega, &
865                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
866                           Plim1, plim2, asupmax, supmax0, &
867                           asupmaxmin, cbmf, plfc, wbeff)
868        if (prt_level >= 10) PRINT *, 'cv3p1_closure-> plfc,wbeff ', plfc(1), wbeff(1)
869      END IF ! iflag_clos==2
870
871      IF (iflag_clos == 3) THEN
872        call driver_log('cv3p2_closure')
873        CALL cv3p2_closure(nloc, ncum, nd, icb, inb, &
874                           pbase, plcl, p, ph, tv, tvp, buoy, &
875                           supmax, ok_inhib, Ale, Alp, omega, &
876                           sig, w0, ptop2, cape, cin, m, iflag, coef_clos, &
877                           Plim1, plim2, asupmax, supmax0, &
878                           asupmaxmin, cbmf, plfc, wbeff)
879        if (prt_level >= 10) PRINT *, 'cv3p2_closure-> plfc,wbeff ', plfc(1), wbeff(1)
880      END IF ! iflag_clos==3
881
882! -------------------------------------------------------------------
883! --- MIXING(2)
884! -------------------------------------------------------------------
885
886      IF (iflag_mix == 0) THEN
887        call driver_log('cv3_mixing')
888        CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb, &             ! na->nd
889                        ph, t, q, qs, u, v, tra, h, lv, lf, frac_s, qnk, &
890                        unk, vnk, hp, tv, tvp, ep, clw, m, sig, &
891                        ment, qent, uent, vent, nent, sigij, elij, ments, qents, traent)
892        CALL zilch(hent, nloc*nd*nd)
893      ELSE
894        mm(:, :) = m(:, :)
895        CALL cv3_mixscale(nloc, ncum, nd, ment, mm)
896        IF (debut) PRINT *, ' cv3_mixscale-> '
897      END IF
898
899      IF (debut) PRINT *, ' cv_mixing ->'
900
901! -------------------------------------------------------------------
902! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
903! -------------------------------------------------------------------
904      IF (debut) PRINT *, ' cva_driver -> cv3_unsat '
905
906      call driver_log('cv3_unsat')
907      CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb, iflag, &
908                     t_wake, q_wake, qs_wake, gz, u, v, tra, p, ph, &
909                     th_wake, tv_wake, lv_wake, lf_wake, cpn_wake, &
910                     ep, sigp, clw, frac_s, qpreca, frac_a, qta, &
911                     m, ment, elij, delt, plcl, coef_clos, &
912                     mp, qp, up, vp, trap, wt, water, evap, fondue, ice, &
913                     faci, b, sigd, &
914                     wdtrainA, wdtrainS, wdtrainM)
915      IF (prt_level >= 10) THEN
916        Print *, 'cva_driver after cv3_unsat:mp , water, ice, evap, fondue '
917        DO k = 1, nd
918          write (6, '(i4,5(1x,e13.6))'), &
919            k, mp(igout, k), water(igout, k), ice(igout, k), &
920            evap(igout, k), fondue(igout, k)
921        ENDDO
922        Print *, 'cva_driver after cv3_unsat: wdtrainA, wdtrainS, wdtrainM '
923        DO k = 1, nd
924          write (6, '(i4,3(1x,e13.6))'), &
925            k, wdtrainA(igout, k), wdtrainS(igout, k), wdtrainM(igout, k)
926        ENDDO
927      ENDIF
928
929      IF (debut) PRINT *, 'cv_unsat-> '
930! -------------------------------------------------------------------
931! YIELD
932! (tendencies, precipitation, variables of interface with other processes, etc)
933! -------------------------------------------------------------------
934
935      call driver_log('cv3_yield')
936      call enter_profile("cv3_yield")
937      CALL cv3_yield(nloc, ncum, nd, nd, ntra, ok_conserv_q, &
938                     icb, inb, delt, &
939                     t, q, t_wake, q_wake, s_wake, u, v, tra, &
940                     gz, p, ph, h, hp, lv, lf, cpn, th, th_wake, &
941                     ep, clw, qpreca, m, tp, mp, qp, up, vp, trap, &
942                     wt, water, ice, evap, fondue, faci, b, sigd, &
943                     ment, qent, hent, iflag_mix, uent, vent, &
944                     nent, elij, traent, sig, &
945                     tv, tvp, wghti, &
946                     iflag, precip, Vprecip, Vprecipi, ft, fq, fu, fv, ftra, &
947                     cbmf, upwd, dnwd, dnwd0, ma, mip, &
948                     qcondc, wd, &
949                     ftd, fqd, qta, qtc, sigt, tau_cld_cv, coefw_cld_cv)
950      call exit_profile("cv3_yield")
951      ! Test conseravtion de l'eau
952      IF (debut) PRINT *, ' cv3_yield -> fqd(1) = ', fqd(igout, 1)
953      IF (prt_level >= 10) THEN
954        Print *, 'cva_driver after cv3_yield:ft(1) , ftd(1) ', &
955          ft(igout, 1), ftd(igout, 1)
956        Print *, 'cva_driver after cv3_yield:fq(1) , fqd(1) ', &
957          fq(igout, 1), fqd(igout, 1)
958      ENDIF
959
960!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
961!--- passive tracers
962!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
963
964      call driver_log('cv3_tracer')
965      call enter_profile("cv3_tracer")
966      CALL cv3_tracer(nloc, -1, ncum, nd, nd, &
967                      ment, sigij, da, phi, phi2, d1a, dam, &
968                      ep, vprecip, elij, clw, epmlmMm, eplaMm, &
969                      icb, inb)
970      call exit_profile("cv3_tracer")
971    END IF ! ncum>0
972    call exit_profile("cv3a_compressed")
973  end subroutine
974
975  subroutine driver_log(message)
976    use print_control_mod, only: prt_level
977    character(*) :: message
978    if (prt_level >= 9) PRINT *, 'cva_driver ->', message
979  end subroutine
980end module
Note: See TracBrowser for help on using the repository browser.