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