195 & IM,KM,A,B,C,U1,V1,T1,Q1,KPBL, &
196 & PRSI,DEL,PRSL,PRSLK,PHII, PHIL,DELTIM,KDT, &
197 & HPRIME,OC,OA4,CLX4,THETA,SIGMA,GAMMA,ELVMAX, &
198 & DUSFC,DVSFC,dtaux2d_ms,dtauy2d_ms,dtaux2d_bl, &
199 & dtauy2d_bl,dusfc_ms,dvsfc_ms,dusfc_bl,dvsfc_bl, &
200 & G, CP, RD, RV, IMX, &
201 & nmtvr, cdmbgwd, me, lprnt, ipr, rdxzb, ldiag_ugwp, &
296 USE machine ,
ONLY : kind_phys
300 integer,
intent(in) :: im, km, imx, kdt, ipr, me
301 integer,
intent(in) :: KPBL(:)
302 real(kind=kind_phys),
intent(in) :: &
303 & deltim, g, cp, rd, rv, cdmbgwd(:)
304 real(kind=kind_phys),
intent(inout) :: &
305 & a(:,:), b(:,:), c(:,:)
306 real(kind=kind_phys),
intent(in) :: &
307 & u1(:,:), v1(:,:), t1(:,:), &
308 & q1(:,:), prsi(:,:), del(:,:), &
309 & prsl(:,:), prslk(:,:), phil(:,:), &
311 real(kind=kind_phys),
intent(in) :: &
312 & oc(:), oa4(:,:), clx4(:,:), hprime(:)
313 real(kind=kind_phys),
intent(inout) :: elvmax(:)
314 real(kind=kind_phys),
intent(in) :: &
315 & theta(:), sigma(:), gamma(:)
316 real(kind=kind_phys),
intent(inout) :: dusfc(:), dvsfc(:), &
318 real(kind=kind_phys),
intent(inout),
optional :: dtaux2d_ms(:,:), &
319 & dtauy2d_ms(:,:), dtaux2d_bl(:,:), &
321 real(kind=kind_phys),
intent(inout),
optional :: dusfc_ms(:), &
322 & dvsfc_ms(:), dusfc_bl(:), dvsfc_bl(:)
323 integer,
intent(in) :: nmtvr
324 logical,
intent(in) :: lprnt
325 logical,
intent(in) :: ldiag_ugwp
326 character(len=*),
intent(out) :: errmsg
327 integer,
intent(out) :: errflg
331 real(kind=kind_phys) wk(im)
332 real(kind=kind_phys) bnv2lm(im,km),pe(im),ek(im),zbk(im),up(im)
333 real(kind=kind_phys) db(im,km),ang(im,km),uds(im,km)
334 real(kind=kind_phys) zlen, rtrm, phiang, cdmb, dbim, zr, cdmbo4
335 real(kind=kind_phys) eng0, eng1
339 real(kind=kind_phys) pi, dw2min, rimin, ric, bnv2min, efmin
340 &, efmax,hpmax,hpmin, rad_to_deg, deg_to_rad
341 parameter(pi=3.1415926535897931)
342 parameter(rad_to_deg=180.0/pi, deg_to_rad=pi/180.0)
343 parameter(dw2min=1., rimin=-100., ric=0.25, bnv2min=1.0e-5)
345 parameter(efmin=0.0, efmax=10.0, hpmax=2400.0, hpmin=1.0)
347 real(kind=kind_phys) frc, ce, ceofrc, frmax, cg, gmax
348 &, veleps, factop, rlolev, rdi
350 parameter(frc=1.0, ce=0.8, ceofrc=ce/frc, frmax=100., cg=0.5)
351 parameter(gmax=1.0, veleps=1.0, factop=0.5)
353 parameter(rlolev=50000.0)
357 real(kind=kind_phys) dpmin,hminmt,hncrit,minwnd,sigfac
360 parameter(hncrit=8000.)
362 parameter(sigfac=4.0)
363 parameter(hminmt=50.)
364 parameter(minwnd=0.1)
370 parameter(dpmin=5000.0)
373 real(kind=kind_phys) fdir
375 parameter(mdir=8, fdir=mdir/(pi+pi))
377 data nwdir/6,7,5,8,2,3,1,4/
384 real(kind=kind_phys) taub(im), xn(im), yn(im), ubar(im) &
385 &, vbar(im), ulow(im), oa(im), clx(im) &
386 &, roll(im), uloi(im) &
387 &, dtfac(im), xlinv(im), delks(im)
390 real(kind=kind_phys) bnv2(im,km), taup(im,km+1), ri_n(im,km) &
391 &, taud(im,km), ro(im,km), vtk(im,km) &
392 &, vtj(im,km), scor(im), velco(im,km-1) &
393 &, bnv2bar(im), cdsigohp(im)
396 integer kref(IM), kint(im), iwk(im), ipt(im)
400 integer idxzb(im), ktrial, klevm1
402 real(kind=kind_phys) gor, gocp, fv, gr2, bnv, fr &
403 &, brvf, cleff, tem, tem1, tem2, temc, temv&
404 &, wdir, ti, rdz, dw2, shr2, bvf2 &
405 &, rdelks, efact, coefm, gfobnv, onebg &
406 &, scork, rscor, hd, fro, rim, sira &
407 &, dtaux, dtauy, pkp1log, pklog &
408 &, cosang, sinang, cos2a, sin2a, oneocpdt
410 integer kmm1, kmm2, lcap, lcapp1, kbps, kbpsp1,kbpsm1 &
411 &, kmps, idir, nwd, i, j, k, klcap, kp1, kmpbl, npt, npr, kmll
422 cdmb = 4.0 * 192.0/float(imx)
423 if (cdmbgwd(1) >= 0.0) cdmb = cdmb * cdmbgwd(1)
442 oneocpdt = 1.0 / (cp*deltim)
452 IF ( nmtvr == 14)
then
457 IF (elvmax(i) > hminmt .and. hprime(i) > hpmin)
then
497 elvmax(j) = min(elvmax(j) + sigfac * hprime(j), hncrit)
498 cdsigohp(i) = cdmbo4 * sigma(j) / hprime(j)
506 pkp1log = phil(j,k+1) * onebg
507 pklog = phil(j,k) * onebg
510 if (elvmax(j) <= pkp1log .and. elvmax(j) >= pklog)
THEN
514 wk(i) = g * elvmax(j) / (phil(j,k+1) - phil(j,k))
515 iwklm(i) = max(iwklm(i), k+1)
521 vtj(i,k) = t1(j,k) * (1.0+fv*q1(j,k))
522 vtk(i,k) = vtj(i,k) / prslk(j,k)
523 ro(i,k) = rdi * prsl(j,k) / vtj(i,k)
546 rdz = g / (phil(j,kp1) - phil(j,k))
549 bnv2lm(i,k) = (g+g) * rdz * (vtk(i,kp1) - vtk(i,k))
550 & / (vtk(i,kp1) + vtk(i,k))
551 bnv2lm(i,k) = max( bnv2lm(i,k), bnv2min )
558 delks(i) = 1.0 / (prsi(j,1) - prsi(j,iwklm(i)))
565 bnv2bar(i) = (prsi(j,1)-prsl(j,1)) * delks(i) * bnv2lm(i,1)
590 rdelks = del(j,k) * delks(i)
591 ubar(i) = ubar(i) + rdelks * u1(j,k)
592 vbar(i) = vbar(i) + rdelks * v1(j,k)
593 roll(i) = roll(i) + rdelks * ro(i,k)
594 if (k < iwklm(i)-1)
then
595 rdelks = (prsl(j,k)-prsl(j,k+1)) * delks(i)
597 rdelks = (prsl(j,k)-prsi(j,k+1)) * delks(i)
599 bnv2bar(i) = bnv2bar(i) + bnv2lm(i,k) * rdelks
611 DO k = iwklm(i), 1, -1
612 phiang = atan2(v1(j,k),u1(j,k))*rad_to_deg
613 ang(i,k) = ( theta(j) - phiang )
614 if ( ang(i,k) > 90. ) ang(i,k) = ang(i,k) - 180.
615 if ( ang(i,k) < -90. ) ang(i,k) = ang(i,k) + 180.
616 ang(i,k) = ang(i,k) * deg_to_rad
626 & max(sqrt(u1(j,k)*u1(j,k) + v1(j,k)*v1(j,k)), minwnd)
628 IF (idxzb(i) == 0 )
then
629 pe(i) = pe(i) + bnv2lm(i,k) * (g*elvmax(j) - phil(j,k))
630 & * (phii(j,k+1) - phii(j,k))
636 up(i) = uds(i,k) * cos(ang(i,k))
637 ek(i) = 0.5 * up(i) * up(i)
640 IF (pe(i) >= ek(i))
THEN
642 rdxzb(j) = real(k,kind=kind_phys)
676 & - sqrt(ubar(i)*ubar(i) + vbar(i)*vbar(i))/bnv2bar(i)
691 IF (idxzb(i) > 0)
then
692 DO k = idxzb(i), 1, -1
693 IF (phil(j,idxzb(i)) > phil(j,k))
then
702 zlen = sqrt( (phil(j,idxzb(i)) - phil(j,k))
703 & / (phil(j,k ) + g*hprime(j)) )
715 cosang = cos(ang(i,k))
716 sinang = sin(ang(i,k))
717 cos2a = cosang * cosang
718 sin2a = sinang * sinang
719 tem = cos2a + gamma(j)*sin2a
722 if (abs(tem) > 1.e-06)
then
723 rtrm = (gamma(j)*cos2a + sin2a) / tem
724 elseif (tem > 0.0)
then
725 rtrm = (gamma(j)*cos2a + sin2a) * 1.0e6
727 rtrm = - (gamma(j)*cos2a + sin2a) * 1.0e6
729 zr = max( 2.0 - rtrm, 0. )
742 db(i,k) = cdsigohp(i) * zr * ro(i,k) * zlen
743 & * max(cosang, gamma(j)*sinang) * uds(i,k)
761 ELSEIF ( nmtvr /= 14)
then
766 IF ( hprime(i) > hpmin )
then
769 if (ipr == i) npr = npt
796 cleff = 0.5e-5 / sqrt(float(imx)/192.0)
802 if (cdmbgwd(2) >= 0.0) cleff = cleff * cdmbgwd(2)
807 vtj(i,k) = t1(j,k) * (1.0+fv*q1(j,k))
808 vtk(i,k) = vtj(i,k) / prslk(j,k)
809 ro(i,k) = rdi * prsl(j,k) / vtj(i,k)
817 ti = 2.0 / (t1(j,k)+t1(j,kp1))
818 tem = ti / (prsl(j,k)-prsl(j,kp1))
819 rdz = g / (phil(j,kp1) - phil(j,k))
820 tem1 = u1(j,k) - u1(j,kp1)
821 tem2 = v1(j,k) - v1(j,kp1)
822 dw2 = tem1*tem1 + tem2*tem2
823 shr2 = max(dw2,dw2min) * rdz * rdz
824 bvf2 = g*(gocp+rdz*(vtj(i,kp1)-vtj(i,k))) * ti
825 ri_n(i,k) = max(bvf2/shr2,rimin)
829 bnv2(i,k) = (g+g) * rdz * (vtk(i,kp1)-vtk(i,k))
830 & / (vtk(i,kp1)+vtk(i,k))
831 bnv2(i,k) = max( bnv2(i,k), bnv2min )
857 tem = (prsi(j,1) - prsi(j,k))
858 if (tem < dpmin) iwk(i) = k
868 kref(i) = max(iwk(i), kpbl(j)+1 )
869 delks(i) = 1.0 / (prsi(j,1) - prsi(j,kref(i)))
874 kbps = max(kbps, kref(i))
875 kmps = min(kmps, kref(i))
877 bnv2bar(i) = (prsi(j,1)-prsl(j,1)) * delks(i) * bnv2(i,1)
884 IF (k < kref(i))
THEN
886 rdelks = del(j,k) * delks(i)
887 ubar(i) = ubar(i) + rdelks * u1(j,k)
888 vbar(i) = vbar(i) + rdelks * v1(j,k)
890 roll(i) = roll(i) + rdelks * ro(i,k)
891 if (k < kref(i)-1)
then
892 rdelks = (prsl(j,k)-prsl(j,k+1)) * delks(i)
894 rdelks = (prsl(j,k)-prsi(j,k+1)) * delks(i)
896 bnv2bar(i) = bnv2bar(i) + bnv2(i,k) * rdelks
911 wdir = atan2(ubar(i),vbar(i)) + pi
912 idir = mod(nint(fdir*wdir),mdir) + 1
914 oa(i) = (1-2*int( (nwd-1)/4 )) * oa4(j,mod(nwd-1,4)+1)
915 clx(i) = clx4(j,mod(nwd-1,4)+1)
941 ulow(i) = max(sqrt(ubar(i)*ubar(i) + vbar(i)*vbar(i)), 1.0)
942 uloi(i) = 1.0 / ulow(i)
949 velco(i,k) = 0.5 * ((u1(j,k)+u1(j,kp1))*ubar(i)
950 & + (v1(j,k)+v1(j,kp1))*vbar(i))
951 velco(i,k) = velco(i,k) * uloi(i)
984 bnv = sqrt( bnv2bar(i) )
985 fr = bnv * uloi(i) * min(hprime(j),hpmax)
987 xn(i) = ubar(i) * uloi(i)
988 yn(i) = vbar(i) * uloi(i)
1033 efact = (oa(i) + 2.) ** (ceofrc*fr)
1034 efact = min( max(efact,efmin), efmax )
1036 coefm = (1. + clx(i)) ** (oa(i)+1.)
1038 xlinv(i) = coefm * cleff
1040 tem = fr * fr * oc(j)
1041 gfobnv = gmax * tem / ((tem + cg)*bnv)
1043 taub(i) = xlinv(i) * roll(i) * ulow(i) * ulow(i)
1044 & * ulow(i) * gfobnv * efact
1049 k = max(1, kref(i)-1)
1050 tem = max(velco(i,k)*velco(i,k), 0.1)
1051 scor(i) = bnv2(i,k) / tem
1059 IF (k <= kref(i)) taup(i,k) = taub(i)
1073 IF (k >= kref(i))
THEN
1074 icrilv(i) = icrilv(i) .OR. ( ri_n(i,k) < ric)
1075 & .OR. (velco(i,k) <= 0.0)
1092 IF (k >= kref(i))
THEN
1093 IF (.NOT.icrilv(i) .AND. taup(i,k) > 0.0 )
THEN
1094 temv = 1.0 / max(velco(i,k), 0.01)
1096 IF (oa(i) > 0. .AND. kp1 < kint(i))
THEN
1097 scork = bnv2(i,k) * temv * temv
1098 rscor = min(1.0, scork / scor(i))
1116 brvf = sqrt(bnv2(i,k))
1118 tem1 = xlinv(i)*(ro(i,kp1)+ro(i,k))*brvf*0.5
1119 & * max(velco(i,k),0.01)
1120 hd = sqrt(taup(i,k) / tem1)
1121 fro = brvf * hd * temv
1134 tem2 = sqrt(ri_n(i,k))
1135 tem = 1. + tem2 * fro
1136 rim = ri_n(i,k) * (1.-fro) / (tem * tem)
1157 IF (rim <= ric .AND.
1159 & (oa(i) <= 0. .OR. kp1 >= kint(i) ))
THEN
1160 temc = 2.0 + 1.0 / tem2
1161 hd = velco(i,k) * (2.*sqrt(temc)-temc) / brvf
1162 taup(i,kp1) = tem1 * hd * hd
1164 taup(i,kp1) = taup(i,k) * rscor
1166 taup(i,kp1) = min(taup(i,kp1), taup(i,k))
1177 DO klcap = lcapp1, km+1
1179 sira = prsi(ipt(i),klcap) / prsi(ipt(i),lcap)
1180 taup(i,klcap) = sira * taup(i,lcap)
1189 taud(i,k) = g * (taup(i,k+1) - taup(i,k)) / del(ipt(i),k)
1198 taud(i,klcap) = taud(i,klcap) * factop
1208 IF (k > kref(i) .and. prsi(ipt(i),k) >= rlolev)
THEN
1209 IF(taud(i,k) /= 0.)
THEN
1210 tem = deltim * taud(i,k)
1211 dtfac(i) = min(dtfac(i),abs(velco(i,k)/tem))
1230 taud(i,k) = taud(i,k) * dtfac(i)
1231 dtaux = taud(i,k) * xn(i)
1232 dtauy = taud(i,k) * yn(i)
1233 eng0 = 0.5*(u1(j,k)*u1(j,k)+v1(j,k)*v1(j,k))
1235 if (k < idxzb(i))
then
1237 dbim = db(i,k) / (1.+db(i,k)*deltim)
1238 a(j,k) = - dbim * v1(j,k) + a(j,k)
1239 b(j,k) = - dbim * u1(j,k) + b(j,k)
1240 eng1 = eng0*(1.0-dbim*deltim)*(1.0-dbim*deltim)
1246 tem1 = dbim * del(j,k)
1247 dusfc(j) = dusfc(j) + onebg * tem1 * u1(j,k)
1248 dvsfc(j) = dvsfc(j) + onebg * tem1 * v1(j,k)
1250 if (ldiag_ugwp)
then
1251 dtaux2d_bl(j,k) = - dbim * u1(j,k)
1252 dtauy2d_bl(j,k) = - dbim * v1(j,k)
1253 dusfc_bl(j) = dusfc_bl(j) + dtaux2d_bl(j,k) * del(j,k)
1254 dvsfc_bl(j) = dvsfc_bl(j) + dtauy2d_bl(j,k) * del(j,k)
1258 a(j,k) = dtauy + a(j,k)
1259 b(j,k) = dtaux + b(j,k)
1260 tem1 = u1(j,k) + dtaux*deltim
1261 tem2 = v1(j,k) + dtauy*deltim
1262 eng1 = 0.5 * (tem1*tem1+tem2*tem2)
1263 dusfc(j) = dusfc(j) - onebg * dtaux * del(j,k)
1264 dvsfc(j) = dvsfc(j) - onebg * dtauy * del(j,k)
1266 if (ldiag_ugwp)
then
1267 dtaux2d_ms(j,k) = dtaux
1268 dtauy2d_ms(j,k) = dtauy
1269 dusfc_ms(j) = dusfc_ms(j) + dtaux * del(j,k)
1270 dvsfc_ms(j) = dvsfc_ms(j) + dtauy * del(j,k)
1273 c(j,k) = c(j,k) + max(eng0-eng1,0.) * oneocpdt
1283 if (ldiag_ugwp)
then
1287 dusfc_ms(j) = - onebg * dusfc_ms(j)
1288 dvsfc_ms(j) = - onebg * dvsfc_ms(j)
1289 dusfc_bl(j) = - onebg * dusfc_bl(j)
1290 dvsfc_bl(j) = - onebg * dvsfc_bl(j)