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

Diff of /trunk/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 198 by guez, Tue May 31 16:17:35 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  
   
     real, intent(inout):: cbmf1(klon) ! cloud base mass flux  
     real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft  
   
     real, intent(inout):: w01(klon, klev)  
     ! 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  
   
     real, intent(out):: upwd1(klon, klev)  
     ! total upward mass flux (adiab + mixed)  
48    
49      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)      integer, intent(out):: iflag1(:) ! (klon)
50      real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux      ! Flag for Emanuel conditions.
51    
52      real qcondc1(klon, klev) ! cld      ! 0: Moist convection occurs.
     ! 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  
53    
54      real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)      ! 1: Moist convection occurs, but a CFL condition on the
55      real, intent(inout):: mp1(klon, klev)      ! subsidence warming is violated. This does not cause the scheme
56        ! to terminate.
57    
58      ! ARGUMENTS      ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
59    
60      ! On input:      ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
61    
62      ! t: Array of absolute temperature (K) of dimension KLEV, with first      ! 4: No moist convection; atmosphere is not unstable.
     ! 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.  
63    
64      ! nl: The maximum number of levels to which convection can penetrate, plus 1      ! 6: No moist convection because ihmin <= minorig.
     ! NL MUST be less than or equal to KLEV-1.  
65    
66      ! delt: The model time step (sec) between calls to CONVECT      ! 7: No moist convection because unreasonable parcel level
67        ! temperature or specific humidity.
68    
69      ! On Output:      ! 8: No moist convection: lifted condensation level is above the
70        ! 200 mbar level.
71    
72      ! iflag: An output integer whose value denotes the following:      ! 9: No moist convection: cloud base is higher than the level NL-1.
     ! 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.  
73    
74      ! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same      real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
75      ! grid levels as T, Q, QS and P.      real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
76    
77      ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension KLEV,      real, intent(out):: fu1(klon, klev), fv1(klon, klev)
78      ! defined at same grid levels as T, Q, QS and P.      ! forcing (tendency) of zonal and meridional velocity (m/s^2)
79    
80      ! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV,      real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
     ! defined at same grid levels as T.  
81    
82      ! fv: Same as FU, but for forcing of meridional velocity.      real, intent(out):: VPrecip1(klon, klev + 1)
83        ! vertical profile of convective precipitation (kg/m2/s)
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
     integer icbmax  
     integer nk1(klon)  
111      integer icbs1(klon)      integer icbs1(klon)
   
112      real plcl1(klon)      real plcl1(klon)
113      real tnk1(klon)      real tnk1(klon)
114      real qnk1(klon)      real qnk1(klon)
115      real gznk1(klon)      real gznk1(klon)
116      real pbase1(klon)      real pbase1(klon)
117      real buoybase1(klon)      real buoybase1(klon)
   
118      real lv1(klon, klev)      real lv1(klon, klev)
119      real cpn1(klon, klev)      real cpn1(klon, klev)
120      real tv1(klon, klev)      real tv1(klon, klev)
# Line 192  contains Line 125  contains
125      real tvp1(klon, klev)      real tvp1(klon, klev)
126      real clw1(klon, klev)      real clw1(klon, klev)
127      real th1(klon, klev)      real th1(klon, klev)
   
128      integer ncum      integer ncum
129    
130      ! (local) compressed fields:      ! Compressed fields:
131        integer, allocatable:: idcum(:), iflag(:) ! (ncum)
132      integer idcum(klon)      integer, allocatable:: icb(:) ! (ncum)
     integer iflag(klon), nk(klon), icb(klon)  
133      integer nent(klon, klev)      integer nent(klon, klev)
134      integer icbs(klon)      integer icbs(klon)
135      integer inb(klon), inbis(klon)      integer inb(klon)
136        real, allocatable:: plcl(:) ! (ncum)
137      real cbmf(klon), plcl(klon), tnk(klon), qnk(klon), gznk(klon)      real 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) ! pressure at full level, in hPa
142        real ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
143      real clw(klon, klev)      real clw(klon, klev)
     real dph(klon, klev)  
