/[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 188 - (show annotations)
Tue Mar 22 16:31:39 2016 UTC (8 years, 2 months ago) by guez
File size: 11268 byte(s)
Removed argument ncum of cv30_unsat, arguments nloc, ncum, nd, na of cv30_yield.

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, nl
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, allocatable:: b(:, :) ! (ncum, nl)
186 real ft(klon, klev), fq(klon, klev)
187 real fu(klon, klev), fv(klon, klev)
188 real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
189 real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
190 real tps(klon, klev)
191 real precip(klon)
192 real VPrecip(klon, klev + 1)
193 real qcondc(klon, klev) ! cld
194 real wd(klon) ! gust
195
196 !-------------------------------------------------------------------
197
198 ! SET CONSTANTS AND PARAMETERS
199
200 ! set thermodynamical constants:
201 ! (common cvthermo)
202 CALL cv_thermo
203
204 ! set convect parameters
205 ! includes microphysical parameters and parameters that
206 ! control the rate of approach to quasi-equilibrium)
207 ! (common cvparam)
208 CALL cv30_param(delt)
209
210 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
211
212 do k = 1, klev
213 do i = 1, klon
214 ft1(i, k) = 0.
215 fq1(i, k) = 0.
216 fu1(i, k) = 0.
217 fv1(i, k) = 0.
218 tvp1(i, k) = 0.
219 tp1(i, k) = 0.
220 clw1(i, k) = 0.
221 clw(i, k) = 0.
222 gz1(i, k) = 0.
223 VPrecip1(i, k) = 0.
224 Ma1(i, k) = 0.
225 upwd1(i, k) = 0.
226 dnwd1(i, k) = 0.
227 dnwd01(i, k) = 0.
228 qcondc1(i, k) = 0.
229 end do
230 end do
231
232 do i = 1, klon
233 precip1(i) = 0.
234 iflag1(i) = 0
235 wd1(i) = 0.
236 cape1(i) = 0.
237 VPrecip1(i, klev + 1) = 0.
238 end do
239
240 do il = 1, klon
241 sig1(il, klev) = sig1(il, klev) + 1.
242 sig1(il, klev) = min(sig1(il, klev), 12.1)
243 enddo
244
245 ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
246 CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
247 gz1, h1, hm1, th1)
248
249 ! CONVECTIVE FEED
250 CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &
251 icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na
252
253 ! UNDILUTE (ADIABATIC) UPDRAFT / 1st part
254 ! (up through ICB for convect4, up through ICB + 1 for convect3)
255 ! Calculates the lifted parcel virtual temperature at nk, the
256 ! actual temperature, and the adiabatic liquid water content.
257 CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &
258 tp1, tvp1, clw1, icbs1) ! klev->na
259
260 ! TRIGGERING
261 CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
262 buoybase1, iflag1, sig1, w01) ! klev->na
263
264 ! Moist convective adjustment is necessary
265
266 ncum = 0
267 do i = 1, klon
268 if (iflag1(i) == 0) then
269 ncum = ncum + 1
270 idcum(ncum) = i
271 endif
272 end do
273
274 IF (ncum > 0) THEN
275 allocate(b(ncum, nl))
276
277 ! COMPRESS THE FIELDS
278 ! (-> vectorization over convective gridpoints)
279 CALL cv30_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, &
280 plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, &
281 v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, &
282 sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, &
283 buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, &
284 tvp, clw, sig, w0)
285
286 CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, &
287 h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, &
288 ep, sigp, buoy)
289
290 ! CLOSURE
291 CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &
292 buoy, sig, w0, cape, m) ! na->klev
293
294 ! MIXING
295 CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &
296 v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &
297 sij, elij, ments, qents)
298
299 ! Unsaturated (precipitating) downdrafts
300 CALL cv30_unsat(icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, ph, th, &
301 tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, plcl, mp, qp, &
302 up, vp, wt, water, evap, b)
303
304 ! Yield (tendencies, precipitation, variables of interface with
305 ! other processes, etc)
306 CALL cv30_yield(icb(:ncum), inb(:ncum), delt, t, q, u, v, gz, p, ph, &
307 h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, wt, &
308 water(:ncum, :nl), evap(:ncum, :nl), b, ment, qent, uent, vent, &
309 nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, &
310 upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, wd)
311
312 ! passive tracers
313 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
314
315 ! UNCOMPRESS THE FIELDS
316
317 ! set iflag1 = 42 for non convective points
318 iflag1 = 42
319
320 CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
321 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, &
322 da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
323 fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
324 cape1, da1, phi1, mp1)
325 ENDIF
326
327 end SUBROUTINE cv_driver
328
329 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21