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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 181 - (show annotations)
Tue Mar 15 17:51:30 2016 UTC (8 years, 1 month ago) by guez
File size: 13895 byte(s)
Removed the option iflag_con == 4. This seems to be a coding of the
Emanuel scheme equivalent to and older than the coding for iflag_con
== 3.

1 module cv_driver_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, &
8 fq1, fu1, fv1, precip1, VPrecip1, cbmf1, sig1, w01, icb1, inb1, delt, &
9 Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1, da1, phi1, mp1)
10
11 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12 ! Main driver for convection
13 ! Author: S. Bony, March 2002
14
15 ! Several modules corresponding to different physical processes
16
17 use cv3_compress_m, only: cv3_compress
18 use cv3_feed_m, only: cv3_feed
19 use cv3_mixing_m, only: cv3_mixing
20 use cv3_param_m, only: cv3_param
21 use cv3_prelim_m, only: cv3_prelim
22 use cv3_tracer_m, only: cv3_tracer
23 use cv3_uncompress_m, only: cv3_uncompress
24 use cv3_unsat_m, only: cv3_unsat
25 use cv3_yield_m, only: cv3_yield
26 USE dimphy, ONLY: klev, klon
27
28 real, intent(in):: t1(klon, klev) ! temperature
29 real, intent(in):: q1(klon, klev) ! specific hum
30 real, intent(in):: qs1(klon, klev) ! sat specific hum
31 real, intent(in):: u1(klon, klev) ! u-wind
32 real, intent(in):: v1(klon, klev) ! v-wind
33 real, intent(in):: p1(klon, klev) ! full level pressure
34 real, intent(in):: ph1(klon, klev + 1) ! half level pressure
35 integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions
36 real, intent(out):: ft1(klon, klev) ! temp tend
37 real, intent(out):: fq1(klon, klev) ! spec hum tend
38 real, intent(out):: fu1(klon, klev) ! u-wind tend
39 real, intent(out):: fv1(klon, klev) ! v-wind tend
40 real, intent(out):: precip1(klon) ! precipitation
41
42 real, intent(out):: VPrecip1(klon, klev + 1)
43 ! vertical profile of precipitation
44
45 real, intent(inout):: cbmf1(klon) ! cloud base mass flux
46 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 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
57 real, intent(out):: upwd1(klon, klev)
58 ! total upward mass flux (adiab + mixed)
59
60 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
61 real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
62
63 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
70 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
71 real, intent(inout):: mp1(klon, klev)
72
73 ! ARGUMENTS
74
75 ! On input:
76
77 ! 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
82 ! 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
87 ! 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
92 ! 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
97 ! v: Same as u but for meridional velocity.
98
99 ! 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
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 ! det: Array of detrainment mass flux of dimension KLEV.
168
169 ! Local arrays
170
171 real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
172
173 integer i, k, il
174 integer icbmax
175 integer nk1(klon)
176 integer icbs1(klon)
177
178 real plcl1(klon)
179 real tnk1(klon)
180 real qnk1(klon)
181 real gznk1(klon)
182 real pbase1(klon)
183 real buoybase1(klon)
184
185 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
196 integer ncum
197
198 ! (local) compressed fields:
199
200 integer idcum(klon)
201 integer iflag(klon), nk(klon), icb(klon)
202 integer nent(klon, klev)
203 integer icbs(klon)
204 integer inb(klon), inbis(klon)
205
206 real cbmf(klon), 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 dph(klon, klev)
213 real pbase(klon), buoybase(klon), th(klon, klev)
214 real tvp(klon, klev)
215 real sig(klon, klev), w0(klon, klev)
216 real hp(klon, klev), ep(klon, klev), sigp(klon, klev)
217 real frac(klon), buoy(klon, klev)
218 real cape(klon)
219 real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
220 real uent(klon, klev, klev), vent(klon, klev, klev)
221 real ments(klon, klev, klev), qents(klon, klev, klev)
222 real sij(klon, klev, klev), elij(klon, klev, klev)
223 real qp(klon, klev), up(klon, klev), vp(klon, klev)
224 real wt(klon, klev), water(klon, klev), evap(klon, klev)
225 real b(klon, klev), ft(klon, klev), fq(klon, klev)
226 real fu(klon, klev), fv(klon, klev)
227 real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
228 real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
229 real tps(klon, klev), qprime(klon), tprime(klon)
230 real precip(klon)
231 real VPrecip(klon, klev + 1)
232 real qcondc(klon, klev) ! cld
233 real wd(klon) ! gust
234
235 !-------------------------------------------------------------------
236
237 ! SET CONSTANTS AND PARAMETERS
238
239 ! set simulation flags:
240 ! (common cvflag)
241
242 CALL cv_flag
243
244 ! set thermodynamical constants:
245 ! (common cvthermo)
246
247 CALL cv_thermo
248
249 ! set convect parameters
250
251 ! includes microphysical parameters and parameters that
252 ! control the rate of approach to quasi-equilibrium)
253 ! (common cvparam)
254
255 CALL cv3_param(klev, delt)
256
257 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
258
259 do k = 1, klev
260 do i = 1, klon
261 ft1(i, k) = 0.0
262 fq1(i, k) = 0.0
263 fu1(i, k) = 0.0
264 fv1(i, k) = 0.0
265 tvp1(i, k) = 0.0
266 tp1(i, k) = 0.0
267 clw1(i, k) = 0.0
268 !ym
269 clw(i, k) = 0.0
270 gz1(i, k) = 0.
271 VPrecip1(i, k) = 0.
272 Ma1(i, k) = 0.0
273 upwd1(i, k) = 0.0
274 dnwd1(i, k) = 0.0
275 dnwd01(i, k) = 0.0
276 qcondc1(i, k) = 0.0
277 end do
278 end do
279
280 do i = 1, klon
281 precip1(i) = 0.0
282 iflag1(i) = 0
283 wd1(i) = 0.0
284 cape1(i) = 0.0
285 VPrecip1(i, klev + 1) = 0.0
286 end do
287
288 do il = 1, klon
289 sig1(il, klev) = sig1(il, klev) + 1.
290 sig1(il, klev) = min(sig1(il, klev), 12.1)
291 enddo
292
293 ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
294
295 CALL cv3_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
296 gz1, h1, hm1, th1)
297
298 ! CONVECTIVE FEED
299
300 CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
301 icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
302
303 ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
304 ! (up through ICB for convect4, up through ICB + 1 for convect3)
305 ! Calculates the lifted parcel virtual temperature at nk, the
306 ! actual temperature, and the adiabatic liquid water content.
307
308 CALL cv3_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
309 tp1, tvp1, clw1, icbs1) ! klev->na
310
311 ! TRIGGERING
312
313 CALL cv3_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
314 buoybase1, iflag1, sig1, w01) ! klev->na
315
316 ! Moist convective adjustment is necessary
317
318 ncum = 0
319 do i = 1, klon
320 if (iflag1(i) == 0) then
321 ncum = ncum + 1
322 idcum(ncum) = i
323 endif
324 end do
325
326 IF (ncum > 0) THEN
327 ! COMPRESS THE FIELDS
328 ! (-> vectorization over convective gridpoints)
329
330 CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
331 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
332 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
333 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
334 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
335 tvp, clw, sig, w0)
336
337 ! UNDILUTE (ADIABATIC) UPDRAFT / second part :
338 ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
339 ! &
340 ! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
341 ! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
342 ! &
343 ! FIND THE LEVEL OF NEUTRAL BUOYANCY
344
345 CALL cv3_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, &
346 t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, &
347 tvp, clw, hp, ep, sigp, buoy) !na->klev
348
349 ! CLOSURE
350
351 CALL cv3_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
352 buoy, sig, w0, cape, m) ! na->klev
353
354 ! MIXING
355
356 CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
357 v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
358 sij, elij, ments, qents)
359
360 ! UNSATURATED (PRECIPITATING) DOWNDRAFTS
361
362 CALL cv3_unsat(klon, ncum, klev, klev, icb, inb, t, q, qs, gz, u, &
363 v, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, &
364 plcl, mp, qp, up, vp, wt, water, evap, b)! na->klev
365
366 ! YIELD
367 ! (tendencies, precipitation, variables of interface with other
368 ! processes, etc)
369
370 CALL cv3_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &
371 gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &
372 wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &
373 tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &
374 dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev
375
376 ! passive tracers
377
378 CALL cv3_tracer(klon, ncum, klev, ment, sij, da, phi)
379
380 ! UNCOMPRESS THE FIELDS
381
382 ! set iflag1 = 42 for non convective points
383 do i = 1, klon
384 iflag1(i) = 42
385 end do
386
387 CALL cv3_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
388 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
389 da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
390 fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
391 cape1, da1, phi1, mp1)
392 ENDIF ! ncum>0
393
394 end SUBROUTINE cv_driver
395
396 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21