/[lmdze]/trunk/Sources/phylmd/cv_driver.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 180 by guez, Tue Mar 15 17:07:47 2016 UTC revision 187 by guez, Mon Mar 21 18:01:02 2016 UTC
# Line 4  module cv_driver_m Line 4  module cv_driver_m
4    
5  contains  contains
6    
7    SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, &    SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8         fq1, fu1, fv1, precip1, VPrecip1, cbmf1, sig1, w01, icb1, inb1, delt, &         fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, &
9         Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)         dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
10    
11      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12      ! Main driver for convection      ! Main driver for convection
# Line 14  contains Line 14  contains
14    
15      ! Several modules corresponding to different physical processes      ! Several modules corresponding to different physical processes
16    
17      ! Several versions of convect may be used:      use cv30_closure_m, only: cv30_closure
18      ! - iflag_con = 3: version lmd      use cv30_compress_m, only: cv30_compress
19      ! - iflag_con = 4: version 4.3b      use cv30_feed_m, only: cv30_feed
20        use cv30_mixing_m, only: cv30_mixing
21      use clesphys2, only: iflag_con      use cv30_param_m, only: cv30_param
22      use cv3_compress_m, only: cv3_compress      use cv30_prelim_m, only: cv30_prelim
23      use cv3_feed_m, only: cv3_feed      use cv30_tracer_m, only: cv30_tracer
24      use cv3_mixing_m, only: cv3_mixing      use cv30_uncompress_m, only: cv30_uncompress
25      use cv3_param_m, only: cv3_param      use cv30_undilute2_m, only: cv30_undilute2
26      use cv3_prelim_m, only: cv3_prelim      use cv30_unsat_m, only: cv30_unsat
27      use cv3_tracer_m, only: cv3_tracer      use cv30_yield_m, only: cv30_yield
     use cv3_uncompress_m, only: cv3_uncompress  
     use cv3_unsat_m, only: cv3_unsat  
     use cv3_yield_m, only: cv3_yield  
     use cv_feed_m, only: cv_feed  
     use cv_uncompress_m, only: cv_uncompress  
28      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
29    
30      real, intent(in):: t1(klon, klev) ! temperature      real, intent(in):: t1(klon, klev)
31      real, intent(in):: q1(klon, klev) ! specific hum      ! temperature (K), with first index corresponding to lowest model
32      real, intent(in):: qs1(klon, klev) ! sat specific hum      ! level
33      real, intent(in):: u1(klon, klev) ! u-wind  
34      real, intent(in):: v1(klon, klev) ! v-wind      real, intent(in):: q1(klon, klev)
35      real, intent(in):: p1(klon, klev) ! full level pressure      ! Specific humidity, with first index corresponding to lowest
36      real, intent(in):: ph1(klon, klev + 1) ! half level pressure      ! model level. Must be defined at same grid levels as T1.
37      integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions  
38      real, intent(out):: ft1(klon, klev) ! temp tend      real, intent(in):: qs1(klon, klev)
39      real, intent(out):: fq1(klon, klev) ! spec hum tend      ! Saturation specific humidity, with first index corresponding to
40      real, intent(out):: fu1(klon, klev) ! u-wind tend      ! lowest model level. Must be defined at same grid levels as
41      real, intent(out):: fv1(klon, klev) ! v-wind tend      ! T1.
42      real, intent(out):: precip1(klon) ! precipitation  
43        real, intent(in):: u1(klon, klev), v1(klon, klev)
44        ! Zonal wind and meridional velocity (m/s), witth first index
45        ! corresponding with the lowest model level. Defined at same
46        ! levels as T1.
47    
48        real, intent(in):: p1(klon, klev)
49        ! Full level pressure (mb) of dimension KLEV, with first index
50        ! corresponding to lowest model level. Must be defined at same
51        ! grid levels as T1.
52    
53        real, intent(in):: ph1(klon, klev + 1)
54        ! Half level pressure (mb), with first index corresponding to
55        ! lowest level. These pressures are defined at levels intermediate
56        ! between those of P1, T1, Q1 and QS1. The first value of PH
57        ! should be greater than (i.e. at a lower level than) the first
58        ! value of the array P1.
59    
60        integer, intent(out):: iflag1(klon)
61        ! Flag for Emanuel conditions.
62    
63        ! 0: Moist convection occurs.
64    
65        ! 1: Moist convection occurs, but a CFL condition on the
66        ! subsidence warming is violated. This does not cause the scheme
67        ! to terminate.
68    
69        ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
70    
71        ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
72    
73        ! 4: No moist convection; atmosphere is not unstable
74    
75        ! 6: No moist convection because ihmin le minorig.
76    
77        ! 7: No moist convection because unreasonable parcel level
78        ! temperature or specific humidity.
79    
80        ! 8: No moist convection: lifted condensation level is above the
81        ! 200 mb level.
82    
83        ! 9: No moist convection: cloud base is higher then the level NL-1.
84    
85        real, intent(out):: ft1(klon, klev)
86        ! Temperature tendency (K/s), defined at same grid levels as T1,
87        ! Q1, QS1 and P1.
88    
89        real, intent(out):: fq1(klon, klev)
90        ! Specific humidity tendencies (s-1), defined at same grid levels
91        ! as T1, Q1, QS1 and P1.
92    
93        real, intent(out):: fu1(klon, klev), fv1(klon, klev)
94        ! Forcing (tendency) of zonal and meridional velocity (m/s^2),
95        ! defined at same grid levels as T1.
96    
97        real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
98    
99      real, intent(out):: VPrecip1(klon, klev + 1)      real, intent(out):: VPrecip1(klon, klev + 1)
100      ! vertical profile of precipitation      ! vertical profile of convective precipitation (kg/m2/s)
101    
     real, intent(inout):: cbmf1(klon) ! cloud base mass flux  
