/[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 189 by guez, Tue Mar 29 15:20:23 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, dnwd1, &
9         Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)         dnwd01, qcondc1, 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, nl
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_trigger_m, only: cv30_trigger
25      use cv3_param_m, only: cv3_param      use cv30_uncompress_m, only: cv30_uncompress
26      use cv3_prelim_m, only: cv3_prelim      use cv30_undilute2_m, only: cv30_undilute2
27      use cv3_tracer_m, only: cv3_tracer      use cv30_unsat_m, only: cv30_unsat
28      use cv3_uncompress_m, only: cv3_uncompress      use cv30_yield_m, only: cv30_yield
     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  
29      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
30    
31      real, intent(in):: t1(klon, klev) ! temperature      real, intent(in):: t1(klon, klev) ! temperature (K)
32      real, intent(in):: q1(klon, klev) ! specific hum      real, intent(in):: q1(klon, klev) ! specific humidity
33      real, intent(in):: qs1(klon, klev) ! sat specific hum      real, intent(in):: qs1(klon, klev) ! saturation specific humidity
34      real, intent(in):: u1(klon, klev) ! u-wind  
35      real, intent(in):: v1(klon, klev) ! v-wind      real, intent(in):: u1(klon, klev), v1(klon, klev)
36      real, intent(in):: p1(klon, klev) ! full level pressure      ! zonal wind and meridional velocity (m/s)
37      real, intent(in):: ph1(klon, klev + 1) ! half level pressure  
38      integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions      real, intent(in):: p1(klon, klev) ! full level pressure (hPa)
39      real, intent(out):: ft1(klon, klev) ! temp tend  
40      real, intent(out):: fq1(klon, klev) ! spec hum tend      real, intent(in):: ph1(klon, klev + 1)
41      real, intent(out):: fu1(klon, klev) ! u-wind tend      ! Half level pressure (hPa). These pressures are defined at levels
42      real, intent(out):: fv1(klon, klev) ! v-wind tend      ! intermediate between those of P1, T1, Q1 and QS1. The first
43      real, intent(out):: precip1(klon) ! precipitation      ! value of PH should be greater than (i.e. at a lower level than)
44        ! the first value of the array P1.
45    
46      real, intent(out):: VPrecip1(klon, klev + 1)      integer, intent(out):: iflag1(klon)
47      ! vertical profile of precipitation      ! Flag for Emanuel conditions.
48    
49      real, intent(inout):: cbmf1(klon) ! cloud base mass flux      ! 0: Moist convection occurs.
     real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft  
50    
51      real, intent(inout):: w01(klon, klev)      ! 1: Moist convection occurs, but a CFL condition on the
52      ! vertical velocity within adiabatic updraft      ! subsidence warming is violated. This does not cause the scheme
53        ! to terminate.
54    
55      integer, intent(out):: icb1(klon)      ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
     integer, intent(inout):: inb1(klon)  
     real, intent(in):: delt ! time step  
     real Ma1(klon, klev)  
     ! Ma1 Real Output mass flux adiabatic updraft  
56    
57      real, intent(out):: upwd1(klon, klev)      ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
     ! total upward mass flux (adiab + mixed)  
58    
59      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)      ! 4: No moist convection; atmosphere is not unstable
     real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux  
   
     real qcondc1(klon, klev) ! cld  
     ! qcondc1 Real Output in-cld mixing ratio of condensed water  
     real wd1(klon) ! gust  
     ! wd1 Real Output downdraft velocity scale for sfc fluxes  
     real cape1(klon)  
     ! cape1 Real Output CAPE  
60    
61      real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)      ! 6: No moist convection because ihmin le minorig.
     real, intent(inout):: mp1(klon, klev)  
62    
63      ! ARGUMENTS      ! 7: No moist convection because unreasonable parcel level
64        ! temperature or specific humidity.
     ! 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.  
