CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mfpbltq.f
1
5
10 contains
16 subroutine mfpbltq(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,wush,tkemean,vez0fun,xmf,
19 & tcko,qcko,ucko,vcko,xlamueq,a1)
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), hpbl(im), vpert(im),
40 & buo(im,km), wush(im,km),
41 & tkemean(im),vez0fun(im),xmf(im,km),
42 & tcko(im,km),qcko(im,km,ntrac1),
43 & ucko(im,km),vcko(im,km),
44 & xlamueq(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,
52 & cm, cq, tkcrt,
53 & factor, gocp, cmxfac,
54 & g, b1, f1,
55 & bb1, bb2,
56 & alp, vpertmax,a1, pgcon,
57 & qmin, qlmin, xmmx, rbint,
58 & tem, tem1, tem2,
59 & ptem, ptem1, ptem2
60!
61 real(kind=kind_phys) elocp, el2orc, qs, es,
62 & tlu, gamma, qlu,
63 & thup, thvu, dq
64!
65 real(kind=kind_phys) rbdn(im), rbup(im), hpblx(im),
66 & xlamue(im,km-1), xlamuem(im,km-1)
67 real(kind=kind_phys) delz(im), xlamax(im), ce0t(im)
68!
69 real(kind=kind_phys) wu2(im,km), thlu(im,km),
70 & qtx(im,km), qtu(im,km)
71!
72 real(kind=kind_phys) xlamavg(im), sigma(im),
73 & scaldfunc(im), sumx(im)
74!
75 logical totflg, flg(im)
76!
77! physical parameters
78 parameter(g=grav)
79 parameter(gocp=g/cp)
80 parameter(elocp=hvap/cp,el2orc=hvap*hvap/(rv*cp))
81 parameter(ce0=0.4,cm=1.0,cq=1.0,tkcrt=2.,cmxfac=5.)
82 parameter(qmin=1.e-8,qlmin=1.e-12)
83 parameter(alp=1.5,vpertmax=3.0,pgcon=0.55)
84 parameter(b1=0.5,f1=0.15)
85!
86!************************************************************************
87!!
88 totflg = .true.
89 do i=1,im
90 totflg = totflg .and. (.not. cnvflg(i))
91 enddo
92 if(totflg) return
93!!
94!
95 dt2 = delt
96!
97 do k = 1, km
98 do i=1,im
99 if (cnvflg(i)) then
100 buo(i,k) = 0.
101 wu2(i,k) = 0.
102 qtx(i,k) = q1(i,k,1) + q1(i,k,ntcw)
103 endif
104 enddo
105 enddo
106!
108!
109 do i=1,im
110 if(cnvflg(i)) then
111 ptem = alp * vpert(i)
112 ptem = min(ptem, vpertmax)
113 thlu(i,1)= thlx(i,1) + ptem
114 qtu(i,1) = qtx(i,1)
115 buo(i,1) = g * ptem / thvx(i,1)
116 endif
117 enddo
118!
120!
121! if tkemean>tkcrt, ce0t=sqrt(tkemean/tkcrt)*ce0
122!
123 do i=1,im
124 if(cnvflg(i)) then
125 ce0t(i) = ce0 * vez0fun(i)
126 if(tkemean(i) > tkcrt) then
127 tem = sqrt(tkemean(i)/tkcrt)
128 tem1 = min(tem, cmxfac)
129 tem2 = tem1 * ce0
130 ce0t(i) = max(ce0t(i), tem2)
131 endif
132 endif
133 enddo
134!
135 do i=1,im
136 if(cnvflg(i)) then
137 k = kpbl(i) / 2
138 k = max(k, 1)
139 delz(i) = zl(i,k+1) - zl(i,k)
140 xlamax(i) = ce0t(i) / delz(i)
141 endif
142 enddo
143!
144 do k = 1, kmpbl
145 do i=1,im
146 if(cnvflg(i)) then
147 if(k < kpbl(i)) then
148 ptem = 1./(zm(i,k)+delz(i))
149 tem = max((hpbl(i)-zm(i,k)+delz(i)) ,delz(i))
150 ptem1 = 1./tem
151 xlamue(i,k) = ce0t(i) * (ptem+ptem1)
152 else
153 xlamue(i,k) = xlamax(i)
154 endif
155!
156 xlamueq(i,k) = cq * xlamue(i,k)
157 xlamuem(i,k) = cm * xlamue(i,k)
158 endif
159 enddo
160 enddo
161!
163!
164 do k = 2, kmpbl
165 do i=1,im
166 if(cnvflg(i)) then
167 dz = zl(i,k) - zl(i,k-1)
168 tem = 0.5 * xlamue(i,k-1) * dz
169 factor = 1. + tem
170!
171 thlu(i,k) = ((1.-tem)*thlu(i,k-1)+tem*
172 & (thlx(i,k-1)+thlx(i,k)))/factor
173!
174 tem = 0.5 * xlamueq(i,k-1) * dz
175 factor = 1. + tem
176 qtu(i,k) = ((1.-tem)*qtu(i,k-1)+tem*
177 & (qtx(i,k-1)+qtx(i,k)))/factor
178!
179 tlu = thlu(i,k) / pix(i,k)
180 es = 0.01 * fpvs(tlu) ! fpvs in pa
181 qs = max(qmin, eps * es / (plyr(i,k)+epsm1*es))
182 dq = qtu(i,k) - qs
183!
184 if (dq > 0.) then
185 gamma = el2orc * qs / (tlu**2)
186 qlu = dq / (1. + gamma)
187 qtu(i,k) = qs + qlu
188 tem1 = 1. + fv * qs - qlu
189 thup = thlu(i,k) + pix(i,k) * elocp * qlu
190 thvu = thup * tem1
191 else
192 tem1 = 1. + fv * qtu(i,k)
193 thvu = thlu(i,k) * tem1
194 endif
195 buo(i,k) = g * (thvu / thvx(i,k) - 1.)
196!
197 endif
198 enddo
199 enddo
200!
203!
204! tem = 1.-2.*f1
205! bb1 = 2. * b1 / tem
206! bb2 = 2. / tem
207! from Soares et al. (2004,QJRMS)
208! bb1 = 2.
209! bb2 = 4.
210!
211! from Bretherton et al. (2004, MWR)
212! bb1 = 4.
213! bb2 = 2.
214!
215! from our tuning
216 bb1 = 2.0
217 bb2 = 4.0
218!
219 do i = 1, im
220 if(cnvflg(i)) then
221 dz = zm(i,1)
222 tem = 0.5*bb1*xlamue(i,1)*dz
223 tem1 = bb2 * buo(i,1) * dz
224 ptem1 = 1. + tem
225 wu2(i,1) = tem1 / ptem1
226 endif
227 enddo
228 do k = 2, kmpbl
229 do i = 1, im
230 if(cnvflg(i)) then
231 dz = zm(i,k) - zm(i,k-1)
232 tem = 0.25*bb1*(xlamue(i,k-1)+xlamue(i,k))*dz
233 tem1 = max(wu2(i,k-1), 0.)
234 tem1 = bb2 * buo(i,k) - wush(i,k) * sqrt(tem1)
235 tem2 = tem1 * dz
236 ptem = (1. - tem) * wu2(i,k-1)
237 ptem1 = 1. + tem
238 wu2(i,k) = (ptem + tem2) / ptem1
239 endif
240 enddo
241 enddo
242!
244!
245 do i=1,im
246 flg(i) = .true.
247 kpblx(i) = 1
248 kpbly(i) = kpbl(i)
249 if(cnvflg(i)) then
250 flg(i) = .false.
251 rbup(i) = wu2(i,1)
252 endif
253 enddo
254 do k = 2, kmpbl
255 do i = 1, im
256 if(.not.flg(i)) then
257 rbdn(i) = rbup(i)
258 rbup(i) = wu2(i,k)
259 kpblx(i)= k
260 flg(i) = rbup(i).le.0.
261 endif
262 enddo
263 enddo
264 do i = 1,im
265 if(cnvflg(i)) then
266 k = kpblx(i)
267 if(rbdn(i) <= 0.) then
268 rbint = 0.
269 elseif(rbup(i) >= 0.) then
270 rbint = 1.
271 else
272 rbint = rbdn(i)/(rbdn(i)-rbup(i))
273 endif
274 hpblx(i) = zm(i,k-1) + rbint*(zm(i,k)-zm(i,k-1))
275 endif
276 enddo
277!
278 do i = 1,im
279 if(cnvflg(i)) then
280 if(kpblx(i) < kpbl(i)) then
281 kpbl(i) = kpblx(i)
282 hpbl(i) = hpblx(i)
283 endif
284 if(kpbl(i) <= 1) cnvflg(i)=.false.
285 endif
286 enddo
287!
289!
290 do i=1,im
291 if(cnvflg(i)) then
292 k = kpbl(i) / 2
293 k = max(k, 1)
294 delz(i) = zl(i,k+1) - zl(i,k)
295 xlamax(i) = ce0t(i) / delz(i)
296 endif
297 enddo
298!
299 do k = 1, kmpbl
300 do i=1,im
301 if(cnvflg(i) .and. kpblx(i) < kpbly(i)) then
302! if(cnvflg(i)) then
303 if(k < kpbl(i)) then
304 ptem = 1./(zm(i,k)+delz(i))
305 tem = max((hpbl(i)-zm(i,k)+delz(i)) ,delz(i))
306 ptem1 = 1./tem
307 xlamue(i,k) = ce0t(i) * (ptem+ptem1)
308 else
309 xlamue(i,k) = xlamax(i)
310 endif
311!
312 xlamueq(i,k) = cq * xlamue(i,k)
313 xlamuem(i,k) = cm * xlamue(i,k)
314 endif
315 enddo
316 enddo
317!
319!
320 do i = 1, im
321 xlamavg(i) = 0.
322 sumx(i) = 0.
323 enddo
324 do k = 1, kmpbl
325 do i = 1, im
326 if (cnvflg(i) .and. k < kpbl(i)) then
327 dz = zl(i,k+1) - zl(i,k)
328 xlamavg(i) = xlamavg(i) + xlamue(i,k) * dz
329 sumx(i) = sumx(i) + dz
330 endif
331 enddo
332 enddo
333 do i = 1, im
334 if(cnvflg(i)) then
335 xlamavg(i) = xlamavg(i) / sumx(i)
336 endif
337 enddo
338!
340!
341 do k = 1, kmpbl
342 do i = 1, im
343 if (cnvflg(i) .and. k < kpbl(i)) then
344 xmf(i,k) = a1 * sqrt(wu2(i,k))
345 endif
346 enddo
347 enddo
348!
351!
352 do i = 1, im
353 if(cnvflg(i)) then
354 tem = 0.2 / xlamavg(i)
355 tem1 = 3.14 * tem * tem
356 sigma(i) = tem1 / (gdx(i) * gdx(i))
357 sigma(i) = max(sigma(i), 0.001)
358 sigma(i) = min(sigma(i), 0.999)
359 endif
360 enddo
361!
364!
365 do i = 1, im
366 if(cnvflg(i)) then
367 if (sigma(i) > a1) then
368 scaldfunc(i) = (1.-sigma(i)) * (1.-sigma(i))
369 scaldfunc(i) = max(min(scaldfunc(i), 1.0), 0.)
370 else
371 scaldfunc(i) = 1.0
372 endif
373 endif
374 enddo
375!
377!
378 do k = 1, kmpbl
379 do i = 1, im
380 if (cnvflg(i) .and. k < kpbl(i)) then
381 xmf(i,k) = scaldfunc(i) * xmf(i,k)
382 dz = zl(i,k+1) - zl(i,k)
383 xmmx = dz / dt2
384 xmf(i,k) = min(xmf(i,k),xmmx)
385 endif
386 enddo
387 enddo
388!
389!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
391!
392 do i=1,im
393 if(cnvflg(i)) then
394 thlu(i,1)= thlx(i,1)
395 endif
396 enddo
397!
398! do i=1,im
399! if(cnvflg(i)) then
400! ptem1 = max(qcko(i,1,ntcw), 0.)
401! tlu = thlu(i,1) / pix(i,1)
402! tcko(i,1) = tlu + elocp * ptem1
403! endif
404! enddo
405!
406 do k = 2, kmpbl
407 do i=1,im
408 if(cnvflg(i) .and. k <= kpbl(i)) then
409 dz = zl(i,k) - zl(i,k-1)
410 tem = 0.5 * xlamue(i,k-1) * dz
411 factor = 1. + tem
412!
413 thlu(i,k) = ((1.-tem)*thlu(i,k-1)+tem*
414 & (thlx(i,k-1)+thlx(i,k)))/factor
415!
416 tem = 0.5 * xlamueq(i,k-1) * dz
417 factor = 1. + tem
418 qtu(i,k) = ((1.-tem)*qtu(i,k-1)+tem*
419 & (qtx(i,k-1)+qtx(i,k)))/factor
420!
421 tlu = thlu(i,k) / pix(i,k)
422 es = 0.01 * fpvs(tlu) ! fpvs in pa
423 qs = max(qmin, eps * es / (plyr(i,k)+epsm1*es))
424 dq = qtu(i,k) - qs
425!
426 if (dq > 0.) then
427 gamma = el2orc * qs / (tlu**2)
428 qlu = dq / (1. + gamma)
429 qtu(i,k) = qs + qlu
430 qcko(i,k,1) = qs
431 qcko(i,k,ntcw) = qlu
432 tcko(i,k) = tlu + elocp * qlu
433 else
434 qcko(i,k,1) = qtu(i,k)
435 qcko(i,k,ntcw) = 0.
436 tcko(i,k) = tlu
437 endif
438!
439 endif
440 enddo
441 enddo
442!
443 do k = 2, kmpbl
444 do i = 1, im
445 if (cnvflg(i) .and. k <= kpbl(i)) then
446 dz = zl(i,k) - zl(i,k-1)
447 tem = 0.5 * xlamuem(i,k-1) * dz
448 factor = 1. + tem
449 ptem = tem + pgcon
450 ptem1= tem - pgcon
451 ucko(i,k) = ((1.-tem)*ucko(i,k-1)+ptem*u1(i,k)
452 & +ptem1*u1(i,k-1))/factor
453 vcko(i,k) = ((1.-tem)*vcko(i,k-1)+ptem*v1(i,k)
454 & +ptem1*v1(i,k-1))/factor
455 endif
456 enddo
457 enddo
458!
459 if(ntcw > 2) then
460!
461 do n = 2, ntcw-1
462 do k = 2, kmpbl
463 do i = 1, im
464 if (cnvflg(i) .and. k <= kpbl(i)) then
465 dz = zl(i,k) - zl(i,k-1)
466 tem = 0.5 * xlamueq(i,k-1) * dz
467 factor = 1. + tem
468!
469 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
470 & (q1(i,k,n)+q1(i,k-1,n)))/factor
471 endif
472 enddo
473 enddo
474 enddo
475!
476 endif
477!
478 ndc = ntrac1 - ntcw
479!
480 if(ndc > 0) then
481!
482 do n = ntcw+1, ntrac1
483 do k = 2, kmpbl
484 do i = 1, im
485 if (cnvflg(i) .and. k <= kpbl(i)) then
486 dz = zl(i,k) - zl(i,k-1)
487 tem = 0.5 * xlamueq(i,k-1) * dz
488 factor = 1. + tem
489!
490 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
491 & (q1(i,k,n)+q1(i,k-1,n)))/factor
492 endif
493 enddo
494 enddo
495 enddo
496!
497 endif
498!
499 return
500 end
501
502 end module mfpbltq_mod
subroutine mfpbltq(im, ix, km, kmpbl, ntcw, ntrac1, delt,
This subroutine computes mass flux and updraft parcel properties for thermals driven by surface heati...
Definition mfpbltq.f:17
This module contains the subroutine that calculates mass flux and updraft parcel properties for therm...
Definition mfpbltq.f:9