102      real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft      real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
103    
104      real, intent(inout):: w01(klon, klev)      real, intent(inout):: w01(klon, klev)
# Line 57  contains Line 106  contains
106    
107      integer, intent(out):: icb1(klon)      integer, intent(out):: icb1(klon)
108      integer, intent(inout):: inb1(klon)      integer, intent(inout):: inb1(klon)
109      real, intent(in):: delt ! time step      real, intent(in):: delt ! the model time step (sec) between calls
110      real Ma1(klon, klev)  
111      ! Ma1 Real Output mass flux adiabatic updraft      real Ma1(klon, klev) ! Output mass flux adiabatic updraft
112    
113      real, intent(out):: upwd1(klon, klev)      real, intent(out):: upwd1(klon, klev)
114      ! total upward mass flux (adiab + mixed)      ! total upward mass flux (adiab + mixed)
# Line 67  contains Line 116  contains
116      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
117      real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux      real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
118    
119      real qcondc1(klon, klev) ! cld      real qcondc1(klon, klev) ! Output in-cld mixing ratio of condensed water
120      ! qcondc1 Real Output in-cld mixing ratio of condensed water  
121      real wd1(klon) ! gust      real wd1(klon) ! gust
122      ! wd1 Real Output downdraft velocity scale for sfc fluxes      ! Output downdraft velocity scale for surface fluxes
123      real cape1(klon)      ! A convective downdraft velocity scale. For use in surface
124      ! cape1 Real Output CAPE      ! flux parameterizations. See convect.ps file for details.
125    
126        real cape1(klon) ! Output
127      real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)      real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
128      real, intent(inout):: mp1(klon, klev)      real, intent(inout):: mp1(klon, klev)
129    
130      ! ARGUMENTS      ! Local:
   
     ! On input:  
   
     ! t: Array of absolute temperature (K) of dimension KLEV, with first  
     ! index corresponding to lowest model level. Note that this array  
     ! will be altered by the subroutine if dry convective adjustment  
     ! occurs and if IPBL is not equal to 0.  
   
     ! q: Array of specific humidity (gm/gm) of dimension KLEV, with first  
     ! index corresponding to lowest model level. Must be defined  
     ! at same grid levels as T. Note that this array will be altered  
     ! if dry convective adjustment occurs and if IPBL is not equal to 0.  
   
     ! qs: Array of saturation specific humidity of dimension KLEV, with first  
     ! index corresponding to lowest model level. Must be defined  
     ! at same grid levels as T. Note that this array will be altered  
     ! if dry convective adjustment occurs and if IPBL is not equal to 0.  
   
     ! u: Array of zonal wind velocity (m/s) of dimension KLEV, witth first  
     ! index corresponding with the lowest model level. Defined at  
     ! same levels as T. Note that this array will be altered if  
     ! dry convective adjustment occurs and if IPBL is not equal to 0.  
   
     ! v: Same as u but for meridional velocity.  
   
     ! p: Array of pressure (mb) of dimension KLEV, with first  
     ! index corresponding to lowest model level. Must be defined  
     ! at same grid levels as T.  
   
     ! ph: Array of pressure (mb) of dimension KLEV + 1, with first index  
     ! corresponding to lowest level. These pressures are defined at  
     ! levels intermediate between those of P, T, Q and QS. The first  
     ! value of PH should be greater than (i.e. at a lower level than)  
     ! the first value of the array P.  
   
     ! nl: The maximum number of levels to which convection can penetrate, plus 1  
     ! NL MUST be less than or equal to KLEV-1.  
   
     ! delt: The model time step (sec) between calls to CONVECT  
   
     ! On Output:  
   
     ! iflag: An output integer whose value denotes the following:  
     ! VALUE INTERPRETATION  
     ! ----- --------------  
     ! 0 Moist convection occurs.  
     ! 1 Moist convection occurs, but a CFL condition  
     ! on the subsidence warming is violated. This  
     ! does not cause the scheme to terminate.  
     ! 2 Moist convection, but no precip because ep(inb) lt 0.0001  
     ! 3 No moist convection because new cbmf is 0 and old cbmf is 0.  
     ! 4 No moist convection; atmosphere is not  
     ! unstable  
     ! 6 No moist convection because ihmin le minorig.  
     ! 7 No moist convection because unreasonable  
     ! parcel level temperature or specific humidity.  
     ! 8 No moist convection: lifted condensation  
     ! level is above the 200 mb level.  
     ! 9 No moist convection: cloud base is higher  
     ! then the level NL-1.  
   
     ! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same  
     ! grid levels as T, Q, QS and P.  
   
     ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension KLEV,  
     ! defined at same grid levels as T, Q, QS and P.  
   
     ! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV,  
     ! defined at same grid levels as T.  
   
     ! fv: Same as FU, but for forcing of meridional velocity.  
   
     ! precip: Scalar convective precipitation rate (mm/day).  
   
     ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).  
   
     ! wd: A convective downdraft velocity scale. For use in surface  
     ! flux parameterizations. See convect.ps file for details.  
   
     ! tprime: A convective downdraft temperature perturbation scale (K).  
     ! For use in surface flux parameterizations. See convect.ps  
     ! file for details.  
   
     ! qprime: A convective downdraft specific humidity  
     ! perturbation scale (gm/gm).  
     ! For use in surface flux parameterizations. See convect.ps  
     ! file for details.  
   
     ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST  
     ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT  
     ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"  
     ! by the calling program between calls to CONVECT.  
   
     ! det: Array of detrainment mass flux of dimension KLEV.  
   
     ! Local arrays  
