!WRF:model_layer:physics ! module module_bl_shinhong ! contains ! !------------------------------------------------------------------------------- ! subroutine shinhong(u3d,v3d,th3d,t3d,qv3d,qc3d,qi3d,p3d,p3di,pi3d, & rublten,rvblten,rthblten, & rqvblten,rqcblten,rqiblten,flag_qi, & cp,g,rovcp,rd,rovg,ep1,ep2,karman,xlv,rv, & dz8w,psfc, & znu,znw,p_top, & znt,ust,hpbl,psim,psih, & xland,hfx,qfx,wspd,br, & dt,kpbl2d, & exch_h, & u10,v10, & shinhong_tke_diag,tke_pbl,el_pbl,corf, & dx,dy, & ids,ide, jds,jde, kds,kde, & ims,ime, jms,jme, kms,kme, & its,ite, jts,jte, kts,kte, & !optional ctopo,ctopo2, & wstar,delta, & regime ) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- !-- u3d 3d u-velocity interpolated to theta points (m/s) !-- v3d 3d v-velocity interpolated to theta points (m/s) !-- th3d 3d potential temperature (k) !-- t3d temperature (k) !-- qv3d 3d water vapor mixing ratio (kg/kg) !-- qc3d 3d cloud mixing ratio (kg/kg) !-- qi3d 3d ice mixing ratio (kg/kg) ! (note: if P_QI0..and.dth(k-1)<=0.)then dth(k)=dth(k)+ct exit endif enddo ! ! compute local gradient richardson number ! do k = kte,kts+1,-1 rdz=2./(z(k+1)-z(k-1)) s2l=((u(k)-u(k-1))**2+(v(k)-v(k-1))**2)*rdz*rdz ! s**2 if(pblflg.and.k.le.lpbl)then suk=(u(k)-u(k-1))*rdz svk=(v(k)-v(k-1))*rdz s2l=(suk-hgamu/hpbl-ufxpbl(k))*suk+(svk-hgamv/hpbl-vfxpbl(k))*svk endif s2l=max(s2l,epsgm) s2(k)=s2l ! tem=(t(k)+t(k-1))*0.5 thm=(the(k)+the(k-1))*0.5 a=thm*p608 b=(elocp/tem-1.-p608)*thm ghl=(dth(k)*((q(k)+q(k-1)+cwm(k)+cwm(k-1))*(0.5*p608)+1.) & +(q(k)-q(k-1)+cwm(k)-cwm(k-1))*a & +(cwm(k)-cwm(k-1))*b)*rdz ! dtheta/dz if(pblflg.and.k.le.lpbl)then ghl=ghl-mf(k)-(hgamq/hpbl+qfxpbl(k))*a endif if(abs(ghl)<=epsgh)ghl=epsgh ! en2(k)=ghl*g/thm ! n**2 gh(k)=ghl ri(k)=en2(k)/s2l enddo ! ! find maximum mixing lengths and the level of the pbl top ! do k = kte,kts+1,-1 s2l=s2(k) ghl=gh(k) if(ghl>=epsgh)then if(s2l/ghl<=requ)then elm(k)=epsl else aubr=(aubm*s2l+aubh*ghl)*ghl bubr= bubm*s2l+bubh*ghl qol2st=(-0.5*bubr+sqrt(bubr*bubr*0.25-aubr*cubr))*rcubr eloq2x=1./qol2st elm(k)=max(sqrt(eloq2x*q2(k)),epsl) endif else aden=(adnm*s2l+adnh*ghl)*ghl bden= bdnm*s2l+bdnh*ghl qol2un=-0.5*bden+sqrt(bden*bden*0.25-aden) eloq2x=1./(qol2un+epsru) ! repsr1/qol2un elm(k)=max(sqrt(eloq2x*q2(k)),epsl) endif enddo ! do k = lpbl,lmh,-1 q1(k)=sqrt(q2(k)) enddo ! szq=0. sq =0. do k = kte,kts+1,-1 qdzl=(q1(k)+q1(k-1))*(z(k)-z(k-1)) szq=(z(k)+z(k-1)-z(lmh)-z(lmh))*qdzl+szq sq=qdzl+sq enddo ! ! computation of asymptotic l in blackadar formula ! el0=min(alph*szq*0.5/sq,el0max) el0=max(el0 ,el0min) ! ! above the pbl top ! lpblm=min(lpbl+1,kte) do k = kte,lpblm,-1 el(k)=(z(k+1)-z(k-1))*elfc rel(k)=el(k)/elm(k) enddo ! ! inside the pbl ! epshol=min(epshol,0.0) ckp=elcbl*((1.0-8.0*epshol)**(1./3.)) if(lpbl>lmh)then do k = lpbl,lmh+1,-1 vkrmz=(z(k)-z(lmh))*vkarman if(pblflg) then vkrmz=ckp*(z(k)-z(lmh))*vkarman el(k)=vkrmz/(vkrmz/el0+1.) else el(k)=vkrmz/(vkrmz/el0+1.) endif rel(k)=el(k)/elm(k) enddo endif ! do k = lpbl-1,lmh+2,-1 srel=min(((rel(k-1)+rel(k+1))*0.5+rel(k))*0.5,rel(k)) el(k)=max(srel*elm(k),epsl) enddo ! ! mixing length for the qnse model in stable case ! f=max(corf,eps1) rlambda=f/(blckdr*ustar) do k = kte,kts+1,-1 if(en2(k)>=0.0)then ! stable case vkrmz=(z(k)-z(lmh))*vkarman rlb=rlambda+1./vkrmz rln=sqrt(2.*en2(k)/q2(k))/cn el(k)=1./(rlb+rln) endif enddo ! end subroutine mixlen !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- subroutine prodq2(lmh,dtturbl,ustar,s2,ri,q2,el,z,akm,akh, & uxk,vxk,thxk,thvxk, & hgamu,hgamv,hgamq,delxy, & hpbl,pblflg,kpbl, & mf,ufxpbl,vfxpbl,qfxpbl, & p608, & ids,ide,jds,jde,kds,kde, & ims,ime,jms,jme,kms,kme, & its,ite,jts,jte,kts,kte) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- ! real,parameter :: epsq2l = 0.01,c0 = 0.55,ceps = 16.6,g = 9.81 ! integer, intent(in ) :: ids,ide, jds,jde, kds,kde, & ims,ime, jms,jme, kms,kme, & its,ite, jts,jte, kts,kte integer, intent(in ) :: lmh,kpbl ! real, intent(in ) :: p608,dtturbl,ustar real, intent(in ) :: hgamu,hgamv,hgamq,delxy,hpbl ! logical, intent(in ) :: pblflg ! real, dimension( kts:kte ) , & intent(in ) :: uxk, & vxk, & thxk, & thvxk real, dimension( kts+1:kte ) , & intent(in ) :: s2, & ri, & akm, & akh, & el, & mf, & ufxpbl, & vfxpbl, & qfxpbl ! real, dimension( kts:kte+1 ) , & intent(in ) :: z ! real, dimension( kts:kte ) , & intent(inout) :: q2 ! ! local vars ! integer :: k ! real :: s2l,q2l,deltaz,akml,akhl,en2,pr,bpr,dis,rc02 real :: suk,svk,gthvk,govrthvk,pru,prv real :: thm,disel ! !------------------------------------------------------------------------------- ! rc02=2.0/(c0*c0) ! ! start of production/dissipation loop ! main_integration: do k = kts+1,kte deltaz=0.5*(z(k+1)-z(k-1)) s2l=s2(k) q2l=q2(k) suk=(uxk(k)-uxk(k-1))/deltaz svk=(vxk(k)-vxk(k-1))/deltaz gthvk=(thvxk(k)-thvxk(k-1))/deltaz govrthvk=g/(0.5*(thvxk(k)+thvxk(k-1))) akml=akm(k) akhl=akh(k) en2=ri(k)*s2l !n**2 thm=(thxk(k)+thxk(k-1))*0.5 ! ! turbulence production term ! if(pblflg.and.k.le.kpbl)then pru=(akml*(suk-hgamu/hpbl-ufxpbl(k)))*suk prv=(akml*(svk-hgamv/hpbl-vfxpbl(k)))*svk else pru=akml*suk*suk prv=akml*svk*svk endif pr=pru+prv ! ! buoyancy production ! if(pblflg.and.k.le.kpbl)then bpr=(akhl*(gthvk-mf(k)-(hgamq/hpbl+qfxpbl(k))*p608*thm))*govrthvk else bpr=akhl*gthvk*govrthvk endif ! ! dissipation ! disel=min(delxy,ceps*el(k)) dis=(q2l)**1.5/disel ! q2l=q2l+2.0*(pr-bpr-dis)*dtturbl q2(k)=amax1(q2l,epsq2l) ! ! end of production/dissipation loop ! enddo main_integration ! ! lower boundary condition for q2 ! q2(kts)=amax1(rc02*ustar*ustar,epsq2l) ! end subroutine prodq2 !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- subroutine vdifq(lmh,dtdif,q2,el,z, & akhk,ptke1, & hgame,hpbl,pblflg,kpbl, & efxpbl, & ids,ide,jds,jde,kds,kde, & ims,ime,jms,jme,kms,kme, & its,ite,jts,jte,kts,kte) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- ! real,parameter :: c_k=1.0,esq=5.0 ! integer, intent(in ) :: ids,ide, jds,jde, kds,kde, & ims,ime, jms,jme, kms,kme, & its,ite, jts,jte, kts,kte integer, intent(in ) :: lmh,kpbl ! real, intent(in ) :: dtdif,hpbl,efxpbl ! logical, intent(in ) :: pblflg ! real, dimension( kts:kte ) , & intent(in ) :: hgame, & ptke1 real, dimension( kts+1:kte ) , & intent(in ) :: el, & akhk real, dimension( kts:kte+1 ) , & intent(in ) :: z ! real, dimension( kts:kte ) , & intent(inout) :: q2 ! ! local vars ! integer :: k ! real :: aden,akqs,bden,besh,besm,cden,cf,dtozs,ell,eloq2,eloq4 real :: elqdz,esh,esm,esqhf,ghl,gml,q1l,rden,rdz real :: zak ! real, dimension( kts+1:kte ) :: zfacentk real, dimension( kts+2:kte ) :: akq, & cm, & cr, & dtoz, & rsq2 ! !------------------------------------------------------------------------------- ! ! vertical turbulent diffusion ! esqhf=0.5*esq do k = kts+1,kte zak=0.5*(z(k)+z(k-1)) !zak of vdifq = za(k-1) of shinhong2d zfacentk(k)=(zak/hpbl)**3.0 enddo ! do k = kte,kts+2,-1 dtoz(k)=(dtdif+dtdif)/(z(k+1)-z(k-1)) akq(k)=c_k*(akhk(k)/(z(k+1)-z(k-1))+akhk(k-1)/(z(k)-z(k-2))) akq(k)=akq(k)*ptke1(k) cr(k)=-dtoz(k)*akq(k) enddo ! akqs=c_k*akhk(kts+1)/(z(kts+2)-z(kts)) akqs=akqs*ptke1(kts+1) cm(kte)=dtoz(kte)*akq(kte)+1. rsq2(kte)=q2(kte) ! do k = kte-1,kts+2,-1 cf=-dtoz(k)*akq(k+1)/cm(k+1) cm(k)=-cr(k+1)*cf+(akq(k+1)+akq(k))*dtoz(k)+1. rsq2(k)=-rsq2(k+1)*cf+q2(k) if(pblflg.and.k.lt.kpbl) then rsq2(k)=rsq2(k)-dtoz(k)*(2.0*hgame(k)/hpbl)*akq(k+1)*(z(k+1)-z(k)) & +dtoz(k)*(2.0*hgame(k-1)/hpbl)*akq(k)*(z(k)-z(k-1)) rsq2(k)=rsq2(k)-dtoz(k)*2.0*efxpbl*zfacentk(k+1) & +dtoz(k)*2.0*efxpbl*zfacentk(k) endif enddo ! dtozs=(dtdif+dtdif)/(z(kts+2)-z(kts)) cf=-dtozs*akq(lmh+2)/cm(lmh+2) ! if(pblflg.and.((lmh+1).lt.kpbl)) then q2(lmh+1)=(dtozs*akqs*q2(lmh)-rsq2(lmh+2)*cf+q2(lmh+1) & -dtozs*(2.0*hgame(lmh+1)/hpbl)*akq(lmh+2)*(z(lmh+2)-z(lmh+1)) & +dtozs*(2.0*hgame(lmh)/hpbl)*akqs*(z(lmh+1)-z(lmh))) q2(lmh+1)=q2(lmh+1)-dtozs*2.0*efxpbl*zfacentk(lmh+2) & +dtozs*2.0*efxpbl*zfacentk(lmh+1) q2(lmh+1)=q2(lmh+1)/((akq(lmh+2)+akqs)*dtozs-cr(lmh+2)*cf+1.) else q2(lmh+1)=(dtozs*akqs*q2(lmh)-rsq2(lmh+2)*cf+q2(lmh+1)) & /((akq(lmh+2)+akqs)*dtozs-cr(lmh+2)*cf+1.) endif ! do k = lmh+2,kte q2(k)=(-cr(k)*q2(k-1)+rsq2(k))/cm(k) enddo ! end subroutine vdifq !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- subroutine shinhonginit(rublten,rvblten,rthblten,rqvblten, & rqcblten,rqiblten, & tke_pbl, & p_qi,p_first_scalar, & restart, allowed_to_read, & ids, ide, jds, jde, kds, kde, & ims, ime, jms, jme, kms, kme, & its, ite, jts, jte, kts, kte ) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- ! real,parameter :: epsq2l = 0.01 logical , intent(in) :: restart, allowed_to_read integer , intent(in) :: ids, ide, jds, jde, kds, kde, & ims, ime, jms, jme, kms, kme, & its, ite, jts, jte, kts, kte integer , intent(in) :: p_qi,p_first_scalar real , dimension( ims:ime , kms:kme , jms:jme ), intent(out) :: & rublten, & rvblten, & rthblten, & rqvblten, & rqcblten, & rqiblten real , dimension( ims:ime , kms:kme , jms:jme ), intent(out) :: & tke_pbl integer :: i, j, k, itf, jtf, ktf ! jtf = min0(jte,jde-1) ktf = min0(kte,kde-1) itf = min0(ite,ide-1) ! if(.not.restart)then do j = jts,jtf do k = kts,ktf do i = its,itf rublten(i,k,j) = 0. rvblten(i,k,j) = 0. rthblten(i,k,j) = 0. rqvblten(i,k,j) = 0. rqcblten(i,k,j) = 0. tke_pbl(i,k,j) = epsq2l/2. enddo enddo enddo endif ! if (p_qi .ge. p_first_scalar .and. .not.restart) then do j = jts,jtf do k = kts,ktf do i = its,itf rqiblten(i,k,j) = 0. enddo enddo enddo endif ! end subroutine shinhonginit !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- function pu(d,h) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- real :: pu real,parameter :: pmin = 0.0,pmax = 1.0 real,parameter :: a1 = 1.0, a2 = 0.070, a3 = 1.0, a4 = 0.142, a5 = 0.071 real,parameter :: b1 = 2.0, b2 = 0.6666667 real :: d,h,doh,num,den ! doh=d/h num=a1*(doh)**b1+a2*(doh)**b2 den=a3*(doh)**b1+a4*(doh)**b2+a5 pu=num/den pu=max(pu,pmin) pu=min(pu,pmax) ! return end function !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- function pq(d,h) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- real :: pq real,parameter :: pmin = 0.0,pmax = 1.0 real,parameter :: a1 = 1.0, a2 = -0.098, a3 = 1.0, a4 = 0.106, a5 = 0.5 real,parameter :: b1 = 2.0 real :: d,h,doh,num,den ! doh=d/h num=a1*(doh)**b1+a2 den=a3*(doh)**b1+a4 pq=a5*num/den+(1.-a5) pq=max(pq,pmin) pq=min(pq,pmax) ! return end function !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- function pthnl(d,h) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- real :: pthnl real,parameter :: pmin = 0.0,pmax = 1.0 real,parameter :: a1 = 1.000, a2 = 0.936, a3 = -1.110, & a4 = 1.000, a5 = 0.312, a6 = 0.329, a7 = 0.243 real,parameter :: b1 = 2.0, b2 = 0.875 real :: d,h,doh,num,den ! doh=d/h num=a1*(doh)**b1+a2*(doh)**b2+a3 den=a4*(doh)**b1+a5*(doh)**b2+a6 pthnl=a7*num/den+(1.-a7) pthnl=max(pthnl,pmin) pthnl=min(pthnl,pmax) ! return end function !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- function pthl(d,h) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- real :: pthl real,parameter :: pmin = 0.0,pmax = 1.0 real,parameter :: a1 = 1.000, a2 = 0.870, a3 = -0.913, & a4 = 1.000, a5 = 0.153, a6 = 0.278, a7 = 0.280 real,parameter :: b1 = 2.0, b2 = 0.5 real :: d,h,doh,num,den ! doh=d/h num=a1*(doh)**b1+a2*(doh)**b2+a3 den=a4*(doh)**b1+a5*(doh)**b2+a6 pthl=a7*num/den+(1.-a7) pthl=max(pthl,pmin) pthl=min(pthl,pmax) ! return end function !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- function ptke(d,h) !------------------------------------------------------------------------------- implicit none !------------------------------------------------------------------------------- real :: ptke real,parameter :: pmin = 0.0,pmax = 1.0 real,parameter :: a1 = 1.000, a2 = 0.070, & a3 = 1.000, a4 = 0.142, a5 = 0.071 real,parameter :: b1 = 2.0, b2 = 0.6666667 real :: d,h,doh,num,den ! doh=d/h num=a1*(doh)**b1+a2*(doh)**b2 den=a3*(doh)**b1+a4*(doh)**b2+a5 ptke=num/den ptke=max(ptke,pmin) ptke=min(ptke,pmax) ! return end function !------------------------------------------------------------------------------- ! !------------------------------------------------------------------------------- end module module_bl_shinhong !-------------------------------------------------------------------------------