/[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

trunk/libf/phylmd/cv_driver.f90 revision 52 by guez, Fri Sep 23 12:28:01 2011 UTC trunk/Sources/phylmd/cv_driver.f revision 183 by guez, Wed Mar 16 14:42:58 2016 UTC
# Line 4  module cv_driver_m Line 4  module cv_driver_m
4    
5  contains  contains
6    
7    SUBROUTINE cv_driver(len, nd, ndp1, ntra, iflag_con, t1, q1, qs1, u1, v1, &    SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8         tra1, p1, ph1, iflag1, ft1, fq1, fu1, fv1, ftra1, precip1, VPrecip1, &         fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, &
9         cbmf1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, dnwd1, dnwd01, &         dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
10         qcondc1, wd1, cape1, da1, phi1, mp1)  
11        ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3 2005/04/15 12:36:17      ! Main driver for convection
13        ! Author: S. Bony, March 2002
14      use dimens_m  
15      use dimphy      ! Several modules corresponding to different physical processes
16      !  
17      ! PARAMETERS:      use cv3_compress_m, only: cv3_compress
18      !      Name            Type         Usage            Description      use cv3_feed_m, only: cv3_feed
19      !   ----------      ----------     -------  ----------------------------      use cv3_mixing_m, only: cv3_mixing
20      !      use cv3_param_m, only: cv3_param
21      !      len           Integer        Input        first (i) dimension      use cv3_prelim_m, only: cv3_prelim
22      !      nd            Integer        Input        vertical (k) dimension      use cv3_tracer_m, only: cv3_tracer
23      !      ndp1          Integer        Input        nd + 1      use cv3_uncompress_m, only: cv3_uncompress
24      !      ntra          Integer        Input        number of tracors      use cv3_undilute2_m, only: cv3_undilute2
25      !      iflag_con     Integer        Input        version of convect (3/4)      use cv3_unsat_m, only: cv3_unsat
26      !      t1            Real           Input        temperature      use cv3_yield_m, only: cv3_yield
27      !      q1            Real           Input        specific hum      USE dimphy, ONLY: klev, klon
28      !      qs1           Real           Input        sat specific hum  
29      !      u1            Real           Input        u-wind      real, intent(in):: t1(klon, klev) ! temperature
30      !      v1            Real           Input        v-wind      real, intent(in):: q1(klon, klev) ! specific hum
31      !      tra1          Real           Input        tracors      real, intent(in):: qs1(klon, klev) ! sat specific hum
32      !      p1            Real           Input        full level pressure      real, intent(in):: u1(klon, klev) ! u-wind
33      !      ph1           Real           Input        half level pressure      real, intent(in):: v1(klon, klev) ! v-wind
34      !      iflag1        Integer        Output       flag for Emanuel conditions      real, intent(in):: p1(klon, klev) ! full level pressure
35      !      ft1           Real           Output       temp tend      real, intent(in):: ph1(klon, klev + 1) ! half level pressure
36      !      fq1           Real           Output       spec hum tend      integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions
37      !      fu1           Real           Output       u-wind tend      real, intent(out):: ft1(klon, klev) ! temp tend
38      !      fv1           Real           Output       v-wind tend      real, intent(out):: fq1(klon, klev) ! spec hum tend
39      !      ftra1         Real           Output       tracor tend      real, intent(out):: fu1(klon, klev) ! u-wind tend
40      !      precip1       Real           Output       precipitation      real, intent(out):: fv1(klon, klev) ! v-wind tend
41      !      VPrecip1      Real           Output       vertical profile of precipitations      real, intent(out):: precip1(klon) ! precipitation
42      !      cbmf1         Real           Output       cloud base mass flux  
43      !      sig1          Real           In/Out       section adiabatic updraft      real, intent(out):: VPrecip1(klon, klev + 1)
44      !      w01           Real           In/Out       vertical velocity within adiab updraft      ! vertical profile of precipitation
45      !      delt          Real           Input        time step  
46      !      Ma1           Real           Output       mass flux adiabatic updraft      real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
47      !      upwd1         Real           Output       total upward mass flux (adiab+mixed)  
48      !      dnwd1         Real           Output       saturated downward mass flux (mixed)      real, intent(inout):: w01(klon, klev)
49      !      dnwd01        Real           Output       unsaturated downward mass flux      ! vertical velocity within adiabatic updraft
50      !      qcondc1       Real           Output       in-cld mixing ratio of condensed water  
51      !      wd1           Real           Output       downdraft velocity scale for sfc fluxes      integer, intent(out):: icb1(klon)
52      !      cape1         Real           Output       CAPE      integer, intent(inout):: inb1(klon)
53      !      real, intent(in):: delt ! time step
54      ! S. Bony, Mar 2002:      real Ma1(klon, klev)
55      !     * Several modules corresponding to different physical processes      ! Ma1 Real Output mass flux adiabatic updraft
56      !     * Several versions of convect may be used:  
57      !        - iflag_con=3: version lmd  (previously named convect3)      real, intent(out):: upwd1(klon, klev)
58      !        - iflag_con=4: version 4.3b (vect. version, previously convect1/2)      ! total upward mass flux (adiab + mixed)
59      !   + tard:    - iflag_con=5: version lmd with ice (previously named convectg)  
60      ! S. Bony, Oct 2002:      real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
61      !     * Vectorization of convect3 (ie version lmd)      real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
62      !  
63      !..............................END PROLOGUE.............................      real qcondc1(klon, klev) ! cld
64      !      ! qcondc1 Real Output in-cld mixing ratio of condensed water
65      !      real wd1(klon) ! gust
66        ! wd1 Real Output downdraft velocity scale for sfc fluxes
67      integer len      real cape1(klon)
68      integer nd      ! cape1 Real Output CAPE
69      integer ndp1  
70      integer noff      real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
71      integer, intent(in):: iflag_con      real, intent(inout):: mp1(klon, klev)
72      integer ntra  
73      real, intent(in):: t1(len, nd)      ! ARGUMENTS
74      real q1(len, nd)  
75      real qs1(len, nd)      ! On input:
76      real u1(len, nd)  
77      real v1(len, nd)      ! t: Array of absolute temperature (K) of dimension KLEV, with first
78      real p1(len, nd)      ! index corresponding to lowest model level. Note that this array
79      real ph1(len, ndp1)      ! will be altered by the subroutine if dry convective adjustment
80      integer iflag1(len)      ! occurs and if IPBL is not equal to 0.
81      real ft1(len, nd)  
82      real fq1(len, nd)      ! q: Array of specific humidity (gm/gm) of dimension KLEV, with first
83      real fu1(len, nd)      ! index corresponding to lowest model level. Must be defined
84      real fv1(len, nd)      ! at same grid levels as T. Note that this array will be altered
85      real precip1(len)      ! if dry convective adjustment occurs and if IPBL is not equal to 0.
86      real cbmf1(len)  
87      real VPrecip1(len, nd+1)      ! qs: Array of saturation specific humidity of dimension KLEV, with first
88      real Ma1(len, nd)      ! index corresponding to lowest model level. Must be defined
89      real upwd1(len, nd)      ! at same grid levels as T. Note that this array will be altered
90      real dnwd1(len, nd)      ! if dry convective adjustment occurs and if IPBL is not equal to 0.
91      real dnwd01(len, nd)  
92        ! u: Array of zonal wind velocity (m/s) of dimension KLEV, witth first
93      real qcondc1(len, nd)     ! cld      ! index corresponding with the lowest model level. Defined at
94      real wd1(len)            ! gust      ! same levels as T. Note that this array will be altered if
95      real cape1(len)      ! dry convective adjustment occurs and if IPBL is not equal to 0.
96    
97      real da1(len, nd), phi1(len, nd, nd), mp1(len, nd)      ! v: Same as u but for meridional velocity.
98      real da(len, nd), phi(len, nd, nd), mp(len, nd)  
99      real, intent(in):: tra1(len, nd, ntra)      ! p: Array of pressure (mb) of dimension KLEV, with first
100      real ftra1(len, nd, ntra)      ! index corresponding to lowest model level. Must be defined
101        ! at same grid levels as T.
102    
103        ! ph: Array of pressure (mb) of dimension KLEV + 1, with first index
104        ! corresponding to lowest level. These pressures are defined at
105        ! levels intermediate between those of P, T, Q and QS. The first
106        ! value of PH should be greater than (i.e. at a lower level than)
107        ! the first value of the array P.
108    
109        ! nl: The maximum number of levels to which convection can penetrate, plus 1
110        ! NL MUST be less than or equal to KLEV-1.
111    
112        ! delt: The model time step (sec) between calls to CONVECT
113    
114        ! On Output:
115    
116        ! iflag: An output integer whose value denotes the following:
117        ! VALUE INTERPRETATION
118        ! ----- --------------
119        ! 0 Moist convection occurs.
120        ! 1 Moist convection occurs, but a CFL condition
121        ! on the subsidence warming is violated. This
122        ! does not cause the scheme to terminate.
123        ! 2 Moist convection, but no precip because ep(inb) lt 0.0001
124        ! 3 No moist convection because new cbmf is 0 and old cbmf is 0.
125        ! 4 No moist convection; atmosphere is not
126        ! unstable
127        ! 6 No moist convection because ihmin le minorig.
128        ! 7 No moist convection because unreasonable
129        ! parcel level temperature or specific humidity.
130        ! 8 No moist convection: lifted condensation
131        ! level is above the 200 mb level.
132        ! 9 No moist convection: cloud base is higher
133        ! then the level NL-1.
134    
135        ! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same
136        ! grid levels as T, Q, QS and P.
137    
138        ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension KLEV,
139        ! defined at same grid levels as T, Q, QS and P.
140    
141        ! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV,
142        ! defined at same grid levels as T.
143    
144        ! fv: Same as FU, but for forcing of meridional velocity.
145    
146        ! precip: Scalar convective precipitation rate (mm/day).
147    
148        ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
149    
150        ! wd: A convective downdraft velocity scale. For use in surface
151        ! flux parameterizations. See convect.ps file for details.
152    
153        ! tprime: A convective downdraft temperature perturbation scale (K).
154        ! For use in surface flux parameterizations. See convect.ps
155        ! file for details.
156    
157        ! qprime: A convective downdraft specific humidity
158        ! perturbation scale (gm/gm).
159        ! For use in surface flux parameterizations. See convect.ps
160        ! file for details.
161    
162        ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
163        ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
164        ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"
165        ! by the calling program between calls to CONVECT.
166    
167      real, intent(in):: delt      ! det: Array of detrainment mass flux of dimension KLEV.
168    
169      !-------------------------------------------------------------------      ! Local arrays
     ! --- ARGUMENTS  
     !-------------------------------------------------------------------  
     ! --- On input:  
     !  
     !  t:   Array of absolute temperature (K) of dimension ND, 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 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.  
     !  
     !  tra: Array of passive tracer mixing ratio, of dimensions (ND, NTRA),  
     !       where NTRA is the number of different tracers. If no  
     !       convective tracer transport is needed, define a dummy  
     !       input array of dimension (ND, 1). Tracers are defined at  
     !       same vertical levels as T. Note that this array will be altered  
     !       if dry convective adjustment occurs and if IPBL is not equal to 0.  
     !  
     !  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.  
     !  
     !  ftra: Array of forcing of tracer content, in tracer mixing ratio per  
     !        second, defined at same levels as T. Dimensioned (ND, NTRA).  
     !  
     !  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.  
     !  
     !  det:   Array of detrainment mass flux of dimension ND.  
     !  
     !-------------------------------------------------------------------  
     !  
     !  Local arrays  
     !  
170    
171      integer i, k, n, il, j      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
172    
173        integer i, k, il
174      integer icbmax      integer icbmax
175      integer nk1(klon)      integer nk1(klon)
     integer icb1(klon)  
     integer inb1(klon)  
176      integer icbs1(klon)      integer icbs1(klon)
177    
178      real plcl1(klon)      real plcl1(klon)
179      real tnk1(klon)      real tnk1(klon)
180      real qnk1(klon)      real qnk1(klon)
181      real gznk1(klon)      real gznk1(klon)
     real pnk1(klon)  
     real qsnk1(klon)  
182      real pbase1(klon)      real pbase1(klon)
183      real buoybase1(klon)      real buoybase1(klon)
184    
# Line 239  contains Line 191  contains
191      real tp1(klon, klev)      real tp1(klon, klev)
192      real tvp1(klon, klev)      real tvp1(klon, klev)
193      real clw1(klon, klev)      real clw1(klon, klev)
     real sig1(klon, klev)  
     real w01(klon, klev)  
194      real th1(klon, klev)      real th1(klon, klev)
195      !  
196      integer ncum      integer ncum
197      !  
198      ! (local) compressed fields:      ! (local) compressed fields:
     !  
     integer nloc  
     parameter (nloc=klon) ! pour l'instant  
   
     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 tra(nloc, klev, ntra), trap(nloc, klev, ntra)  
     real ftra(nloc, klev, ntra), traent(nloc, klev, klev, ntra)  
     real qcondc(nloc, klev)  ! cld  
     real wd(nloc)           ! gust  
199    
200      !-------------------------------------------------------------------      integer idcum(klon)
201      ! --- SET CONSTANTS AND PARAMETERS      integer iflag(klon), nk(klon), icb(klon)
202        integer nent(klon, klev)
203        integer icbs(klon)
204        integer inb(klon)
205    
206        real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
207        real t(klon, klev), q(klon, klev), qs(klon, klev)
208        real u(klon, klev), v(klon, klev)
209        real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
210        real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
211        real clw(klon, klev)
212        real pbase(klon), buoybase(klon), th(klon, klev)
213        real tvp(klon, klev)
214        real sig(klon, klev), w0(klon, klev)
215        real hp(klon, klev), ep(klon, klev), sigp(klon, klev)
216        real buoy(klon, klev)
217        real cape(klon)
218        real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
219        real uent(klon, klev, klev), vent(klon, klev, klev)
220        real ments(klon, klev, klev), qents(klon, klev, klev)
221        real sij(klon, klev, klev), elij(klon, klev, klev)
222        real qp(klon, klev), up(klon, klev), vp(klon, klev)
223        real wt(klon, klev), water(klon, klev), evap(klon, klev)
224        real b(klon, klev), ft(klon, klev), fq(klon, klev)
225        real fu(klon, klev), fv(klon, klev)
226        real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
227        real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
228        real tps(klon, klev)
229        real precip(klon)
230        real VPrecip(klon, klev + 1)
231        real qcondc(klon, klev) ! cld
232        real wd(klon) ! gust
233    
234      !-------------------------------------------------------------------      !-------------------------------------------------------------------
235    
236      ! -- set simulation flags:      ! SET CONSTANTS AND PARAMETERS
     !   (common cvflag)  
237    
238        ! set simulation flags:
239        ! (common cvflag)
240      CALL cv_flag      CALL cv_flag
241    
242      ! -- set thermodynamical constants:      ! set thermodynamical constants:
243      !     (common cvthermo)      ! (common cvthermo)
244        CALL cv_thermo
245      CALL cv_thermo(iflag_con)  
246        ! set convect parameters
247      ! -- set convect parameters      ! includes microphysical parameters and parameters that
248      !      ! control the rate of approach to quasi-equilibrium)
249      !     includes microphysical parameters and parameters that      ! (common cvparam)
250      !     control the rate of approach to quasi-equilibrium)  
251      !     (common cvparam)      CALL cv3_param(klev, delt)
252    
253      if (iflag_con.eq.3) then      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
254         CALL cv3_param(nd, delt)  
255      endif      do k = 1, klev
256           do i = 1, klon
257      if (iflag_con.eq.4) then            ft1(i, k) = 0.0
258         CALL cv_param(nd)            fq1(i, k) = 0.0
259      endif            fu1(i, k) = 0.0
260              fv1(i, k) = 0.0
261      !---------------------------------------------------------------------            tvp1(i, k) = 0.0
262      ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS            tp1(i, k) = 0.0
263      !---------------------------------------------------------------------            clw1(i, k) = 0.0
264              clw(i, k) = 0.0
     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  