131    
132      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
133    
# Line 202  contains Line 156  contains
156    
157      integer ncum      integer ncum
158    
159      ! (local) compressed fields:      ! Compressed fields:
160    
161      integer idcum(klon)      integer idcum(klon)
162      integer iflag(klon), nk(klon), icb(klon)      integer iflag(klon), nk(klon), icb(klon)
163      integer nent(klon, klev)      integer nent(klon, klev)
164      integer icbs(klon)      integer icbs(klon)
165      integer inb(klon), inbis(klon)      integer inb(klon)
166    
167      real cbmf(klon), plcl(klon), tnk(klon), qnk(klon), gznk(klon)      real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
168      real t(klon, klev), q(klon, klev), qs(klon, klev)      real t(klon, klev), q(klon, klev), qs(klon, klev)
169      real u(klon, klev), v(klon, klev)      real u(klon, klev), v(klon, klev)
170      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
171      real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)      real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
172      real clw(klon, klev)      real clw(klon, klev)
     real dph(klon, klev)  
173      real pbase(klon), buoybase(klon), th(klon, klev)      real pbase(klon), buoybase(klon), th(klon, klev)
174      real tvp(klon, klev)      real tvp(klon, klev)
175      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
176      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)
177      real frac(klon), buoy(klon, klev)      real buoy(klon, klev)
178      real cape(klon)      real cape(klon)
179      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
180      real uent(klon, klev, klev), vent(klon, klev, klev)      real uent(klon, klev, klev), vent(klon, klev, klev)
# Line 233  contains Line 186  contains
186      real fu(klon, klev), fv(klon, klev)      real fu(klon, klev), fv(klon, klev)
187      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
188      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
189      real tps(klon, klev), qprime(klon), tprime(klon)      real tps(klon, klev)
190      real precip(klon)      real precip(klon)
191      real VPrecip(klon, klev + 1)      real VPrecip(klon, klev + 1)
192      real qcondc(klon, klev) ! cld      real qcondc(klon, klev) ! cld
# Line 243  contains Line 196  contains
196    
197      ! SET CONSTANTS AND PARAMETERS      ! SET CONSTANTS AND PARAMETERS
198    
     ! set simulation flags:  
     ! (common cvflag)  
   
     CALL cv_flag  
   
