/[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 187 - (show annotations)
Mon Mar 21 18:01:02 2016 UTC (8 years, 2 months ago) by guez
File size: 11216 byte(s)
Made variable nl of module cv30_param_m a parameter. There was no
coding allowing it to change.

Removed arguments nloc and nd of cv30_undilute2, arguments nloc, nd
and na of cv30_unsat. Just use klon and klev directly (going for
clarity).

Removed the option cvflag_grav = f. This was a lot of redundant code,
probably obsolete, and cvflag_grav was initialized to true with no
provision for changing it (as in LMDZ).

In cv30_unsat, downdraft_loop started at i = nl + 1, but for i >= nl,
i > inb, so num1 = 0.

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, fq1, fu1, &
8 fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, &
9 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 cv30_closure_m, only: cv30_closure
18 use cv30_compress_m, only: cv30_compress
19 use cv30_feed_m, only: cv30_feed
20 use cv30_mixing_m, only: cv30_mixing
21 use cv30_param_m, only: cv30_param
22 use cv30_prelim_m, only: cv30_prelim
23 use cv30_tracer_m, only: cv30_tracer
24 use cv30_uncompress_m, only: cv30_uncompress
25 use cv30_undilute2_m, only: cv30_undilute2
26 use cv30_unsat_m, only: cv30_unsat
27 use cv30_yield_m, only: cv30_yield
28 USE dimphy, ONLY: klev, klon
29
30 real, intent(in):: t1(klon, klev)
31 ! temperature (K), with first index corresponding to lowest model
32 ! level
33
34 real, intent(in):: q1(klon, klev)
35 ! Specific humidity, with first index corresponding to lowest
36 ! model level. Must be defined at same grid levels as T1.
37
38 real, intent(in):: qs1(klon, klev)
39 ! Saturation specific humidity, with first index corresponding to
40 ! lowest model level. Must be defined at same grid levels as
41 ! T1.
42
43 real, intent(in):: u1(klon, klev), v1(klon, klev)
44 ! Zonal wind and meridional velocity (m/s), witth first index
45 ! corresponding with the lowest model level. Defined at same
46 ! levels as T1.
47
48 real, intent(in):: p1(klon, klev)
49 ! Full level pressure (mb) of dimension KLEV, with first index
50 ! corresponding to lowest model level. Must be defined at same
51 ! grid levels as T1.
52
53 real, intent(in):: ph1(klon, klev + 1)
54 ! Half level pressure (mb), with first index corresponding to
55 ! lowest level. These pressures are defined at levels intermediate
56 ! between those of P1, T1, Q1 and QS1. The first value of PH
57 ! should be greater than (i.e. at a lower level than) the first
58 ! value of the array P1.
59
60 integer, intent(out):: iflag1(klon)
61 ! Flag for Emanuel conditions.
62
63 ! 0: Moist convection occurs.
64
65 ! 1: Moist convection occurs, but a CFL condition on the
66 ! subsidence warming is violated. This does not cause the scheme
67 ! to terminate.
68
69 ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
70
71 ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
72
73 ! 4: No moist convection; atmosphere is not unstable
74
75 ! 6: No moist convection because ihmin le minorig.
76
77 ! 7: No moist convection because unreasonable parcel level
78 ! temperature or specific humidity.
79
80 ! 8: No moist convection: lifted condensation level is above the
81 ! 200 mb level.
82
83 ! 9: No moist convection: cloud base is higher then the level NL-1.
84
85 real, intent(out):: ft1(klon, klev)
86 ! Temperature tendency (K/s), defined at same grid levels as T1,
87 ! Q1, QS1 and P1.
88
89 real, intent(out):: fq1(klon, klev)
90 ! Specific humidity tendencies (s-1), defined at same grid levels
91 ! as T1, Q1, QS1 and P1.
92
93 real, intent(out):: fu1(klon, klev), fv1(klon, klev)
94 ! Forcing (tendency) of zonal and meridional velocity (m/s^2),
95 ! defined at same grid levels as T1.
96
97 real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
98
99 real, intent(out):: VPrecip1(klon, klev + 1)
100 ! vertical profile of convective precipitation (kg/m2/s)
101
102 real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
103
104 real, intent(inout):: w01(klon, klev)
105 ! vertical velocity within adiabatic updraft
106
107 integer, intent(out):: icb1(klon)
108 integer, intent(inout):: inb1(klon)
109 real, intent(in):: delt ! the model time step (sec) between calls
110
111 real Ma1(klon, klev) ! Output mass flux adiabatic updraft
112
113 real, intent(out):: upwd1(klon, klev)
114 ! total upward mass flux (adiab + mixed)
115
116 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
117 real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
118
119 real qcondc1(klon, klev) ! Output in-cld mixing ratio of condensed water
120
121 real wd1(klon) ! gust
122 ! Output downdraft velocity scale for surface fluxes
123 ! A convective downdraft velocity scale. For use in surface
124 ! flux parameterizations. See convect.ps file for details.
125
126 real cape1(klon) ! Output
127 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
128 real, intent(inout):: mp1(klon, klev)
129
130 ! Local:
131
132 real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
133
134 integer i, k, il
135 integer icbmax
136 integer nk1(klon)
137 integer icbs1(klon)
138
139 real plcl1(klon)
140 real tnk1(klon)
141 real qnk1(klon)
142 real gznk1(klon)
143 real pbase1(klon)
144 real buoybase1(klon)
145
146 real lv1(klon, klev)
147 real cpn1(klon, klev)
148 real tv1(klon, klev)
149 real gz1(klon, klev)
150 real hm1(klon, klev)
151 real h1(klon, klev)
152 real tp1(klon, klev)
153 real tvp1(klon, klev)
154 real clw1(klon, klev)
155 real th1(klon, klev)
156
157 integer ncum
158
159 ! Compressed fields:
160
161 integer idcum(klon)
162 integer iflag(klon), nk(klon), icb(klon)
163 integer nent(klon, klev)
164 integer icbs(klon)
165 integer inb(klon)
166
167 real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
168 real t(klon, klev), q(klon, klev), qs(klon, klev)
169 real u(klon, klev), v(klon, klev)
170 real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
171 real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
172 real clw(klon, klev)
173 real pbase(klon), buoybase(klon), th(klon, klev)
174 real tvp(klon, klev)
175 real sig(klon, klev), w0(klon, klev)
176 real hp(klon, klev), ep(klon, klev), sigp(klon, klev)
177 real buoy(klon, klev)
178 real cape(klon)
179 real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
180 real uent(klon, klev, klev), vent(klon, klev, klev)
181 real ments(klon, klev, klev), qents(klon, klev, klev)
182 real sij(klon, klev, klev), elij(klon, klev, klev)
183 real qp(klon, klev), up(klon, klev), vp(klon, klev)
184 real wt(klon, klev), water(klon, klev), evap(klon, klev)
185 real b(klon, klev), ft(klon, klev), fq(klon, klev)
186 real fu(klon, klev), fv(klon, klev)
187 real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
188 real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
189 real tps(klon, klev)
190 real precip(klon)
191 real VPrecip(klon, klev + 1)
192 real qcondc(klon, klev) ! cld
193 real wd(klon) ! gust
194
195 !-------------------------------------------------------------------
196
197 ! SET CONSTANTS AND PARAMETERS
198
199 ! set thermodynamical constants:
200 ! (common cvthermo)
201 CALL cv_thermo
202
203 ! set convect parameters
204 ! includes microphysical parameters and parameters that
205 ! control the rate of approach to quasi-equilibrium)
206 ! (common cvparam)
207 CALL cv30_param(delt)
208
209 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
210
211 do k = 1, klev
212 do i = 1, klon
213 ft1(i, k) = 0.
214 fq1(i, k) = 0.
215 fu1(i, k) = 0.
216 fv1(i, k) = 0.
217 tvp1(i, k) = 0.
218 tp1(i, k) = 0.
219 clw1(i, k) = 0.
220 clw(i, k) = 0.
221 gz1(i, k) = 0.
222 VPrecip1(i, k) = 0.
223 Ma1(i, k) = 0.
224 upwd1(i, k) = 0.
225 dnwd1(i, k) = 0.
226 dnwd01(i, k) = 0.
227 qcondc1(i, k) = 0.
228 end do
229 end do
230
231 do i = 1, klon
232 precip1(i) = 0.
233 iflag1(i) = 0
234 wd1(i) = 0.
235 cape1(i) = 0.
236 VPrecip1(i, klev + 1) = 0.
237 end do
238
239 do il = 1, klon
240 sig1(il, klev) = sig1(il, klev) + 1.
241 sig1(il, klev) = min(sig1(il, klev), 12.1)
242 enddo
243
244 ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
245 CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
246 gz1, h1, hm1, th1)
247
248 ! CONVECTIVE FEED
249 CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
250 icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
251
252 ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
253 ! (up through ICB for convect4, up through ICB + 1 for convect3)
254 ! Calculates the lifted parcel virtual temperature at nk, the
255 ! actual temperature, and the adiabatic liquid water content.
256 CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
257 tp1, tvp1, clw1, icbs1) ! klev->na
258
259 ! TRIGGERING
260 CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
261 buoybase1, iflag1, sig1, w01) ! klev->na
262
263 ! Moist convective adjustment is necessary
264
265 ncum = 0
266 do i = 1, klon
267 if (iflag1(i) == 0) then
268 ncum = ncum + 1
269 idcum(ncum) = i
270 endif
271 end do
272
273 IF (ncum > 0) THEN
274 ! COMPRESS THE FIELDS
275 ! (-> vectorization over convective gridpoints)
276 CALL cv30_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
277 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
278 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
279 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
280 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
281 tvp, clw, sig, w0)
282
283 CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, &
284 h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, &
285 ep, sigp, buoy)
286
287 ! CLOSURE
288 CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
289 buoy, sig, w0, cape, m) ! na->klev
290
291 ! MIXING
292 CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
293 v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
294 sij, elij, ments, qents)
295
296 ! Unsaturated (precipitating) downdrafts
297 CALL cv30_unsat(ncum, icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, &
298 ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, plcl, &
299 mp, qp, up, vp, wt, water, evap, b(:ncum, :))
300
301 ! Yield (tendencies, precipitation, variables of interface with
302 ! other processes, etc)
303 CALL cv30_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, &
304 gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, &
305 wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, &
306 tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, &
307 dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev
308
309 ! passive tracers
310 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
311
312 ! UNCOMPRESS THE FIELDS
313
314 ! set iflag1 = 42 for non convective points
315 iflag1 = 42
316
317 CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
318 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
319 da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
320 fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
321 cape1, da1, phi1, mp1)
322 ENDIF
323
324 end SUBROUTINE cv_driver
325
326 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21