/[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 97 by guez, Fri Apr 25 14:58:31 2014 UTC revision 266 by guez, Thu Apr 19 17:54:55 2018 UTC
# Line 4  module cv_driver_m Line 4  module cv_driver_m
4    
5  contains  contains
6    
7    SUBROUTINE cv_driver(len, nd, 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)         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 comconst, only: dtphys
18      ! - iflag_con = 3: version lmd      use cv30_closure_m, only: cv30_closure
19      ! - iflag_con = 4: version 4.3b      use cv30_compress_m, only: cv30_compress
20        use cv30_feed_m, only: cv30_feed
21      use clesphys2, only: iflag_con      use cv30_mixing_m, only: cv30_mixing
22      use cv3_compress_m, only: cv3_compress      use cv30_param_m, only: cv30_param, nl
23      use cv3_mixing_m, only: cv3_mixing      use cv30_prelim_m, only: cv30_prelim
24      use cv3_param_m, only: cv3_param      use cv30_tracer_m, only: cv30_tracer
25      use cv3_prelim_m, only: cv3_prelim      use cv30_trigger_m, only: cv30_trigger
26      use cv3_tracer_m, only: cv3_tracer      use cv30_uncompress_m, only: cv30_uncompress
27      use cv3_uncompress_m, only: cv3_uncompress      use cv30_undilute1_m, only: cv30_undilute1
28      use cv3_unsat_m, only: cv3_unsat      use cv30_undilute2_m, only: cv30_undilute2
29      use cv3_yield_m, only: cv3_yield      use cv30_unsat_m, only: cv30_unsat
30      use cv_uncompress_m, only: cv_uncompress      use cv30_yield_m, only: cv30_yield
31      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
32    
33      integer, intent(in):: len ! first dimension      real, intent(in):: t1(klon, klev) ! temperature, in K
34      integer, intent(in):: nd ! vertical dimension      real, intent(in):: q1(klon, klev) ! specific humidity
35      real, intent(in):: t1(len, nd) ! temperature      real, intent(in):: qs1(klon, klev) ! saturation specific humidity
36      real q1(len, nd) !           Input        specific hum  
37      real qs1(len, nd)      real, intent(in):: u1(klon, klev), v1(klon, klev)
38      !      qs1           Real           Input        sat specific hum      ! zonal wind and meridional velocity (m/s)
39      real, intent(in):: u1(len, nd)  
40      !      u1            Real           Input        u-wind      real, intent(in):: p1(klon, klev) ! full level pressure, in hPa
41      real, intent(in):: v1(len, nd)  
42      !      v1            Real           Input        v-wind      real, intent(in):: ph1(klon, klev + 1)
43      real p1(len, nd)      ! Half level pressure, in hPa. These pressures are defined at levels
44      !      p1            Real           Input        full level pressure      ! intermediate between those of P1, T1, Q1 and QS1. The first
45      real ph1(len, nd + 1)      ! value of PH should be greater than (i.e. at a lower level than)
46      !      ph1           Real           Input        half level pressure      ! the first value of the array P1.
47      integer iflag1(len)  
48      !      iflag1        Integer        Output       flag for Emanuel conditions      integer, intent(out):: iflag1(:) ! (klon)
49      real ft1(len, nd)      ! Flag for Emanuel conditions.
50      !      ft1           Real           Output       temp tend  
51      real fq1(len, nd)      ! 0: Moist convection occurs.
52      !      fq1           Real           Output       spec hum tend  
53      real fu1(len, nd)      ! 1: Moist convection occurs, but a CFL condition on the
54      !      fu1           Real           Output       u-wind tend      ! subsidence warming is violated. This does not cause the scheme
55      real fv1(len, nd)      ! to terminate.
56      !      fv1           Real           Output       v-wind tend  
57      real precip1(len)      ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
58      !      precip1       Real           Output       precipitation  
59      real VPrecip1(len, nd+1)      ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
60      !      VPrecip1      Real           Output       vertical profile of precipitations  
61      real cbmf1(len)      ! 4: No moist convection; atmosphere is not unstable.
62      !      cbmf1         Real           Output       cloud base mass flux  
63      real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft      ! 6: No moist convection because ihmin <= minorig.
64    
65        ! 7: No moist convection because unreasonable parcel level
66        ! temperature or specific humidity.
67    
68        ! 8: No moist convection: lifted condensation level is above the
69        ! 200 mbar level.
70    
71        ! 9: No moist convection: cloud base is higher than the level NL-1.
72    
73        real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
74        real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
75    
76        real, intent(out):: fu1(klon, klev), fv1(klon, klev)
77        ! forcing (tendency) of zonal and meridional velocity (m/s^2)
78    
79        real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
80    
81        real, intent(out):: VPrecip1(klon, klev + 1)
82        ! vertical profile of convective precipitation (kg/m2/s)
83    
84        real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft
85    
86      real, intent(inout):: w01(klon, klev)      real, intent(inout):: w01(klon, klev)
87      ! vertical velocity within adiabatic updraft      ! vertical velocity within adiabatic updraft
88    
89      integer icb1(klon)      integer, intent(out):: icb1(klon)
90      integer inb1(klon)      integer, intent(inout):: inb1(klon)
91      real, intent(in):: delt      real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
     !      delt          Real           Input        time step  
     real Ma1(len, nd)  
     !      Ma1           Real           Output       mass flux adiabatic updraft  
     real, intent(out):: upwd1(len, nd) ! total upward mass flux (adiab+mixed)  
     real, intent(out):: dnwd1(len, nd) ! saturated downward mass flux (mixed)  
     real, intent(out):: dnwd01(len, nd) ! unsaturated downward mass flux  
   
     real qcondc1(len, nd)     ! cld  
     !      qcondc1       Real           Output       in-cld mixing ratio of condensed water  
     real wd1(len)            ! gust  
     !      wd1           Real           Output       downdraft velocity scale for sfc fluxes  
     real cape1(len)  
     !      cape1         Real           Output       CAPE  
92    
93      real da1(len, nd), phi1(len, nd, nd), mp1(len, nd)      real, intent(out):: upwd1(klon, klev)
94        ! total upward mass flux (adiabatic + mixed)
95    
96      !-------------------------------------------------------------------      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
     ! --- ARGUMENTS  
     !-------------------------------------------------------------------  
     ! --- On input:  
97    
98      !  t:   Array of absolute temperature (K) of dimension ND, with first      real, intent(out):: qcondc1(klon, klev)
99      !       index corresponding to lowest model level. Note that this array      ! in-cloud mixing ratio of condensed water
     !       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 ND, 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 ND, 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 ND, 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 ND, with first  
     !       index corresponding to lowest model level. Must be defined  
     !       at same grid levels as T.  
   
     !  ph:  Array of pressure (mb) of dimension ND+1, with first index  
     !       corresponding to lowest level. These pressures are defined at  
     !       levels intermediate between those of P, T, Q and QS. The first  
     !       value of PH should be greater than (i.e. at a lower level than)  
     !       the first value of the array P.  
   
     !  nl:  The maximum number of levels to which convection can penetrate, plus 1.  
     !       NL MUST be less than or equal to ND-1.  
   
     !  delt: The model time step (sec) between calls to CONVECT  
   
     !----------------------------------------------------------------------------  
     ! ---   On Output:  
   
     !  iflag: An output integer whose value denotes the following:  
     !       VALUE   INTERPRETATION  
     !       -----   --------------  
     !         0     Moist convection occurs.  
     !         1     Moist convection occurs, but a CFL condition  
     !               on the subsidence warming is violated. This  
     !               does not cause the scheme to terminate.  
     !         2     Moist convection, but no precip because ep(inb) lt 0.0001  
     !         3     No moist convection because new cbmf is 0 and old cbmf is 0.  
     !         4     No moist convection; atmosphere is not  
     !               unstable  
     !         6     No moist convection because ihmin le minorig.  
     !         7     No moist convection because unreasonable  
     !               parcel level temperature or specific humidity.  
     !         8     No moist convection: lifted condensation  
     !               level is above the 200 mb level.  
     !         9     No moist convection: cloud base is higher  
     !               then the level NL-1.  
   
     !  ft:   Array of temperature tendency (K/s) of dimension ND, defined at same  
     !        grid levels as T, Q, QS and P.  
   
     !  fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,  
     !        defined at same grid levels as T, Q, QS and P.  
   
     !  fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,  
     !        defined at same grid levels as T.  
   
     !  fv:   Same as FU, but for forcing of meridional velocity.  
   
     !  precip: Scalar convective precipitation rate (mm/day).  
   
     !  VPrecip: Vertical profile of convective precipitation (kg/m2/s).  
   
     !  wd:   A convective downdraft velocity scale. For use in surface  
     !        flux parameterizations. See convect.ps file for details.  
   
     !  tprime: A convective downdraft temperature perturbation scale (K).  
     !          For use in surface flux parameterizations. See convect.ps  
     !          file for details.  
   
     !  qprime: A convective downdraft specific humidity  
     !          perturbation scale (gm/gm).  
     !          For use in surface flux parameterizations. See convect.ps  
     !          file for details.  
   
     !  cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST  
     !        BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT  
     !        ITS NEXT CALL. That is, the value of CBMF must be "remembered"  
     !        by the calling program between calls to CONVECT.  
100    
101      !  det:   Array of detrainment mass flux of dimension ND.      real, intent(out):: cape1(klon)
102        real, intent(out):: da1(:, :) ! (klon, klev)
103        real, intent(out):: phi1(:, :, :) ! (klon, klev, klev)
104    
105      !-------------------------------------------------------------------      real, intent(out):: mp1(:, :) ! (klon, klev) Mass flux of the
106        ! unsaturated downdraft, defined positive downward, in kg m-2
107        ! s-1. M_p in Emanuel (1991 928).
108    
109      !  Local arrays      ! Local:
110    
111      real da(len, nd), phi(len, nd, nd), mp(len, nd)      real da(klon, klev), phi(klon, klev, klev)
112    
113        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 200  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 209  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 nloc      integer, allocatable:: icb(:) ! (ncum)
145      parameter (nloc = klon) ! pour l'instant      integer nent(klon, klev)
146        integer icbs(klon)
147    
148        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)
155        real u(klon, klev), v(klon, klev)
156        real gz(klon, klev), h(klon, klev)
157    
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)
167        real pbase(klon), buoybase(klon)
168        real, allocatable:: th(:, :) ! (ncum, nl)
169        real tvp(klon, klev)
170        real sig(klon, klev), w0(klon, klev)
171        real hp(klon, klev), ep(klon, klev)
172        real buoy(klon, klev)
173        real cape(klon)
174        real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
175        real uent(klon, klev, klev), vent(klon, klev, klev)
176        real ments(klon, klev, klev), qents(klon, klev, klev)
177        real sij(klon, klev, klev), elij(klon, klev, klev)
178        real qp(klon, klev), up(klon, klev), vp(klon, klev)
179        real wt(klon, klev), water(klon, klev)
180        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)
184        real upwd(klon, klev), dnwd(klon, klev)
185        real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
186        real tps(klon, klev)
187        real precip(klon)
188        real VPrecip(klon, klev + 1)
189        real qcondc(klon, klev) ! cloud
190    
     integer idcum(nloc)  
     integer iflag(nloc), nk(nloc), icb(nloc)  
     integer nent(nloc, klev)  
     integer icbs(nloc)  
     integer inb(nloc), inbis(nloc)  
   
     real cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)  
     real t(nloc, klev), q(nloc, klev), qs(nloc, klev)  
     real u(nloc, klev), v(nloc, klev)  
     real gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)  
     real p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)  
     real clw(nloc, klev)  
     real dph(nloc, klev)  
     real pbase(nloc), buoybase(nloc), th(nloc, klev)  
     real tvp(nloc, klev)  
     real sig(nloc, klev), w0(nloc, klev)  
     real hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)  
     real frac(nloc), buoy(nloc, klev)  
     real cape(nloc)  
     real m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)  
     real uent(nloc, klev, klev), vent(nloc, klev, klev)  
     real ments(nloc, klev, klev), qents(nloc, klev, klev)  
     real sij(nloc, klev, klev), elij(nloc, klev, klev)  
     real qp(nloc, klev), up(nloc, klev), vp(nloc, klev)  
     real wt(nloc, klev), water(nloc, klev), evap(nloc, klev)  
     real b(nloc, klev), ft(nloc, klev), fq(nloc, klev)  
     real fu(nloc, klev), fv(nloc, klev)  
     real upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)  
     real Ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)  
     real tps(nloc, klev), qprime(nloc), tprime(nloc)  
     real precip(nloc)  
     real VPrecip(nloc, klev+1)  
     real qcondc(nloc, klev)  ! cld  
     real wd(nloc)           ! gust  
   
     !-------------------------------------------------------------------  
     ! --- SET CONSTANTS AND PARAMETERS  