199      ! set thermodynamical constants:      ! set thermodynamical constants:
200      ! (common cvthermo)      ! (common cvthermo)
   
201      CALL cv_thermo      CALL cv_thermo
202    
203      ! set convect parameters      ! set convect parameters
   
204      ! includes microphysical parameters and parameters that      ! includes microphysical parameters and parameters that
205      ! control the rate of approach to quasi-equilibrium)      ! control the rate of approach to quasi-equilibrium)
206      ! (common cvparam)      ! (common cvparam)
207        CALL cv30_param(delt)
     if (iflag_con == 3) CALL cv3_param(klev, delt)  
208    
209      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
210    
211      do k = 1, klev      do k = 1, klev
212         do i = 1, klon         do i = 1, klon
213            ft1(i, k) = 0.0            ft1(i, k) = 0.
214            fq1(i, k) = 0.0            fq1(i, k) = 0.
215            fu1(i, k) = 0.0            fu1(i, k) = 0.
216            fv1(i, k) = 0.0            fv1(i, k) = 0.
217            tvp1(i, k) = 0.0            tvp1(i, k) = 0.
218            tp1(i, k) = 0.0            tp1(i, k) = 0.
219            clw1(i, k) = 0.0            clw1(i, k) = 0.
220            !ym            clw(i, k) = 0.
           clw(i, k) = 0.0  
221            gz1(i, k) = 0.            gz1(i, k) = 0.
222            VPrecip1(i, k) = 0.            VPrecip1(i, k) = 0.
223            Ma1(i, k) = 0.0            Ma1(i, k) = 0.
224            upwd1(i, k) = 0.0            upwd1(i, k) = 0.
225            dnwd1(i, k) = 0.0            dnwd1(i, k) = 0.
226            dnwd01(i, k) = 0.0            dnwd01(i, k) = 0.
227            qcondc1(i, k) = 0.0            qcondc1(i, k) = 0.
228         end do         end do
229      end do      end do
230    
231      do i = 1, klon      do i = 1, klon
232         precip1(i) = 0.0         precip1(i) = 0.
233         iflag1(i) = 0         iflag1(i) = 0
234         wd1(i) = 0.0         wd1(i) = 0.
235         cape1(i) = 0.0         cape1(i) = 0.
236         VPrecip1(i, klev + 1) = 0.0         VPrecip1(i, klev + 1) = 0.
237      end do      end do
238    
239      if (iflag_con == 3) then      do il = 1, klon
240         do il = 1, klon         sig1(il, klev) = sig1(il, klev) + 1.
241            sig1(il, klev) = sig1(il, klev) + 1.         sig1(il, klev) = min(sig1(il, klev), 12.1)
242            sig1(il, klev) = min(sig1(il, klev), 12.1)      enddo
        enddo  
     endif  