65    
66      ! nl: The maximum number of levels to which convection can penetrate, plus 1      ! 8: No moist convection: lifted condensation level is above the
67      ! NL MUST be less than or equal to KLEV-1.      ! 200 mb level.
68    
69      ! delt: The model time step (sec) between calls to CONVECT      ! 9: No moist convection: cloud base is higher then the level NL-1.
70    
71      ! On Output:      real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
72        real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
73    
74      ! iflag: An output integer whose value denotes the following:      real, intent(out):: fu1(klon, klev), fv1(klon, klev)
75      ! VALUE INTERPRETATION      ! forcing (tendency) of zonal and meridional velocity (m/s^2)
     ! ----- --------------  
     ! 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.  
76    
77      ! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same      real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
     ! grid levels as T, Q, QS and P.  
78    
79      ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension KLEV,      real, intent(out):: VPrecip1(klon, klev + 1)
80      ! defined at same grid levels as T, Q, QS and P.      ! vertical profile of convective precipitation (kg/m2/s)
   
     ! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV,  
     ! defined at same grid levels as T.  
81    
82      ! fv: Same as FU, but for forcing of meridional velocity.      real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft
83    
84      ! precip: Scalar convective precipitation rate (mm/day).      real, intent(inout):: w01(klon, klev)
85        ! vertical velocity within adiabatic updraft
86    
87      ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).      integer, intent(out):: icb1(klon)
88        integer, intent(inout):: inb1(klon)
89        real, intent(in):: delt ! the model time step (sec) between calls
90    
91      ! wd: A convective downdraft velocity scale. For use in surface      real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
     ! flux parameterizations. See convect.ps file for details.  
92    
93      ! tprime: A convective downdraft temperature perturbation scale (K).      real, intent(out):: upwd1(klon, klev)
94      ! For use in surface flux parameterizations. See convect.ps      ! total upward mass flux (adiabatic + mixed)
     ! file for details.  
95    
96      ! qprime: A convective downdraft specific humidity      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
97      ! perturbation scale (gm/gm).      real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
     ! For use in surface flux parameterizations. See convect.ps  
     ! file for details.  
98    
99      ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST      real, intent(out):: qcondc1(klon, klev)
100      ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT      ! in-cloud mixing ratio of condensed water
     ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"  
     ! by the calling program between calls to CONVECT.  
101    
102      ! det: Array of detrainment mass flux of dimension KLEV.      real, intent(out):: cape1(klon)
103        real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
104        real, intent(inout):: mp1(klon, klev)
105    
106      ! Local arrays      ! Local:
107    
108      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
   
109      integer i, k, il      integer i, k, il
110      integer icbmax      integer icbmax
111      integer nk1(klon)      integer nk1(klon)
112      integer icbs1(klon)      integer icbs1(klon)
   
113      real plcl1(klon)      real plcl1(klon)
114      real tnk1(klon)      real tnk1(klon)
115      real qnk1(klon)      real qnk1(klon)
116      real gznk1(klon)      real gznk1(klon)
117      real pbase1(klon)      real pbase1(klon)
118      real buoybase1(klon)      real buoybase1(klon)
   
119      real lv1(klon, klev)      real lv1(klon, klev)
120      real cpn1(klon, klev)      real cpn1(klon, klev)
121      real tv1(klon, klev)      real tv1(klon, klev)
# Line 199  contains Line 126  contains
126      real tvp1(klon, klev)      real tvp1(klon, klev)
127      real clw1(klon, klev)      real clw1(klon, klev)
128      real th1(klon, klev)      real th1(klon, klev)
   
129      integer ncum      integer ncum
130    
131      ! (local) compressed fields:      ! Compressed fields:
   
132      integer idcum(klon)      integer idcum(klon)
133      integer iflag(klon), nk(klon), icb(klon)      integer iflag(klon), nk(klon), icb(klon)
134      integer nent(klon, klev)      integer nent(klon, klev)
135      integer icbs(klon)      integer icbs(klon)
136      integer inb(klon), inbis(klon)      integer inb(klon)
137        real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
     real cbmf(klon), plcl(klon), tnk(klon), qnk(klon), gznk(klon)  
138      real t(klon, klev), q(klon, klev), qs(klon, klev)      real t(klon, klev), q(klon, klev), qs(klon, klev)
139      real u(klon, klev), v(klon, klev)      real u(klon, klev), v(klon, klev)
140      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)
141      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)
142      real clw(klon, klev)      real clw(klon, klev)
     real dph(klon, klev)  