191      !-------------------------------------------------------------------      !-------------------------------------------------------------------
192    
193      ! -- set simulation flags:      ! SET CONSTANTS AND PARAMETERS
194      !   (common cvflag)      CALL cv30_param
   
     CALL cv_flag  
195    
196      ! -- set thermodynamical constants:      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
     !     (common cvthermo)  
197    
198      CALL cv_thermo      da1 = 0.
199        mp1 = 0.
200      ! -- set convect parameters      phi1 = 0.
201    
202      !     includes microphysical parameters and parameters that      do k = 1, klev
203      !     control the rate of approach to quasi-equilibrium)         do i = 1, klon
204      !     (common cvparam)            ft1(i, k) = 0.
205              fq1(i, k) = 0.
206      if (iflag_con.eq.3) then            fu1(i, k) = 0.
207         CALL cv3_param(nd, delt)            fv1(i, k) = 0.
208      endif            tvp1(i, k) = 0.
209              tp1(i, k) = 0.
210      if (iflag_con.eq.4) then            clw1(i, k) = 0.
211         CALL cv_param(nd)            clw(i, k) = 0.
212      endif            gz1(i, k) = 0.
   
     !---------------------------------------------------------------------  
     ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS  
     !---------------------------------------------------------------------  
   
     do k = 1, nd  
        do  i = 1, len  
           ft1(i, k) = 0.0  
           fq1(i, k) = 0.0  
           fu1(i, k) = 0.0  
           fv1(i, k) = 0.0  
           tvp1(i, k) = 0.0  
           tp1(i, k) = 0.0  
           clw1(i, k) = 0.0  
           !ym  
           clw(i, k) = 0.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, len      precip1 = 0.