243    
244      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
245        CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
246      if (iflag_con == 3) then           gz1, h1, hm1, th1)
        CALL cv3_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &  
             gz1, h1, hm1, th1)  
     else  
        ! iflag_con == 4  
        CALL cv_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &  
             gz1, h1, hm1)  
     endif  
247    
248      ! CONVECTIVE FEED      ! CONVECTIVE FEED
249        CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
250      if (iflag_con == 3) then           icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
        CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &  
             icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na  
     else  
        ! iflag_con == 4  
        CALL cv_feed(klon, klev, t1, q1, qs1, p1, hm1, gz1, nk1, icb1, icbmax, &  
             iflag1, tnk1, qnk1, gznk1, plcl1)  
     endif  
251    
252      ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part      ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
253      ! (up through ICB for convect4, up through ICB + 1 for convect3)      ! (up through ICB for convect4, up through ICB + 1 for convect3)
254      ! Calculates the lifted parcel virtual temperature at nk, the      ! Calculates the lifted parcel virtual temperature at nk, the
255      ! actual temperature, and the adiabatic liquid water content.      ! actual temperature, and the adiabatic liquid water content.
256        CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
257      if (iflag_con == 3) then           tp1, tvp1, clw1, icbs1) ! klev->na
        CALL cv3_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &  
             tp1, tvp1, clw1, icbs1) ! klev->na  
     else  
        ! iflag_con == 4  
        CALL cv_undilute1(klon, klev, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, &  
             tp1, tvp1, clw1)  
     endif  
258    
259      ! TRIGGERING      ! TRIGGERING
260        CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
261      if (iflag_con == 3) then           buoybase1, iflag1, sig1, w01) ! klev->na
        CALL cv3_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &  
             buoybase1, iflag1, sig1, w01) ! klev->na  
     else  
        ! iflag_con == 4  
        CALL cv_trigger(klon, klev, icb1, cbmf1, tv1, tvp1, iflag1)  
     end if  
262    
263      ! Moist convective adjustment is necessary      ! Moist convective adjustment is necessary
264    
# Line 358  contains Line 273  contains
273      IF (ncum > 0) THEN      IF (ncum > 0) THEN
274         ! COMPRESS THE FIELDS         ! COMPRESS THE FIELDS
275         ! (-> vectorization over convective gridpoints)         ! (-> vectorization over convective gridpoints)
276           CALL cv30_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
277         if (iflag_con == 3) then              plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
278            CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &              v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
279                 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &              sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
280                 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &              buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
281                 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &              tvp, clw, sig, w0)
282                 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &  
283                 tvp, clw, sig, w0)         CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, &
284         else              h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, &
285            ! iflag_con == 4              ep, sigp, buoy)
           CALL cv_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, cbmf1, &  
                plcl1, tnk1, qnk1, gznk1, t1, q1, qs1, u1, v1, gz1, h1, lv1, &  
                cpn1, p1, ph1, tv1, tp1, tvp1, clw1, iflag, nk, icb, cbmf, &  
                plcl, tnk, qnk, gznk, t, q, qs, u, v, gz, h, lv, cpn, p, ph, &  
                tv, tp, tvp, clw, dph)  
        endif  
   
        ! UNDILUTE (ADIABATIC) UPDRAFT / second part :  
        ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES  
        ! &  
        ! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE  
        ! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD  
        ! &  
        ! FIND THE LEVEL OF NEUTRAL BUOYANCY  
   
        if (iflag_con == 3) then  
           CALL cv3_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, &  
                t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, &  
                tvp, clw, hp, ep, sigp, buoy) !na->klev  
        else  
           ! iflag_con == 4  
           CALL cv_undilute2(klon, ncum, klev, icb, nk, tnk, qnk, gznk, t, &  
                qs, gz, p, dph, h, tv, lv, inb, inbis, tp, tvp, clw, hp, ep, &  
                sigp, frac)  
        endif  