143      real pbase(klon), buoybase(klon), th(klon, klev)      real pbase(klon), buoybase(klon), th(klon, klev)
144      real tvp(klon, klev)      real tvp(klon, klev)
145      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
146      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)
147      real frac(klon), buoy(klon, klev)      real buoy(klon, klev)
148      real cape(klon)      real cape(klon)
149      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
150      real uent(klon, klev, klev), vent(klon, klev, klev)      real uent(klon, klev, klev), vent(klon, klev, klev)
# Line 229  contains Line 152  contains
152      real sij(klon, klev, klev), elij(klon, klev, klev)      real sij(klon, klev, klev), elij(klon, klev, klev)
153      real qp(klon, klev), up(klon, klev), vp(klon, klev)      real qp(klon, klev), up(klon, klev), vp(klon, klev)
154      real wt(klon, klev), water(klon, klev), evap(klon, klev)      real wt(klon, klev), water(klon, klev), evap(klon, klev)
155      real b(klon, klev), ft(klon, klev), fq(klon, klev)      real, allocatable:: b(:, :) ! (ncum, nl - 1)
156        real ft(klon, klev), fq(klon, klev)
157      real fu(klon, klev), fv(klon, klev)      real fu(klon, klev), fv(klon, klev)
158      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
159      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
160      real tps(klon, klev), qprime(klon), tprime(klon)      real tps(klon, klev)
161      real precip(klon)      real precip(klon)
162      real VPrecip(klon, klev + 1)      real VPrecip(klon, klev + 1)
163      real qcondc(klon, klev) ! cld      real qcondc(klon, klev) ! cld
     real wd(klon) ! gust  
164    
165      !-------------------------------------------------------------------      !-------------------------------------------------------------------
166    
167      ! SET CONSTANTS AND PARAMETERS      ! SET CONSTANTS AND PARAMETERS
168    
     ! set simulation flags:  
     ! (common cvflag)  
   
     CALL cv_flag  
   
169      ! set thermodynamical constants:      ! set thermodynamical constants:
170      ! (common cvthermo)      ! (common cvthermo)
   
171      CALL cv_thermo      CALL cv_thermo
172    
173      ! set convect parameters      ! set convect parameters
   
174      ! includes microphysical parameters and parameters that      ! includes microphysical parameters and parameters that
175      ! control the rate of approach to quasi-equilibrium)      ! control the rate of approach to quasi-equilibrium)
176      ! (common cvparam)      ! (common cvparam)
177        CALL cv30_param(delt)
     if (iflag_con == 3) CALL cv3_param(klev, delt)  
178    
179      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
180    
181      do k = 1, klev      do k = 1, klev
182         do i = 1, klon         do i = 1, klon
183            ft1(i, k) = 0.0            ft1(i, k) = 0.
184            fq1(i, k) = 0.0            fq1(i, k) = 0.
185            fu1(i, k) = 0.0            fu1(i, k) = 0.
186            fv1(i, k) = 0.0            fv1(i, k) = 0.
187            tvp1(i, k) = 0.0            tvp1(i, k) = 0.
188            tp1(i, k) = 0.0            tp1(i, k) = 0.
189            clw1(i, k) = 0.0            clw1(i, k) = 0.
190            !ym            clw(i, k) = 0.
           clw(i, k) = 0.0  
191            gz1(i, k) = 0.            gz1(i, k) = 0.
192            VPrecip1(i, k) = 0.            VPrecip1(i, k) = 0.
193            Ma1(i, k) = 0.0            Ma1(i, k) = 0.
194            upwd1(i, k) = 0.0            upwd1(i, k) = 0.
195            dnwd1(i, k) = 0.0            dnwd1(i, k) = 0.
196            dnwd01(i, k) = 0.0            dnwd01(i, k) = 0.
197            qcondc1(i, k) = 0.0            qcondc1(i, k) = 0.
198         end do         end do
199      end do      end do
200    
201      do i = 1, klon      do i = 1, klon
202         precip1(i) = 0.0         precip1(i) = 0.
203         iflag1(i) = 0         iflag1(i) = 0
204         wd1(i) = 0.0         cape1(i) = 0.
205         cape1(i) = 0.0         VPrecip1(i, klev + 1) = 0.
        VPrecip1(i, klev + 1) = 0.0  