222         precip1(i) = 0.0      cape1 = 0.
223         iflag1(i) = 0      VPrecip1(:, klev + 1) = 0.
224         wd1(i) = 0.0  
225         cape1(i) = 0.0      do il = 1, klon
226         VPrecip1(i, nd+1) = 0.0         sig1(il, klev) = sig1(il, klev) + 1.
227      end do         sig1(il, klev) = min(sig1(il, klev), 12.1)
228        enddo
229      if (iflag_con.eq.3) then  
230         do il = 1, len      CALL cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
231            sig1(il, nd) = sig1(il, nd) + 1.      CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, icb1, iflag1, tnk1, qnk1, &
232            sig1(il, nd)  =  min(sig1(il, nd), 12.1)           gznk1, plcl1)
233         enddo      CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, icb1, tp1, tvp1, clw1, &
234      endif           icbs1)
235        CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
236      !--------------------------------------------------------------------           iflag1, sig1, w01)
237      ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY  
238      !--------------------------------------------------------------------      ncum = count(iflag1 == 0)
239    
240      if (iflag_con.eq.3) then      IF (ncum > 0) THEN
241         CALL cv3_prelim(len, nd, nd + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, &         ! Moist convective adjustment is necessary
242              h1, hm1, th1)         allocate(idcum(ncum), plcl(ncum), inb(ncum))
243      endif         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      if (iflag_con.eq.4) then         idcum = pack((/(i, i = 1, klon)/), iflag1 == 0)
246         CALL cv_prelim(len, nd, nd + 1, t1, q1, p1, ph1 &         CALL cv30_compress(idcum, iflag1, icb1, icbs1, plcl1, tnk1, qnk1, &
247              , lv1, cpn1, tv1, gz1, h1, hm1)              gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &
248      endif              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      ! --- CONVECTIVE FEED         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      if (iflag_con.eq.3) then         CALL cv30_closure(icb, inb, pbase, p, ph(:ncum, :), tv, buoy, &
255         CALL cv3_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1            &              sig, w0, cape, m)
256              , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! nd->na         CALL cv30_mixing(icb, inb, t, q, qs, u, v, h, lv, &
257      endif              hp, ep, clw, m, sig, ment, qent, uent, vent, nent, sij, elij, &
258                ments, qents)
259      if (iflag_con.eq.4) then         CALL cv30_unsat(icb, inb, t(:ncum, :nl), q(:ncum, :nl), &
260         CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1 &              qs(:ncum, :nl), gz, u(:ncum, :nl), v(:ncum, :nl), p, &
261              , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)              ph(:ncum, :), th(:ncum, :nl - 1), tv, lv, cpn, ep(:ncum, :), &
262      endif              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      ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part         CALL cv30_yield(icb, inb, dtphys, t, q, u, v, gz, p, ph, h, hp, &
266      ! (up through ICB for convect4, up through ICB+1 for convect3)              lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
267      !     Calculates the lifted parcel virtual temperature at nk, the              wt(:ncum, :nl - 1), water(:ncum, :nl), evap, b, ment, qent, uent, &
268      !     actual temperature, and the adiabatic liquid water content.              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)
271      if (iflag_con.eq.3) then         CALL cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
272         CALL cv3_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1   &              fu, fv, inb, Ma, upwd, dnwd, qcondc, cape, da, phi, mp, iflag1, &
273              , tp1, tvp1, clw1, icbs1) ! nd->na              precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, Ma1, &
274      endif              upwd1, dnwd1, qcondc1, cape1, da1, phi1, mp1)
275        ENDIF
     if (iflag_con.eq.4) then  
        CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax &  
             , tp1, tvp1, clw1)  
     endif  
   
     !-------------------------------------------------------------------  
     ! --- TRIGGERING  
     !-------------------------------------------------------------------  
   
     if (iflag_con.eq.3) then  
        CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &  
             buoybase1, iflag1, sig1, w01) ! nd->na  
     endif  
   
     if (iflag_con.eq.4) then  
        CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)  
     endif  
   
     ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY  
   
     ncum = 0  
     do  i = 1, len  
        if(iflag1(i).eq.0)then  
           ncum = ncum+1  
           idcum(ncum) = i  
        endif  
     end do  
   
     !       print*, 'klon, ncum = ', len, ncum  
   
     IF (ncum.gt.0) THEN  
   
        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^  
        ! --- COMPRESS THE FIELDS  
        !        (-> vectorization over convective gridpoints)  
        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^  
   
        if (iflag_con.eq.3) then  
           CALL cv3_compress(len, nloc, ncum, nd, 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)  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_compress( len, nloc, ncum, nd &  
                , 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.eq.3) then  
           CALL cv3_undilute2(nloc, ncum, nd, icb, icbs, nk         &  
                , tnk, qnk, gznk, t, q, qs, gz &  
                , p, h, tv, lv, pbase, buoybase, plcl &  
                , inb, tp, tvp, clw, hp, ep, sigp, buoy) !na->nd  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_undilute2(nloc, ncum, nd, icb, nk &  
                , tnk, qnk, gznk, t, q, qs, gz &  
                , p, dph, h, tv, lv &  
                , inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)  
        endif  
   
        !-------------------------------------------------------------------  
        ! --- CLOSURE  
        !-------------------------------------------------------------------  
   
        if (iflag_con.eq.3) then  
           CALL cv3_closure(nloc, ncum, nd, icb, inb               &  
                , pbase, p, ph, tv, buoy &  
                , sig, w0, cape, m) ! na->nd  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_closure(nloc, ncum, nd, nk, icb &  
                , tv, tvp, p, ph, dph, plcl, cpn &  
                , iflag, cbmf)  
        endif  
   
        !-------------------------------------------------------------------  
        ! --- MIXING  
        !-------------------------------------------------------------------  
   
        if (iflag_con.eq.3) then  
           CALL cv3_mixing(nloc, ncum, nd, nd, icb, nk, inb, ph, t, q, &  
                qs, u, v, h, lv, qnk, hp, tv, tvp, ep, clw, m, sig, ment, &  
                qent, uent, vent, nent, sij, elij, ments, qents)  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_mixing(nloc, ncum, nd, 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.eq.3) then  
           CALL cv3_unsat(nloc, ncum, nd, nd, 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->nd  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_unsat(nloc, ncum, nd, 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.eq.3) then  
           CALL cv3_yield(nloc, ncum, nd, nd             &  
                , 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->nd  
        endif  
   
        if (iflag_con.eq.4) then  
           CALL cv_yield(nloc, ncum, nd, 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.eq.3) then  
           CALL cv3_tracer(nloc, len, ncum, nd, nd, &  
                ment, sij, da, phi)  
        endif  
   
        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^  
        ! --- UNCOMPRESS THE FIELDS  
        !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^  
        ! set iflag1  = 42 for non convective points  
        do  i = 1, len  
           iflag1(i) = 42  
        end do  
   
        if (iflag_con.eq.3) then  
           CALL cv3_uncompress(nloc, len, ncum, nd, idcum &  
                , 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  
   
        if (iflag_con.eq.4) then  
           CALL cv_uncompress(nloc, len, ncum, nd, idcum &  
                , iflag &  
                , precip, cbmf &  
                , ft, fq, fu, fv &  
                , Ma, qcondc             &  
                , iflag1 &  
                , precip1, cbmf1 &  
                , ft1, fq1, fu1, fv1 &  
                , Ma1, qcondc1 )  
        endif  
     ENDIF ! ncum>0  
276    
277    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver
278    

Legend:
Removed from v.97  
changed lines
  Added in v.266

  ViewVC Help
Powered by ViewVC 1.1.21