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

Legend:
Removed from v.91  
changed lines
  Added in v.180

  ViewVC Help
Powered by ViewVC 1.1.21