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

Legend:
Removed from v.185  
changed lines
  Added in v.205

  ViewVC Help
Powered by ViewVC 1.1.21