265            gz1(i, k) = 0.            gz1(i, k) = 0.
266            VPrecip1(i, k) = 0.            VPrecip1(i, k) = 0.
267            Ma1(i, k)=0.0            Ma1(i, k) = 0.0
268            upwd1(i, k)=0.0            upwd1(i, k) = 0.0
269            dnwd1(i, k)=0.0            dnwd1(i, k) = 0.0
270            dnwd01(i, k)=0.0            dnwd01(i, k) = 0.0
271            qcondc1(i, k)=0.0            qcondc1(i, k) = 0.0
        end do  
     end do  
   
     do  j=1, ntra  
        do  k=1, nd  
           do  i=1, len  
              ftra1(i, k, j)=0.0  
           end do  
272         end do         end do
273      end do      end do
274    
275      do  i=1, len      do i = 1, klon
276         precip1(i)=0.0         precip1(i) = 0.0
277         iflag1(i)=0         iflag1(i) = 0
278         wd1(i)=0.0         wd1(i) = 0.0
279         cape1(i)=0.0         cape1(i) = 0.0
280         VPrecip1(i, nd+1)=0.0         VPrecip1(i, klev + 1) = 0.0
281      end do      end do
282    
283      if (iflag_con.eq.3) then      do il = 1, klon
284         do il=1, len         sig1(il, klev) = sig1(il, klev) + 1.
285            sig1(il, nd)=sig1(il, nd)+1.         sig1(il, klev) = min(sig1(il, klev), 12.1)
286            sig1(il, nd)=amin1(sig1(il, nd), 12.1)      enddo
287         enddo  
288      endif      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
289        CALL cv3_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
290      !--------------------------------------------------------------------           gz1, h1, hm1, th1)
291      ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY  
292      !--------------------------------------------------------------------      ! CONVECTIVE FEED
293        CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
294      if (iflag_con.eq.3) then           icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
295         CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, &  
296              h1, hm1, th1)! nd->na      ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
297      endif      ! (up through ICB for convect4, up through ICB + 1 for convect3)
298        ! Calculates the lifted parcel virtual temperature at nk, the
299      if (iflag_con.eq.4) then      ! actual temperature, and the adiabatic liquid water content.
300         CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1 &      CALL cv3_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
301              , lv1, cpn1, tv1, gz1, h1, hm1)           tp1, tvp1, clw1, icbs1) ! klev->na
302      endif  
303        ! TRIGGERING
304      !--------------------------------------------------------------------      CALL cv3_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
305      ! --- CONVECTIVE FEED           buoybase1, iflag1, sig1, w01) ! klev->na
306      !--------------------------------------------------------------------  
307        ! Moist convective adjustment is necessary
308      if (iflag_con.eq.3) then  
309         CALL cv3_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1            &      ncum = 0
310              , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! nd->na      do i = 1, klon
311      endif         if (iflag1(i) == 0) then
312              ncum = ncum + 1
313      if (iflag_con.eq.4) then            idcum(ncum) = i
        CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1 &  
             , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)  
     endif  
   
     !--------------------------------------------------------------------  
     ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part  
     ! (up through ICB for convect4, up through ICB+1 for convect3)  
     !     Calculates the lifted parcel virtual temperature at nk, the  
     !     actual temperature, and the adiabatic liquid water content.  
     !--------------------------------------------------------------------  
   
     if (iflag_con.eq.3) then  
        CALL cv3_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1   &  
             , tp1, tvp1, clw1, icbs1) ! nd->na  
     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  
