CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mfpblt.f
1
5
9 module mfpblt_mod
10 contains
16 subroutine mfpblt(im,ix,km,kmpbl,ntcw,ntrac1,delt, &
17 & cnvflg,zl,zm,q1,t1,u1,v1,plyr,pix,thlx,thvx, &
18 & gdx,hpbl,kpbl,vpert,buo,xmf, &
19 & tcko,qcko,ucko,vcko,xlamue)
20!
21 use machine , only : kind_phys
22 use funcphys , only : fpvs
23 use physcons, grav => con_g, cp => con_cp &
24 &, rv => con_rv, hvap => con_hvap &
25 &, fv => con_fvirt &
26 &, eps => con_eps, epsm1 => con_epsm1
27!
28 implicit none
29!
30 integer im, ix, km, kmpbl, ntcw, ntrac1
31! &, me
32 integer kpbl(im)
33 logical cnvflg(im)
34 real(kind=kind_phys) delt
35 real(kind=kind_phys) q1(ix,km,ntrac1), &
36 & t1(ix,km), u1(ix,km), v1(ix,km), &
37 & plyr(im,km),pix(im,km),thlx(im,km), &
38 & thvx(im,km),zl(im,km), zm(im,km), &
39 & gdx(im), &
40 & hpbl(im), vpert(im), &
41 & buo(im,km), xmf(im,km), &
42 & tcko(im,km),qcko(im,km,ntrac1), &
43 & ucko(im,km),vcko(im,km), &
44 & xlamue(im,km-1)
45!
46c local variables and arrays
47!
48 integer i, j, k, n, ndc
49 integer kpblx(im), kpbly(im)
50!
51 real(kind=kind_phys) dt2, dz, ce0, cm,
52 & factor, gocp,
53 & g, b1, f1,
54 & bb1, bb2,
55 & alp, a1, pgcon,
56 & qmin, qlmin, xmmx, rbint,
57 & tem, tem1, tem2,
58 & ptem, ptem1, ptem2
59!
60 real(kind=kind_phys) elocp, el2orc, qs, es,
61 & tlu, gamma, qlu,
62 & thup, thvu, dq
63!
64 real(kind=kind_phys) rbdn(im), rbup(im), hpblx(im),
65 & xlamuem(im,km-1)
66!
67 real(kind=kind_phys) wu2(im,km), thlu(im,km),
68 & qtx(im,km), qtu(im,km)
69!
70 real(kind=kind_phys) xlamavg(im), sigma(im),
71 & scaldfunc(im), sumx(im)
72!
73 logical totflg, flg(im)
74!
75! physical parameters
76 parameter(g=grav)
77 parameter(gocp=g/cp)
78 parameter(elocp=hvap/cp,el2orc=hvap*hvap/(rv*cp))
79 parameter(ce0=0.4,cm=1.0)
80 parameter(qmin=1.e-8,qlmin=1.e-12)
81 parameter(alp=1.0,pgcon=0.55)
82 parameter(a1=0.13,b1=0.5,f1=0.15)
83!
84!************************************************************************
85!!
86 totflg = .true.
87 do i=1,im
88 totflg = totflg .and. (.not. cnvflg(i))
89 enddo
90 if(totflg) return
91!
92 dt2 = delt
93!!
94 do k = 1, km
95 do i=1,im
96 if (cnvflg(i)) then
97 buo(i,k) = 0.
98 wu2(i,k) = 0.
99 qtx(i,k) = q1(i,k,1) + q1(i,k,ntcw)
100 endif
101 enddo
102 enddo
103!
105!
106 do i=1,im
107 if(cnvflg(i)) then
108 ptem = alp * vpert(i)
109 ptem = min(ptem, 3.0)
110 thlu(i,1)= thlx(i,1) + ptem
111 qtu(i,1) = qtx(i,1)
112 buo(i,1) = g * ptem / thvx(i,1)
113 endif
114 enddo
115!
117!
118 do k = 1, kmpbl
119 do i=1,im
120 if(cnvflg(i)) then
121 dz = zl(i,k+1) - zl(i,k)
122 if(k < kpbl(i)) then
123 ptem = 1./(zm(i,k)+dz)
124 tem = max((hpbl(i)-zm(i,k)+dz) ,dz)
125 ptem1 = 1./tem
126 xlamue(i,k) = ce0 * (ptem+ptem1)
127 else
128 xlamue(i,k) = ce0 / dz
129 endif
130 xlamuem(i,k) = cm * xlamue(i,k)
131 endif
132 enddo
133 enddo
134!
136!
137 do k = 2, kmpbl
138 do i=1,im
139 if(cnvflg(i)) then
140 dz = zl(i,k) - zl(i,k-1)
141 tem = 0.5 * xlamue(i,k-1) * dz
142 factor = 1. + tem
143!
144 thlu(i,k) = ((1.-tem)*thlu(i,k-1)+tem*
145 & (thlx(i,k-1)+thlx(i,k)))/factor
146 qtu(i,k) = ((1.-tem)*qtu(i,k-1)+tem*
147 & (qtx(i,k-1)+qtx(i,k)))/factor
148!
149 tlu = thlu(i,k) / pix(i,k)
150 es = 0.01 * fpvs(tlu) ! fpvs in pa
151 qs = max(qmin, eps * es / (plyr(i,k)+epsm1*es))
152 dq = qtu(i,k) - qs
153!
154 if (dq > 0.) then
155 gamma = el2orc * qs / (tlu**2)
156 qlu = dq / (1. + gamma)
157 qtu(i,k) = qs + qlu
158 tem1 = 1. + fv * qs - qlu
159 thup = thlu(i,k) + pix(i,k) * elocp * qlu
160 thvu = thup * tem1
161 else
162 tem1 = 1. + fv * qtu(i,k)
163 thvu = thlu(i,k) * tem1
164 endif
165 buo(i,k) = g * (thvu / thvx(i,k) - 1.)
166!
167 endif
168 enddo
169 enddo
170!
173!
174! tem = 1.-2.*f1
175! bb1 = 2. * b1 / tem
176! bb2 = 2. / tem
177! from Soares et al. (2004,QJRMS)
178! bb1 = 2.
179! bb2 = 4.
180!
181! from Bretherton et al. (2004, MWR)
182! bb1 = 4.
183! bb2 = 2.
184!
185! from our tuning
186 bb1 = 2.0
187 bb2 = 4.0
188!
189 do i = 1, im
190 if(cnvflg(i)) then
191 dz = zm(i,1)
192 tem = 0.5*bb1*xlamue(i,1)*dz
193 tem1 = bb2 * buo(i,1) * dz
194 ptem1 = 1. + tem
195 wu2(i,1) = tem1 / ptem1
196 endif
197 enddo
198 do k = 2, kmpbl
199 do i = 1, im
200 if(cnvflg(i)) then
201 dz = zm(i,k) - zm(i,k-1)
202 tem = 0.25*bb1*(xlamue(i,k)+xlamue(i,k-1))*dz
203 tem1 = bb2 * buo(i,k) * dz
204 ptem = (1. - tem) * wu2(i,k-1)
205 ptem1 = 1. + tem
206 wu2(i,k) = (ptem + tem1) / ptem1
207 endif
208 enddo
209 enddo
210!
212!
213 do i=1,im
214 flg(i) = .true.
215 kpbly(i) = kpbl(i)
216 if(cnvflg(i)) then
217 flg(i) = .false.
218 rbup(i) = wu2(i,1)
219 endif
220 enddo
221 do k = 2, kmpbl
222 do i = 1, im
223 if(.not.flg(i)) then
224 rbdn(i) = rbup(i)
225 rbup(i) = wu2(i,k)
226 kpblx(i)= k
227 flg(i) = rbup(i).le.0.
228 endif
229 enddo
230 enddo
231 do i = 1,im
232 if(cnvflg(i)) then
233 k = kpblx(i)
234 if(rbdn(i) <= 0.) then
235 rbint = 0.
236 elseif(rbup(i) >= 0.) then
237 rbint = 1.
238 else
239 rbint = rbdn(i)/(rbdn(i)-rbup(i))
240 endif
241 hpblx(i) = zm(i,k-1) + rbint*(zm(i,k)-zm(i,k-1))
242 endif
243 enddo
244!
245 do i = 1,im
246 if(cnvflg(i)) then
247 if(kpbl(i) > kpblx(i)) then
248 kpbl(i) = kpblx(i)
249 hpbl(i) = hpblx(i)
250 endif
251 endif
252 enddo
253!
255!
256 do k = 1, kmpbl
257 do i=1,im
258 if(cnvflg(i) .and. kpbly(i) > kpblx(i)) then
259 dz = zl(i,k+1) - zl(i,k)
260 if(k < kpbl(i)) then
261 ptem = 1./(zm(i,k)+dz)
262 tem = max((hpbl(i)-zm(i,k)+dz) ,dz)
263 ptem1 = 1./tem
264 xlamue(i,k) = ce0 * (ptem+ptem1)
265 else
266 xlamue(i,k) = ce0 / dz
267 endif
268 xlamuem(i,k) = cm * xlamue(i,k)
269 endif
270 enddo
271 enddo
272!
274!
275 do i = 1, im
276 xlamavg(i) = 0.
277 sumx(i) = 0.
278 enddo
279 do k = 1, kmpbl
280 do i = 1, im
281 if (cnvflg(i) .and. k < kpbl(i)) then
282 dz = zl(i,k+1) - zl(i,k)
283 xlamavg(i) = xlamavg(i) + xlamue(i,k) * dz
284 sumx(i) = sumx(i) + dz
285 endif
286 enddo
287 enddo
288 do i = 1, im
289 if(cnvflg(i)) then
290 xlamavg(i) = xlamavg(i) / sumx(i)
291 endif
292 enddo
293!
296!
297 do k = 1, kmpbl
298 do i = 1, im
299 if (cnvflg(i) .and. k < kpbl(i)) then
300 if(wu2(i,k) > 0.) then
301 tem = sqrt(wu2(i,k))
302 else
303 tem = 0.
304 endif
305 xmf(i,k) = a1 * tem
306 endif
307 enddo
308 enddo
309!
312!
313 do i = 1, im
314 if(cnvflg(i)) then
315 tem = 0.2 / xlamavg(i)
316 tem1 = 3.14 * tem * tem
317 sigma(i) = tem1 / (gdx(i) * gdx(i))
318 sigma(i) = max(sigma(i), 0.001)
319 sigma(i) = min(sigma(i), 0.999)
320 endif
321 enddo
322!
325!
326 do i = 1, im
327 if(cnvflg(i)) then
328 if (sigma(i) > a1) then
329 scaldfunc(i) = (1.-sigma(i)) * (1.-sigma(i))
330 scaldfunc(i) = max(min(scaldfunc(i), 1.0), 0.)
331 else
332 scaldfunc(i) = 1.0
333 endif
334 endif
335 enddo
336!
338!
339 do k = 1, kmpbl
340 do i = 1, im
341 if (cnvflg(i) .and. k < kpbl(i)) then
342 xmf(i,k) = scaldfunc(i) * xmf(i,k)
343 dz = zl(i,k+1) - zl(i,k)
344 xmmx = dz / dt2
345 xmf(i,k) = min(xmf(i,k),xmmx)
346 endif
347 enddo
348 enddo
349!
350!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
352!
353 do i=1,im
354 if(cnvflg(i)) then
355 thlu(i,1)= thlx(i,1)
356 endif
357 enddo
358!
359! do i=1,im
360! if(cnvflg(i)) then
361! ptem1 = max(qcko(i,1,ntcw), 0.)
362! tlu = thlu(i,1) / pix(i,1)
363! tcko(i,1) = tlu + elocp * ptem1
364! endif
365! enddo
366!
367 do k = 2, kmpbl
368 do i=1,im
369 if(cnvflg(i) .and. k <= kpbl(i)) then
370 dz = zl(i,k) - zl(i,k-1)
371 tem = 0.5 * xlamue(i,k-1) * dz
372 factor = 1. + tem
373!
374 thlu(i,k) = ((1.-tem)*thlu(i,k-1)+tem*
375 & (thlx(i,k-1)+thlx(i,k)))/factor
376 qtu(i,k) = ((1.-tem)*qtu(i,k-1)+tem*
377 & (qtx(i,k-1)+qtx(i,k)))/factor
378!
379 tlu = thlu(i,k) / pix(i,k)
380 es = 0.01 * fpvs(tlu) ! fpvs in pa
381 qs = max(qmin, eps * es / (plyr(i,k)+epsm1*es))
382 dq = qtu(i,k) - qs
383!
384 if (dq > 0.) then
385 gamma = el2orc * qs / (tlu**2)
386 qlu = dq / (1. + gamma)
387 qtu(i,k) = qs + qlu
388 qcko(i,k,1) = qs
389 qcko(i,k,ntcw) = qlu
390 tcko(i,k) = tlu + elocp * qlu
391 else
392 qcko(i,k,1) = qtu(i,k)
393 qcko(i,k,ntcw) = 0.
394 tcko(i,k) = tlu
395 endif
396!
397 endif
398 enddo
399 enddo
400!
401 do k = 2, kmpbl
402 do i = 1, im
403 if (cnvflg(i) .and. k <= kpbl(i)) then
404 dz = zl(i,k) - zl(i,k-1)
405 tem = 0.5 * xlamuem(i,k-1) * dz
406 factor = 1. + tem
407 ptem = tem + pgcon
408 ptem1= tem - pgcon
409 ucko(i,k) = ((1.-tem)*ucko(i,k-1)+ptem*u1(i,k)
410 & +ptem1*u1(i,k-1))/factor
411 vcko(i,k) = ((1.-tem)*vcko(i,k-1)+ptem*v1(i,k)
412 & +ptem1*v1(i,k-1))/factor
413 endif
414 enddo
415 enddo
416!
417 if(ntcw > 2) then
418!
419 do n = 2, ntcw-1
420 do k = 2, kmpbl
421 do i = 1, im
422 if (cnvflg(i) .and. k <= kpbl(i)) then
423 dz = zl(i,k) - zl(i,k-1)
424 tem = 0.5 * xlamue(i,k-1) * dz
425 factor = 1. + tem
426!
427 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
428 & (q1(i,k,n)+q1(i,k-1,n)))/factor
429 endif
430 enddo
431 enddo
432 enddo
433!
434 endif
435!
436 ndc = ntrac1 - ntcw
437!
438 if(ndc > 0) then
439!
440 do n = ntcw+1, ntrac1
441 do k = 2, kmpbl
442 do i = 1, im
443 if (cnvflg(i) .and. k <= kpbl(i)) then
444 dz = zl(i,k) - zl(i,k-1)
445 tem = 0.5 * xlamue(i,k-1) * dz
446 factor = 1. + tem
447!
448 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
449 & (q1(i,k,n)+q1(i,k-1,n)))/factor
450 endif
451 enddo
452 enddo
453 enddo
454!
455 endif
456!
457 return
458 end
459
460 end module mfpblt_mod
This module contains the subroutine that calculates mass flux and updraft parcel properties for therm...
Definition mfpblt.f:9