206      end do      end do
207    
208      if (iflag_con == 3) then      do il = 1, klon
209         do il = 1, klon         sig1(il, klev) = sig1(il, klev) + 1.
210            sig1(il, klev) = sig1(il, klev) + 1.         sig1(il, klev) = min(sig1(il, klev), 12.1)
211            sig1(il, klev) = min(sig1(il, klev), 12.1)      enddo
        enddo  
     endif  
212    
213      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
214        CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
215      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  
216    
217      ! CONVECTIVE FEED      ! CONVECTIVE FEED
218        CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
219             icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
220    
221      if (iflag_con == 3) then      CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
222         CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &           tp1, tvp1, clw1, icbs1) ! klev->na
             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  
   
     ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part  
     ! (up through ICB for convect4, up through ICB + 1 for convect3)  
     ! Calculates the lifted parcel virtual temperature at nk, the  
     ! actual temperature, and the adiabatic liquid water content.  
   
     if (iflag_con == 3) then  
        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  
223    
224      ! TRIGGERING      ! TRIGGERING
225        CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
226      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  
227    
228      ! Moist convective adjustment is necessary      ! Moist convective adjustment is necessary
229    
# Line 356  contains Line 236  contains
236      end do      end do
237    
238      IF (ncum > 0) THEN      IF (ncum > 0) THEN
239         ! COMPRESS THE FIELDS         allocate(b(ncum, nl - 1))
240         ! (-> vectorization over convective gridpoints)         CALL cv30_compress(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
241                gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &
242         if (iflag_con == 3) then              cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, &
243            CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &              icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, &
244                 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &              th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)
245                 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &         CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, &
246                 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &              h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, &
247                 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &              ep, sigp, buoy)
                tvp, clw, sig, w0)  
        else  
           ! iflag_con == 4  
           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  
248    
249         ! CLOSURE         ! CLOSURE
250           CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
251         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  
252    
253         ! MIXING         ! MIXING
254           CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
255                v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
256                sij, elij, ments, qents)
257    
258           ! Unsaturated (precipitating) downdrafts
259           CALL cv30_unsat(icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, ph, th, &
260                tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, plcl, mp, &
261                qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), wt, water, evap, b)
262    
263           ! Yield (tendencies, precipitation, variables of interface with
264           ! other processes, etc)
265           CALL cv30_yield(icb(:ncum), inb(:ncum), delt, t, q, u, v, gz, p, ph, &
266                h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, wt, &
267                water(:ncum, :nl), evap(:ncum, :nl), b, ment, qent, uent, vent, &
268                nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, &
269                upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
270    
271         if (iflag_con == 3) then         CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
           CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &  
                v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &  
                sij, elij, ments, qents)  
        else  
           ! iflag_con == 4  
           CALL cv_mixing(klon, ncum, klev, icb, nk, inb, inbis, ph, t, q, qs, &  
                u, v, h, lv, qnk, hp, tv, tvp, ep, clw, cbmf, m, ment, qent, &  
                uent, vent, nent, sij, elij)  
        endif  
   
        ! UNSATURATED (PRECIPITATING) DOWNDRAFTS  
   
        if (iflag_con == 3) then  
           CALL cv3_unsat(klon, ncum, klev, klev, icb, inb, t, q, qs, gz, u, &  
                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  
   
        ! passive tracers  
   
        if (iflag_con == 3) CALL cv3_tracer(klon, ncum, klev, ment, sij, da, phi)  
272    
273         ! UNCOMPRESS THE FIELDS         ! UNCOMPRESS THE FIELDS
274           iflag1 = 42 ! for non convective points
275         ! set iflag1 = 42 for non convective points         CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
276         do i = 1, klon              ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, &
277            iflag1(i) = 42              da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
278         end do              fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, &
279                phi1, mp1)
280         if (iflag_con == 3) then      ENDIF
           CALL cv3_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &  
                ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &  
                da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &  
                fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &  
                cape1, da1, phi1, mp1)  
        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  
281    
282    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver
283    

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

  ViewVC Help
Powered by ViewVC 1.1.21