144      real pbase(klon), buoybase(klon), th(klon, klev)      real pbase(klon), buoybase(klon), th(klon, klev)
145      real tvp(klon, klev)      real tvp(klon, klev)
146      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
147      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)      real hp(klon, klev), ep(klon, klev)
148      real frac(klon), buoy(klon, klev)      real buoy(klon, klev)
149      real cape(klon)      real cape(klon)
150      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
151      real uent(klon, klev, klev), vent(klon, klev, klev)      real uent(klon, klev, klev), vent(klon, klev, klev)
152      real ments(klon, klev, klev), qents(klon, klev, klev)      real ments(klon, klev, klev), qents(klon, klev, klev)
153      real sij(klon, klev, klev), elij(klon, klev, klev)      real sij(klon, klev, klev), elij(klon, klev, klev)
154      real qp(klon, klev), up(klon, klev), vp(klon, klev)      real qp(klon, klev), up(klon, klev), vp(klon, klev)
155      real wt(klon, klev), water(klon, klev), evap(klon, klev)      real wt(klon, klev), water(klon, klev)
156      real b(klon, klev), ft(klon, klev), fq(klon, klev)      real, allocatable:: evap(:, :) ! (ncum, nl)
157        real, allocatable:: b(:, :) ! (ncum, nl - 1)
158        real ft(klon, klev), fq(klon, klev)
159      real fu(klon, klev), fv(klon, klev)      real fu(klon, klev), fv(klon, klev)
160      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)      real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
161      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)      real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
162      real tps(klon, klev), qprime(klon), tprime(klon)      real tps(klon, klev)
163      real precip(klon)      real precip(klon)
164      real VPrecip(klon, klev + 1)      real VPrecip(klon, klev + 1)
165      real qcondc(klon, klev) ! cld      real qcondc(klon, klev) ! cld
     real wd(klon) ! gust  
166    
167      !-------------------------------------------------------------------      !-------------------------------------------------------------------
168    
169      ! SET CONSTANTS AND PARAMETERS      ! SET CONSTANTS AND PARAMETERS
   
     ! set simulation flags:  
     ! (common cvflag)  
   
     CALL cv_flag  
   
     ! set thermodynamical constants:  
     ! (common cvthermo)  
   
170      CALL cv_thermo      CALL cv_thermo
171        CALL cv30_param
     ! set convect parameters  
   
     ! includes microphysical parameters and parameters that  
     ! control the rate of approach to quasi-equilibrium)  
     ! (common cvparam)  
   
     CALL cv3_param(klev, delt)  
172    
173      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
174    
175      do k = 1, klev      do k = 1, klev
176         do i = 1, klon         do i = 1, klon
177            ft1(i, k) = 0.0            ft1(i, k) = 0.
178            fq1(i, k) = 0.0            fq1(i, k) = 0.
179            fu1(i, k) = 0.0            fu1(i, k) = 0.
180            fv1(i, k) = 0.0            fv1(i, k) = 0.
181            tvp1(i, k) = 0.0            tvp1(i, k) = 0.
182            tp1(i, k) = 0.0            tp1(i, k) = 0.
183            clw1(i, k) = 0.0            clw1(i, k) = 0.
184            !ym            clw(i, k) = 0.
           clw(i, k) = 0.0  
185            gz1(i, k) = 0.            gz1(i, k) = 0.
186            VPrecip1(i, k) = 0.            VPrecip1(i, k) = 0.
187            Ma1(i, k) = 0.0            Ma1(i, k) = 0.
188            upwd1(i, k) = 0.0            upwd1(i, k) = 0.
189            dnwd1(i, k) = 0.0            dnwd1(i, k) = 0.
190            dnwd01(i, k) = 0.0            dnwd01(i, k) = 0.
191            qcondc1(i, k) = 0.0            qcondc1(i, k) = 0.
192         end do         end do
193      end do      end do
194    
195      do i = 1, klon      precip1 = 0.
196         precip1(i) = 0.0      cape1 = 0.
197         iflag1(i) = 0      VPrecip1(:, klev + 1) = 0.
        wd1(i) = 0.0  
        cape1(i) = 0.0  
        VPrecip1(i, klev + 1) = 0.0  
     end do  
