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

Last change on this file since 5840 was 5840, checked in by jyg, 2 months ago

Getting rid of tracer arrays within cva_driver.
Lot of comments to be cleared later.

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