38 subroutine mfpbl(im,ix,km,ntrac,delt,cnvflg, &
39 & zl,zm,thvx,q1,t1,u1,v1,hpbl,kpbl, &
40 & sflx,ustar,wstar,xmf,tcko,qcko,ucko,vcko)
42 use machine
, only : kind_phys
47 integer im, ix, km, ntrac
51 real(kind=kind_phys) delt
52 real(kind=kind_phys) q1(ix,km,ntrac), t1(ix,km), &
53 & u1(ix,km), v1(ix,km), &
55 & zl(im,km), zm(im,km+1), &
56 & hpbl(im), sflx(im), ustar(im), &
57 & wstar(im), xmf(im,km), &
58 & tcko(im,km),qcko(im,km,ntrac), &
59 & ucko(im,km),vcko(im,km)
61 c local variables and arrays
63 integer i, j, k, n, kmpbl
65 real(kind=kind_phys) dt2, dz, ce0,
69 & alp, a1, qmin, zfmin,
76 real(kind=kind_phys) sigw1(im), usws3(im), xlamax(im),
77 & rbdn(im), rbup(im), delz(im)
79 real(kind=kind_phys) wu2(im,km), xlamue(im,km),
80 & thvu(im,km), zi(im,km),
83 logical totflg, flg(im)
89 parameter(ce0=0.38,qmin=1.e-8,alp=1.0,pgcon=0.55)
90 parameter(a1=0.08,b1=0.5,f1=0.15,c1=0.3,d1=2.58,tau=500.)
91 parameter(zfmin=1.e-8,h1=0.33333333)
93 c-----------------------------------------------------------------------
102 totflg = totflg .and. (.not. cnvflg(i))
119 delz(i) = zl(i,k+1) - zl(i,k)
120 xlamax(i) = ce0 / delz(i)
127 ptem = 1./(zi(i,k)+delz(i))
128 tem = max((hpbl(i)-zi(i,k)+delz(i)) ,delz(i))
130 xlamue(i,k) = ce0 * (ptem+ptem1)
132 xlamue(i,k) = xlamax(i)
138 c compute thermal excess
143 tem = zl(i,1)/hpbl(i)
144 usws3(i) = (ustar(i)/wstar(i))**3.
145 tem1 = usws3(i) + 0.6*tem
146 tem2 = max((1.-tem), zfmin)
147 ptem = (tem1**h1) * sqrt(tem2)
148 sigw1(i) = 1.3 * ptem * wstar(i)
149 ptem1 = alp * sflx(i) / sigw1(i)
150 thvu(i,1) = thvx(i,1) + ptem1
151 buo(i,1) = g * (thvu(i,1)/thvx(i,1)-1.)
155 c compute potential temperature and buoyancy for updraft air parcel
165 dz = zl(i,k) - zl(i,k-1)
166 tem = xlamue(i,k-1) * dz
168 ptem1 = (2. - tem) / ptem
169 tem1 = tem * (thvx(i,k)+thvx(i,k-1)) / ptem
170 thvu(i,k) = ptem1 * thvu(i,k-1) + tem1
171 buo(i,k) = g * (thvu(i,k)/thvx(i,k)-1.)
176 c compute updraft velocity square(wu2)
209 tem = 0.5*bb1*xlamue(i,1)*dz
210 tem1 = bb2 * buo(i,1) * dz
212 wu2(i,1) = tem1 / ptem1
219 dz = zi(i,k) - zi(i,k-1)
220 tem = 0.25*bb1*(xlamue(i,k)+xlamue(i,k-1))*dz
221 tem1 = bb2 * buo(i,k) * dz
222 ptem = (1. - tem) * wu2(i,k-1)
224 wu2(i,k) = (ptem + tem1) / ptem1
229 c update pbl height as the height
where updraft velocity vanishes
246 flg(i) = rbup(i).le.0.
253 if(rbdn(i) <= 0.)
then
255 elseif(rbup(i) >= 0.)
then
258 rbint = rbdn(i)/(rbdn(i)-rbup(i))
260 hpbl(i) = zi(i,k-1) + rbint*(zi(i,k)-zi(i,k-1))
269 delz(i) = zl(i,k+1) - zl(i,k)
270 xlamax(i) = ce0 / delz(i)
274 c update entrainment rate
295 ptem = 1./(zi(i,k)+delz(i))
296 tem = max((hpbl(i)-zi(i,k)+delz(i)) ,delz(i))
298 xlamue(i,k) = ce0 * (ptem+ptem1)
300 xlamue(i,k) = xlamax(i)
306 c updraft mass flux as a
function of sigmaw
307 c (0.3*sigmaw[square root of vertical turbulence variance])
323 c updraft mass flux as a
function of updraft velocity profile
332 if (cnvflg(i) .and. k < kpbl(i))
then
333 xmf(i,k) = a1 * sqrt(wu2(i,k))
334 dz = zl(i,k+1) - zl(i,k)
336 xmf(i,k) = min(xmf(i,k),xmmx)
342 c compute updraft property
359 if (cnvflg(i) .and. k <= kpbl(i))
then
360 dz = zl(i,k) - zl(i,k-1)
361 tem = 0.5 * xlamue(i,k-1) * dz
366 tcko(i,k) = ((1.-tem)*tcko(i,k-1)+tem*
367 & (t1(i,k)+t1(i,k-1))-gocp*dz)/factor
368 ucko(i,k) = ((1.-tem)*ucko(i,k-1)+ptem*u1(i,k)
369 & +ptem1*u1(i,k-1))/factor
370 vcko(i,k) = ((1.-tem)*vcko(i,k-1)+ptem*v1(i,k)
371 & +ptem1*v1(i,k-1))/factor
378 if (cnvflg(i) .and. k <= kpbl(i))
then
379 dz = zl(i,k) - zl(i,k-1)
380 tem = 0.5 * xlamue(i,k-1) * dz
383 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
384 & (q1(i,k,n)+q1(i,k-1,n)))/factor
real(kind=kind_phys), parameter con_g
gravity ( )
subroutine mfpbl(im, ix, km, ntrac, delt, cnvflg, zl, zm, thvx, q1, t1, u1, v1, hpbl, kpbl, sflx, ustar, wstar, xmf, tcko, qcko, ucko, vcko)
This subroutine is used for calculating the mass flux and updraft properties.
real(kind=kind_phys), parameter con_cp
spec heat air at p ( )