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

Annotation of /trunk/Sources/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 186 - (hide annotations)
Mon Mar 21 15:36:26 2016 UTC (8 years, 1 month ago) by guez
File size: 13637 byte(s)
Removed variables nlm and nlp of module cv30_param_m. We do not
believe much in the benefit of these intermediary variables so we go
for clarity.

Removed variable noff of module cv30_param_m. Never used anywhere
else. Just set the value of nl explicitly in cv30_param.

Removed argument nd of cv30_param. Only called with nd = klev.

Replaced calls to zilch by array assignments. There was a strange
double call to zilch with the same arguments in cv30_mixing.

Removed procedure cv_flag. Just set the value of variable cvflag_grav
of module cvflag at declaration.

1 guez 52 module cv_driver_m
2 guez 3
3 guez 52 implicit none
4 guez 3
5 guez 52 contains
6 guez 3
7 guez 183 SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8     fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, &
9     dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
10 guez 3
11 guez 91 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12 guez 69 ! Main driver for convection
13 guez 97 ! Author: S. Bony, March 2002
14 guez 69
15 guez 72 ! Several modules corresponding to different physical processes
16    
17 guez 185 use cv30_compress_m, only: cv30_compress
18     use cv30_feed_m, only: cv30_feed
19     use cv30_mixing_m, only: cv30_mixing
20     use cv30_param_m, only: cv30_param
21     use cv30_prelim_m, only: cv30_prelim
22     use cv30_tracer_m, only: cv30_tracer
23     use cv30_uncompress_m, only: cv30_uncompress
24     use cv30_undilute2_m, only: cv30_undilute2
25     use cv30_unsat_m, only: cv30_unsat
26     use cv30_yield_m, only: cv30_yield
27 guez 62 USE dimphy, ONLY: klev, klon
28    
29 guez 103 real, intent(in):: t1(klon, klev) ! temperature
30     real, intent(in):: q1(klon, klev) ! specific hum
31     real, intent(in):: qs1(klon, klev) ! sat specific hum
32     real, intent(in):: u1(klon, klev) ! u-wind
33     real, intent(in):: v1(klon, klev) ! v-wind
34     real, intent(in):: p1(klon, klev) ! full level pressure
35     real, intent(in):: ph1(klon, klev + 1) ! half level pressure
36     integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions
37     real, intent(out):: ft1(klon, klev) ! temp tend
38     real, intent(out):: fq1(klon, klev) ! spec hum tend
39     real, intent(out):: fu1(klon, klev) ! u-wind tend
40     real, intent(out):: fv1(klon, klev) ! v-wind tend
41     real, intent(out):: precip1(klon) ! precipitation
42    
43 guez 180 real, intent(out):: VPrecip1(klon, klev + 1)
44 guez 103 ! vertical profile of precipitation
45    
46 guez 72 real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
47    
48     real, intent(inout):: w01(klon, klev)
49     ! vertical velocity within adiabatic updraft
50    
51 guez 103 integer, intent(out):: icb1(klon)
52     integer, intent(inout):: inb1(klon)
53     real, intent(in):: delt ! time step
54     real Ma1(klon, klev)
55     ! Ma1 Real Output mass flux adiabatic updraft
56 guez 180
57     real, intent(out):: upwd1(klon, klev)
58     ! total upward mass flux (adiab + mixed)
59    
60 guez 103 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
61     real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
62 guez 3
63 guez 103 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     real cape1(klon)
68     ! cape1 Real Output CAPE
69 guez 3
70 guez 103 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
71     real, intent(inout):: mp1(klon, klev)
72 guez 3
73 guez 180 ! ARGUMENTS
74 guez 103
75 guez 180 ! On input:
76 guez 62
77 guez 103 ! t: Array of absolute temperature (K) of dimension KLEV, with first
78     ! index corresponding to lowest model level. Note that this array
79     ! will be altered by the subroutine if dry convective adjustment
80     ! occurs and if IPBL is not equal to 0.
81 guez 62
82 guez 103 ! q: Array of specific humidity (gm/gm) of dimension KLEV, with first
83     ! index corresponding to lowest model level. Must be defined
84     ! at same grid levels as T. Note that this array will be altered
85     ! if dry convective adjustment occurs and if IPBL is not equal to 0.
86 guez 62
87 guez 103 ! qs: Array of saturation specific humidity of dimension KLEV, with first
88     ! index corresponding to lowest model level. Must be defined
89     ! at same grid levels as T. Note that this array will be altered
90     ! if dry convective adjustment occurs and if IPBL is not equal to 0.
91 guez 62
92 guez 103 ! u: Array of zonal wind velocity (m/s) of dimension KLEV, witth first
93     ! index corresponding with the lowest model level. Defined at
94     ! same levels as T. Note that this array will be altered if
95     ! dry convective adjustment occurs and if IPBL is not equal to 0.
96 guez 62
97 guez 103 ! v: Same as u but for meridional velocity.
98 guez 62
99 guez 103 ! p: Array of pressure (mb) of dimension KLEV, with first
100     ! index corresponding to lowest model level. Must be defined
101     ! at same grid levels as T.
102 guez 62
103 guez 180 ! ph: Array of pressure (mb) of dimension KLEV + 1, with first index
104 guez 103 ! 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 guez 62
109 guez 103 ! 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 guez 62
112 guez 103 ! delt: The model time step (sec) between calls to CONVECT
113 guez 62
114 guez 180 ! On Output:
115 guez 62
116 guez 103 ! 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 guez 62
135 guez 103 ! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same
136     ! grid levels as T, Q, QS and P.
137 guez 62
138 guez 103 ! 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 guez 62
141 guez 103 ! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV,
142     ! defined at same grid levels as T.
143 guez 62
144 guez 103 ! fv: Same as FU, but for forcing of meridional velocity.
145 guez 62
146 guez 103 ! precip: Scalar convective precipitation rate (mm/day).
147 guez 62
148 guez 103 ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
149 guez 62
150 guez 103 ! wd: A convective downdraft velocity scale. For use in surface
151     ! flux parameterizations. See convect.ps file for details.
152 guez 62
153 guez 103 ! tprime: A convective downdraft temperature perturbation scale (K).
154     ! For use in surface flux parameterizations. See convect.ps
155     ! file for details.
156 guez 62
157 guez 103 ! 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 guez 62
162 guez 103 ! 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 guez 62
167 guez 103 ! det: Array of detrainment mass flux of dimension KLEV.
168 guez 62
169 guez 103 ! Local arrays
170 guez 62
171 guez 103 real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
172 guez 3
173 guez 97 integer i, k, il
174 guez 52 integer icbmax
175     integer nk1(klon)
176     integer icbs1(klon)
177 guez 3
178 guez 52 real plcl1(klon)
179     real tnk1(klon)
180     real qnk1(klon)
181     real gznk1(klon)
182     real pbase1(klon)
183     real buoybase1(klon)
184 guez 3
185 guez 52 real lv1(klon, klev)
186     real cpn1(klon, klev)
187     real tv1(klon, klev)
188     real gz1(klon, klev)
189     real hm1(klon, klev)
190     real h1(klon, klev)
191     real tp1(klon, klev)
192     real tvp1(klon, klev)
193     real clw1(klon, klev)
194     real th1(klon, klev)
195 guez 62
196 guez 52 integer ncum
197 guez 62
198 guez 52 ! (local) compressed fields:
199 guez 62
200 guez 103 integer idcum(klon)
201     integer iflag(klon), nk(klon), icb(klon)
202     integer nent(klon, klev)
203     integer icbs(klon)
204 guez 183 integer inb(klon)
205 guez 3
206 guez 182 real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
207 guez 103 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 guez 180 real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
211 guez 103 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 guez 183 real buoy(klon, klev)
217 guez 103 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 guez 183 real tps(klon, klev)
229 guez 103 real precip(klon)
230 guez 180 real VPrecip(klon, klev + 1)
231 guez 103 real qcondc(klon, klev) ! cld
232     real wd(klon) ! gust
233 guez 3
234 guez 52 !-------------------------------------------------------------------
235 guez 3
236 guez 180 ! SET CONSTANTS AND PARAMETERS
237    
238     ! set thermodynamical constants:
239 guez 103 ! (common cvthermo)
240 guez 69 CALL cv_thermo
241 guez 52
242 guez 180 ! set convect parameters
243 guez 103 ! includes microphysical parameters and parameters that
244     ! control the rate of approach to quasi-equilibrium)
245     ! (common cvparam)
246 guez 52
247 guez 186 CALL cv30_param(delt)
248 guez 52
249 guez 180 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
250 guez 3
251 guez 103 do k = 1, klev
252     do i = 1, klon
253 guez 91 ft1(i, k) = 0.0
254     fq1(i, k) = 0.0
255     fu1(i, k) = 0.0
256     fv1(i, k) = 0.0
257     tvp1(i, k) = 0.0
258     tp1(i, k) = 0.0
259     clw1(i, k) = 0.0
260     clw(i, k) = 0.0
261 guez 103 gz1(i, k) = 0.
262 guez 52 VPrecip1(i, k) = 0.
263 guez 91 Ma1(i, k) = 0.0
264     upwd1(i, k) = 0.0
265     dnwd1(i, k) = 0.0
266     dnwd01(i, k) = 0.0
267     qcondc1(i, k) = 0.0
268 guez 52 end do
269     end do
270 guez 3
271 guez 103 do i = 1, klon
272 guez 91 precip1(i) = 0.0
273     iflag1(i) = 0
274     wd1(i) = 0.0
275     cape1(i) = 0.0
276 guez 180 VPrecip1(i, klev + 1) = 0.0
277 guez 52 end do
278 guez 3
279 guez 181 do il = 1, klon
280     sig1(il, klev) = sig1(il, klev) + 1.
281     sig1(il, klev) = min(sig1(il, klev), 12.1)
282     enddo
283 guez 3
284 guez 180 ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
285 guez 185 CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
286 guez 181 gz1, h1, hm1, th1)
287 guez 3
288 guez 180 ! CONVECTIVE FEED
289 guez 185 CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
290 guez 181 icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
291 guez 3
292 guez 180 ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
293     ! (up through ICB for convect4, up through ICB + 1 for convect3)
294 guez 103 ! Calculates the lifted parcel virtual temperature at nk, the
295     ! actual temperature, and the adiabatic liquid water content.
296 guez 185 CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
297 guez 181 tp1, tvp1, clw1, icbs1) ! klev->na
298 guez 3
299 guez 180 ! TRIGGERING
300 guez 185 CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
301 guez 181 buoybase1, iflag1, sig1, w01) ! klev->na
302 guez 3
303 guez 180 ! Moist convective adjustment is necessary
304 guez 3
305 guez 91 ncum = 0
306 guez 103 do i = 1, klon
307 guez 180 if (iflag1(i) == 0) then
308     ncum = ncum + 1
309 guez 91 idcum(ncum) = i
310 guez 52 endif
311     end do
312 guez 3
313 guez 103 IF (ncum > 0) THEN
314 guez 180 ! COMPRESS THE FIELDS
315 guez 103 ! (-> vectorization over convective gridpoints)
316 guez 185 CALL cv30_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
317 guez 181 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
318     v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
319     sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
320     buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
321     tvp, clw, sig, w0)
322 guez 3
323 guez 186 ! Undilute (adiabatic) updraft, second part: find the rest of
324     ! the lifted parcel temperatures; compute the precipitation
325     ! efficiencies and the fraction of precipitation falling
326     ! outside of cloud; find the level of neutral buoyancy.
327 guez 185 CALL cv30_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, &
328 guez 181 t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, &
329     tvp, clw, hp, ep, sigp, buoy) !na->klev
330 guez 3
331 guez 180 ! CLOSURE
332 guez 185 CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
333 guez 181 buoy, sig, w0, cape, m) ! na->klev
334 guez 3
335 guez 180 ! MIXING
336 guez 185 CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
337 guez 181 v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
338     sij, elij, ments, qents)
339 guez 3
340 guez 186 ! Unsaturated (precipitating) downdrafts
341     CALL cv30_unsat(klon, ncum, klev, klev, icb(:ncum), inb(:ncum), t, q, &
342     qs, gz, u, v, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, &
343     elij, delt, plcl, mp, qp, up, vp, wt, water, evap, b)! na->klev
344 guez 3
345 guez 186 ! Yield (tendencies, precipitation, variables of interface with
346     ! other processes, etc)
347 guez 185 CALL cv30_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &
348 guez 181 gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &
349     wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &
350     tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &
351     dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev
352 guez 3
353 guez 180 ! passive tracers
354 guez 185 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
355 guez 3
356 guez 180 ! UNCOMPRESS THE FIELDS
357 guez 103
358     ! set iflag1 = 42 for non convective points
359 guez 186 iflag1 = 42
360 guez 62
361 guez 185 CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
362 guez 181 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
363     da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
364     fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
365     cape1, da1, phi1, mp1)
366 guez 183 ENDIF
367 guez 3
368 guez 52 end SUBROUTINE cv_driver
369 guez 3
370 guez 52 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21