286    
287         ! CLOSURE         ! CLOSURE
288           CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
289         if (iflag_con == 3) then              buoy, sig, w0, cape, m) ! na->klev
           CALL cv3_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &  
                buoy, sig, w0, cape, m) ! na->klev  
        else  
           ! iflag_con == 4  
           CALL cv_closure(klon, ncum, klev, nk, icb, tv, tvp, p, ph, dph, &  
                plcl, cpn, iflag, cbmf)  
        endif  
290    
291         ! MIXING         ! MIXING
292           CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
293         if (iflag_con == 3) then              v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
294            CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &              sij, elij, ments, qents)
295                 v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &  
296                 sij, elij, ments, qents)         ! Unsaturated (precipitating) downdrafts
297         else         CALL cv30_unsat(ncum, icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, &
298            ! iflag_con == 4              ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, plcl, &
299            CALL cv_mixing(klon, ncum, klev, icb, nk, inb, inbis, ph, t, q, qs, &              mp, qp, up, vp, wt, water, evap, b(:ncum, :))
300                 u, v, h, lv, qnk, hp, tv, tvp, ep, clw, cbmf, m, ment, qent, &  
301                 uent, vent, nent, sij, elij)         ! Yield (tendencies, precipitation, variables of interface with
302         endif         ! other processes, etc)
303           CALL cv30_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &
304         ! UNSATURATED (PRECIPITATING) DOWNDRAFTS              gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &
305                wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &
306         if (iflag_con == 3) then              tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &
307            CALL cv3_unsat(klon, ncum, klev, klev, icb, inb, t, q, qs, gz, u, &              dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev
                v, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, &  
                plcl, mp, qp, up, vp, wt, water, evap, b)! na->klev  
        else  
           ! iflag_con == 4  
           CALL cv_unsat(klon, ncum, klev, inb, t, q, qs, gz, u, v, p, ph, h, &  
                lv, ep, sigp, clw, m, ment, elij, iflag, mp, qp, up, vp, wt, &  
                water, evap)  
        endif  
   
        ! YIELD  
        ! (tendencies, precipitation, variables of interface with other  
        ! processes, etc)  
   
        if (iflag_con == 3) then  
           CALL cv3_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &  
                gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &  
                wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &  
                tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &  
                dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev  
        else  
           ! iflag_con == 4  
           CALL cv_yield(klon, ncum, klev, nk, icb, inb, delt, t, q, u, v, gz, &  
                p, ph, h, hp, lv, cpn, ep, clw, frac, m, mp, qp, up, vp, wt, &  
                water, evap, ment, qent, uent, vent, nent, elij, tv, tvp, &  
                iflag, wd, qprime, tprime, precip, cbmf, ft, fq, fu, fv, Ma, &  
                qcondc)  
        endif  
308    
309         ! passive tracers         ! passive tracers
310           CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
        if (iflag_con == 3) CALL cv3_tracer(klon, ncum, klev, ment, sij, da, phi)  
311    
312         ! UNCOMPRESS THE FIELDS         ! UNCOMPRESS THE FIELDS
313    
314         ! set iflag1 = 42 for non convective points         ! set iflag1 = 42 for non convective points
315         do i = 1, klon         iflag1 = 42
           iflag1(i) = 42  
        end do  
316    
317         if (iflag_con == 3) then         CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
318            CALL cv3_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &              ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
319                 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &              da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
320                 da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &              fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
321                 fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &              cape1, da1, phi1, mp1)
322                 cape1, da1, phi1, mp1)      ENDIF
        else  
           ! iflag_con == 4  
           CALL cv_uncompress(idcum(:ncum), iflag, precip, cbmf, ft, fq, fu, &  
                fv, Ma, qcondc, iflag1, precip1, cbmf1, ft1, fq1, fu1, fv1, &  
                Ma1, qcondc1)  
        endif  
     ENDIF ! ncum>0  
323    
324    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver
325    

Legend:
Removed from v.180  
changed lines
  Added in v.187

  ViewVC Help
Powered by ViewVC 1.1.21