CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mfpbl.f
1
3
6 module mfpbl_mod
7 contains
48 subroutine mfpbl(im,ix,km,ntrac,delt,cnvflg, &
49 & zl,zm,thvx,q1,t1,u1,v1,hpbl,kpbl, &
50 & sflx,ustar,wstar,xmf,tcko,qcko,ucko,vcko)
51!
52 use machine , only : kind_phys
53 use physcons, grav => con_g, cp => con_cp
54!
55 implicit none
56!
57 integer im, ix, km, ntrac
58! &, me
59 integer kpbl(im)
60 logical cnvflg(im)
61 real(kind=kind_phys) delt
62 real(kind=kind_phys) q1(ix,km,ntrac), t1(ix,km), &
63 & u1(ix,km), v1(ix,km), &
64 & thvx(im,km), &
65 & zl(im,km), zm(im,km+1), &
66 & hpbl(im), sflx(im), ustar(im), &
67 & wstar(im), xmf(im,km), &
68 & tcko(im,km),qcko(im,km,ntrac), &
69 & ucko(im,km),vcko(im,km)
70!
71c local variables and arrays
72!
73 integer i, j, k, n, kmpbl
74!
75 real(kind=kind_phys) dt2, dz, ce0,
76 & h1, factor, gocp,
77 & g, c1, d1,
78 & b1, f1, bb1, bb2,
79 & alp, a1, qmin, zfmin,
80 & xmmx, rbint, tau,
81! & rbint, tau,
82 & tem, tem1, tem2,
83 & ptem, ptem1, ptem2,
84 & pgcon
85!
86 real(kind=kind_phys) sigw1(im), usws3(im), xlamax(im),
87 & rbdn(im), rbup(im), delz(im)
88!
89 real(kind=kind_phys) wu2(im,km), xlamue(im,km),
90 & thvu(im,km), zi(im,km),
91 & buo(im,km)
92!
93 logical totflg, flg(im)
94!
95c physical parameters
96 parameter(g=grav)
97 parameter(gocp=g/cp)
98! parameter(ce0=0.37,qmin=1.e-8,alp=1.0,pgcon=0.55)
99 parameter(ce0=0.38,qmin=1.e-8,alp=1.0,pgcon=0.55)
100 parameter(a1=0.08,b1=0.5,f1=0.15,c1=0.3,d1=2.58,tau=500.)
101 parameter(zfmin=1.e-8,h1=0.33333333)
102!
103c-----------------------------------------------------------------------
104!
105!************************************************************************
106!
107 kmpbl = km/2 + 1
108 dt2 = delt
110 totflg = .true.
111 do i=1,im
112 totflg = totflg .and. (.not. cnvflg(i))
113 enddo
114 if(totflg) return
115!!
116 do k = 1, km
117 do i=1,im
118 if (cnvflg(i)) then
119 zi(i,k) = zm(i,k+1)
120 endif
121 enddo
122 enddo
125 do i=1,im
126 if(cnvflg(i)) then
127 k = kpbl(i) / 2
128 k = max(k, 1)
129 delz(i) = zl(i,k+1) - zl(i,k)
130 xlamax(i) = ce0 / delz(i)
131 endif
132 enddo
133 do k = 1, kmpbl
134 do i=1,im
135 if(cnvflg(i)) then
136 if(k < kpbl(i)) then
137 ptem = 1./(zi(i,k)+delz(i))
138 tem = max((hpbl(i)-zi(i,k)+delz(i)) ,delz(i))
139 ptem1 = 1./tem
140 xlamue(i,k) = ce0 * (ptem+ptem1)
141 else
142 xlamue(i,k) = xlamax(i)
143 endif
144 endif
145 enddo
146 enddo
147c
148c compute thermal excess
149c
151 do i=1,im
152 if(cnvflg(i)) then
153 tem = zl(i,1)/hpbl(i)
154 usws3(i) = (ustar(i)/wstar(i))**3.
155 tem1 = usws3(i) + 0.6*tem
156 tem2 = max((1.-tem), zfmin)
157 ptem = (tem1**h1) * sqrt(tem2)
158 sigw1(i) = 1.3 * ptem * wstar(i)
159 ptem1 = alp * sflx(i) / sigw1(i)
160 thvu(i,1) = thvx(i,1) + ptem1
161 buo(i,1) = g * (thvu(i,1)/thvx(i,1)-1.)
162 endif
163 enddo
164c
165c compute potential temperature and buoyancy for updraft air parcel
166c
172 do k = 2, kmpbl
173 do i=1,im
174 if(cnvflg(i)) then
175 dz = zl(i,k) - zl(i,k-1)
176 tem = xlamue(i,k-1) * dz
177 ptem = 2. + tem
178 ptem1 = (2. - tem) / ptem
179 tem1 = tem * (thvx(i,k)+thvx(i,k-1)) / ptem
180 thvu(i,k) = ptem1 * thvu(i,k-1) + tem1
181 buo(i,k) = g * (thvu(i,k)/thvx(i,k)-1.)
182 endif
183 enddo
184 enddo
185c
186c compute updraft velocity square(wu2)
187c
193! tem = 1.-2.*f1
194! bb1 = 2. * b1 / tem
195! bb2 = 2. / tem
196! from soares et al. (2004,qjrms)
197! bb1 = 2.
198! bb2 = 4.
199!
200! from bretherton et al. (2004, mwr)
201! bb1 = 4.
202! bb2 = 2.
203!
204! from our tuning
205 bb1 = 1.8
206 bb2 = 3.5
207!
208 do i = 1, im
209 if(cnvflg(i)) then
210!
211! tem = zi(i,1)/hpbl(i)
212! tem1 = usws3(i) + 0.6*tem
213! tem2 = max((1.-tem), zfmin)
214! ptem = (tem1**h1) * sqrt(tem2)
215! ptem1 = 1.3 * ptem * wstar(i)
216! wu2(i,1) = d1*d1*ptem1*ptem1
217!
218 dz = zi(i,1)
219 tem = 0.5*bb1*xlamue(i,1)*dz
220 tem1 = bb2 * buo(i,1) * dz
221 ptem1 = 1. + tem
222 wu2(i,1) = tem1 / ptem1
223!
224 endif
225 enddo
226 do k = 2, kmpbl
227 do i = 1, im
228 if(cnvflg(i)) then
229 dz = zi(i,k) - zi(i,k-1)
230 tem = 0.25*bb1*(xlamue(i,k)+xlamue(i,k-1))*dz
231 tem1 = bb2 * buo(i,k) * dz
232 ptem = (1. - tem) * wu2(i,k-1)
233 ptem1 = 1. + tem
234 wu2(i,k) = (ptem + tem1) / ptem1
235 endif
236 enddo
237 enddo
238c
239c update pbl height as the height where updraft velocity vanishes
240c
243 do i=1,im
244 flg(i) = .true.
245 if(cnvflg(i)) then
246 flg(i) = .false.
247 rbup(i) = wu2(i,1)
248 endif
249 enddo
250 do k = 2, kmpbl
251 do i = 1, im
252 if(.not.flg(i)) then
253 rbdn(i) = rbup(i)
254 rbup(i) = wu2(i,k)
255 kpbl(i) = k
256 flg(i) = rbup(i).le.0.
257 endif
258 enddo
259 enddo
260 do i = 1,im
261 if(cnvflg(i)) then
262 k = kpbl(i)
263 if(rbdn(i) <= 0.) then
264 rbint = 0.
265 elseif(rbup(i) >= 0.) then
266 rbint = 1.
267 else
268 rbint = rbdn(i)/(rbdn(i)-rbup(i))
269 endif
270 hpbl(i) = zi(i,k-1) + rbint*(zi(i,k)-zi(i,k-1))
271 endif
272 enddo
273c
275 do i=1,im
276 if(cnvflg(i)) then
277 k = kpbl(i) / 2
278 k = max(k, 1)
279 delz(i) = zl(i,k+1) - zl(i,k)
280 xlamax(i) = ce0 / delz(i)
281 endif
282 enddo
283c
284c update entrainment rate
285c
286! do k = 1, kmpbl
287! do i=1,im
288! if(cnvflg(i)) then
289! if(k < kpbl(i)) then
290! tem = tau * sqrt(wu2(i,k))
291! tem1 = 1. / tem
292! ptem = ce0 / zi(i,k)
293! xlamue(i,k) = max(tem1, ptem)
294! else
295! xlamue(i,k) = xlamax(i)
296! endif
297! endif
298! enddo
299! enddo
300!
301 do k = 1, kmpbl
302 do i=1,im
303 if(cnvflg(i)) then
304 if(k < kpbl(i)) then
305 ptem = 1./(zi(i,k)+delz(i))
306 tem = max((hpbl(i)-zi(i,k)+delz(i)) ,delz(i))
307 ptem1 = 1./tem
308 xlamue(i,k) = ce0 * (ptem+ptem1)
309 else
310 xlamue(i,k) = xlamax(i)
311 endif
312 endif
313 enddo
314 enddo
315c
316c updraft mass flux as a function of sigmaw
317c(0.3*sigmaw[square root of vertical turbulence variance])
318c
320! do k = 1, kmpbl
321! do i=1,im
322! if(cnvflg(i) .and. k < kpbl(i)) then
323! tem = zi(i,k)/hpbl(i)
324! tem1 = usws3(i) + 0.6*tem
325! tem2 = max((1.-tem), zfmin)
326! ptem = (tem1**h1) * sqrt(tem2)
327! ptem1 = 1.3 * ptem * wstar(i)
328! xmf(i,k) = c1 * ptem1
329! endif
330! enddo
331! enddo
332c
333c updraft mass flux as a function of updraft velocity profile
334c
340 do k = 1, kmpbl
341 do i = 1, im
342 if (cnvflg(i) .and. k < kpbl(i)) then
343 xmf(i,k) = a1 * sqrt(wu2(i,k))
344 dz = zl(i,k+1) - zl(i,k)
345 xmmx = dz / dt2
346 xmf(i,k) = min(xmf(i,k),xmmx)
347 endif
348 enddo
349 enddo
350c
351c!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
352c compute updraft property
353c
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 ptem = tem + pgcon
374 ptem1= tem - pgcon
375!
376 tcko(i,k) = ((1.-tem)*tcko(i,k-1)+tem*
377 & (t1(i,k)+t1(i,k-1))-gocp*dz)/factor
378 ucko(i,k) = ((1.-tem)*ucko(i,k-1)+ptem*u1(i,k)
379 & +ptem1*u1(i,k-1))/factor
380 vcko(i,k) = ((1.-tem)*vcko(i,k-1)+ptem*v1(i,k)
381 & +ptem1*v1(i,k-1))/factor
382 endif
383 enddo
384 enddo
385 do n = 1, ntrac
386 do k = 2, kmpbl
387 do i = 1, im
388 if (cnvflg(i) .and. k <= kpbl(i)) then
389 dz = zl(i,k) - zl(i,k-1)
390 tem = 0.5 * xlamue(i,k-1) * dz
391 factor = 1. + tem
392
393 qcko(i,k,n) = ((1.-tem)*qcko(i,k-1,n)+tem*
394 & (q1(i,k,n)+q1(i,k-1,n)))/factor
395 endif
396 enddo
397 enddo
398 enddo
399!
400 return
401 end
402
403 end module mfpbl_mod
This module contains the subroutine that calculates the updraft properties and mass flux for use in t...
Definition mfpbl.f:6
This module contains the entity of GFS Noah LSM Model(Version 2.7).
Definition sflx.f:5