198    
199      do il = 1, klon      do il = 1, klon
200         sig1(il, klev) = sig1(il, klev) + 1.         sig1(il, klev) = sig1(il, klev) + 1.
201         sig1(il, klev) = min(sig1(il, klev), 12.1)         sig1(il, klev) = min(sig1(il, klev), 12.1)
202      enddo      enddo
203    
204      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY      CALL cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
205        CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, icb1, iflag1, tnk1, qnk1, &
206      CALL cv3_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &           gznk1, plcl1)
207           gz1, h1, hm1, th1)      CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, icb1, tp1, tvp1, clw1, &
208             icbs1)
209      ! CONVECTIVE FEED      CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
210             iflag1, sig1, w01)
211    
212      CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &      ncum = count(iflag1 == 0)
          icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na  
   
     ! 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  
   
     ! Moist convective adjustment is necessary  
   
     ncum = 0  
     do i = 1, klon  
        if (iflag1(i) == 0) then  
           ncum = ncum + 1  
           idcum(ncum) = i  
        endif  
     end do  
213    
214      IF (ncum > 0) THEN      IF (ncum > 0) THEN
215         ! COMPRESS THE FIELDS         ! Moist convective adjustment is necessary
216         ! (-> vectorization over convective gridpoints)         allocate(idcum(ncum), plcl(ncum))
217           allocate(b(ncum, nl - 1), evap(ncum, nl), icb(ncum), iflag(ncum))
218         CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &         idcum = pack((/(i, i = 1, klon)/), iflag1 == 0)
219              plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &         CALL cv30_compress(iflag1, icb1, icbs1, plcl1, tnk1, qnk1, gznk1, &
220              v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &              pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
221              sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &              p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, icb, icbs, plcl, tnk, &
222              buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &              qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, &
223              tvp, clw, sig, w0)              p, ph, tv, tp, tvp, clw, sig, w0)
224           CALL cv30_undilute2(icb, icbs(:ncum), tnk, qnk, gznk, t, qs, gz, p, h, &
225         ! UNDILUTE (ADIABATIC) UPDRAFT / second part :              tv, lv, pbase(:ncum), buoybase(:ncum), plcl, inb(:ncum), tp, tvp, &
226         ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES              clw, hp, ep, buoy)
227         ! &         CALL cv30_closure(icb, inb(:ncum), pbase, p, ph(:ncum, :), tv, buoy, &
228         ! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE              sig, w0, cape, m)
229         ! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD         CALL cv30_mixing(icb, inb(:ncum), t, q, qs, u, v, h, lv, hp, ep, clw, &
230         ! &              m, sig, ment, qent, uent, vent, nent, sij, elij, ments, qents)
231         ! FIND THE LEVEL OF NEUTRAL BUOYANCY         CALL cv30_unsat(icb, inb(:ncum), t(:ncum, :nl), q(:ncum, :nl), &
232                qs(:ncum, :nl), gz, u(:ncum, :nl), v(:ncum, :nl), p, &
233         CALL cv3_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, &              ph(:ncum, :), th(:ncum, :nl - 1), tv, lv(:ncum, :), &
234              t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, &              cpn(:ncum, :nl), ep(:ncum, :), clw(:ncum, :), m(:ncum, :), &
235              tvp, clw, hp, ep, sigp, buoy) !na->klev              ment(:ncum, :, :), elij(:ncum, :, :), dtphys, plcl, mp, &
236                qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), wt(:ncum, :nl), &
237         ! CLOSURE              water(:ncum, :nl), evap, b)
238           CALL cv30_yield(icb, inb(:ncum), dtphys, t, q, u, v, gz, p, ph, h, hp, &
239         CALL cv3_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &              lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
240              buoy, sig, w0, cape, m) ! na->klev              wt(:ncum, :nl - 1), water(:ncum, :nl), evap, b, ment, qent, uent, &
241                vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, &
242         ! MIXING              fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
243           CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
244         CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &         CALL cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
245              v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &              fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, da, phi, mp, &
246              sij, elij, ments, qents)              iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, &
247                Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, phi1, mp1)
248         ! UNSATURATED (PRECIPITATING) DOWNDRAFTS      ENDIF
   
        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)  
   
        ! UNCOMPRESS THE FIELDS  
   
        ! set iflag1 = 42 for non convective points  
        do i = 1, klon  
           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, &  
             da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &  
             fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &  
             cape1, da1, phi1, mp1)  
     ENDIF ! ncum>0  
249    
250    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver
251    

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

  ViewVC Help
Powered by ViewVC 1.1.21