314         endif         endif
315      end do      end do
316    
317      !       print*, 'klon, ncum = ', len, ncum      IF (ncum > 0) THEN
318           ! COMPRESS THE FIELDS
319      IF (ncum.gt.0) THEN         ! (-> vectorization over convective gridpoints)
320           CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
321         !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^              plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
322         ! --- COMPRESS THE FIELDS              v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
323         !        (-> vectorization over convective gridpoints)              sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
324         !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^              buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
325                tvp, clw, sig, w0)
326         if (iflag_con.eq.3) then  
327            CALL cv3_compress( len, nloc, ncum, nd, ntra &         ! UNDILUTE (ADIABATIC) UPDRAFT / second part :
328                 , iflag1, nk1, icb1, icbs1 &         ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
329                 , plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1 &         ! &
330                 , t1, q1, qs1, u1, v1, gz1, th1 &         ! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
331                 , tra1 &         ! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
332                 , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1  &         ! &
333                 , sig1, w01 &         ! FIND THE LEVEL OF NEUTRAL BUOYANCY
334                 , iflag, nk, icb, icbs &         CALL cv3_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, &
335                 , plcl, tnk, qnk, gznk, pbase, buoybase &              t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, &
336                 , t, q, qs, u, v, gz, th &              tvp, clw, hp, ep, sigp, buoy) !na->klev
337                 , tra &  
338                 , h, lv, cpn, p, ph, tv, tp, tvp, clw  &         ! CLOSURE
339                 , sig, w0  )         CALL cv3_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
340         endif              buoy, sig, w0, cape, m) ! na->klev
341    
342         if (iflag_con.eq.4) then         ! MIXING
343            CALL cv_compress( len, nloc, ncum, nd &         CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
344                 , iflag1, nk1, icb1 &              v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
345                 , cbmf1, plcl1, tnk1, qnk1, gznk1 &              sij, elij, ments, qents)
346                 , t1, q1, qs1, u1, v1, gz1 &  
347                 , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1 &         ! UNSATURATED (PRECIPITATING) DOWNDRAFTS
348                 , iflag, nk, icb &         CALL cv3_unsat(klon, ncum, klev, klev, icb, inb, t, q, qs, gz, u, &
349                 , cbmf, plcl, tnk, qnk, gznk &              v, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, &
350                 , t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw  &              plcl, mp, qp, up, vp, wt, water, evap, b)! na->klev
351                 , dph )  
352         endif         ! YIELD
353           ! (tendencies, precipitation, variables of interface with other
354         !-------------------------------------------------------------------         ! processes, etc)
355         ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :         CALL cv3_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &
356         ! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES              gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &
357         ! ---   &              wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &
358         ! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE              tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &
359         ! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD              dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev
360         ! ---   &  
361         ! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY         ! passive tracers
362         !-------------------------------------------------------------------         CALL cv3_tracer(klon, ncum, klev, ment, sij, da, phi)
363    
364         if (iflag_con.eq.3) then         ! UNCOMPRESS THE FIELDS
365            CALL cv3_undilute2(nloc, ncum, nd, icb, icbs, nk         &  
366                 , tnk, qnk, gznk, t, q, qs, gz &         ! set iflag1 = 42 for non convective points
367                 , p, h, tv, lv, pbase, buoybase, plcl &         do i = 1, klon
368                 , inb, tp, tvp, clw, hp, ep, sigp, buoy) !na->nd            iflag1(i) = 42
        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, ntra, icb, nk, inb     &  
                , ph, t, q, qs, u, v, tra, h, lv, qnk &  
                , hp, tv, tvp, ep, clw, m, sig &  
                , ment, qent, uent, vent, nent, sij, elij, ments, qents, traent)! na->nd  
        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, ntra, icb, inb     &  
                , t, q, qs, gz, u, v, tra, p, ph &  
                , th, tv, lv, cpn, ep, sigp, clw &  
                , m, ment, elij, delt, plcl &  
                , mp, qp, up, vp, trap, 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, ntra             &  
                , icb, inb, delt &  
                , t, q, u, v, tra, gz, p, ph, h, hp, lv, cpn, th &  
                , ep, clw, m, tp, mp, qp, up, vp, trap &  
                , wt, water, evap, b &  
                , ment, qent, uent, vent, nent, elij, traent, sig &  
                , tv, tvp &  
                , iflag, precip, VPrecip, ft, fq, fu, fv, ftra &  
                , 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  
