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

Legend:
Removed from v.181  
changed lines
  Added in v.195

  ViewVC Help
Powered by ViewVC 1.1.21