369         end do         end do
        !  
        if (iflag_con.eq.3) then  
           CALL cv3_uncompress(nloc, len, ncum, nd, ntra, idcum &  
                , iflag &  
                , precip, VPrecip, sig, w0 &  
                , ft, fq, fu, fv, ftra &  
                , inb  &  
                , Ma, upwd, dnwd, dnwd0, qcondc, wd, cape &  
                , da, phi, mp &  
                , iflag1 &  
                , precip1, VPrecip1, sig1, w01 &  
                , ft1, fq1, fu1, fv1, ftra1 &  
                , inb1 &  
                , Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1  &  
                , da1, phi1, mp1)  
        endif  
370    
371         if (iflag_con.eq.4) then         CALL cv3_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
372            CALL cv_uncompress(nloc, len, ncum, nd, idcum &              ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
373                 , iflag &              da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
374                 , precip, cbmf &              fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
375                 , ft, fq, fu, fv &              cape1, da1, phi1, mp1)
376                 , Ma, qcondc             &      ENDIF
                , iflag1 &  
                , precip1, cbmf1 &  
                , ft1, fq1, fu1, fv1 &  
                , Ma1, qcondc1 )  
        endif  
     ENDIF ! ncum>0  
377    
378    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver
379    

Legend:
Removed from v.52  
changed lines
  Added in v.183

  ViewVC Help
Powered by ViewVC 1.1.21