!----------------------------------------------------------------------- ! module module_cu_bmj ! !----------------------------------------------------------------------- ! !*** the convection drivers and packages ! !----------------------------------------------------------------------- ! use module_include ! !----------------------------------------------------------------------- ! implicit none ! !----------------------------------------------------------------------- ! private ! public :: bmj_init,bmjdrv ! !----------------------------------------------------------------------- !----------------------------------------------------------------------- !*** for bmj convection !----------------------------------------------------------------------- !----------------------------------------------------------------------- ! real(kind=kfpt),parameter:: & dttop=0.,efifc=5.0,efimn=0.20 & ,efmntl=0.70,efmnts=0.70 & ,eliwv=2.683e6,enplo=98500.,enpup=95000. & ,epsdn=1.05,epsdt=0. & ,epsntp=.0001,epsntt=.0001,epspr=1.e-7 & ,epsup=1.00 & ,fup=1./200000. & ,pbm=13000.,pfrz=15000.,pno=1000. & ,pone=2500.,pqm=20000. & ,psh=20000.,pshu=45000. & ,rendp=1./(enplo-enpup) & ,rhlsc=0.00,rhhsc=1.00 & ,row=1.e3 & ,stabdf=0.90,stabds=0.90 & ,stabs=1.0,stresh=1.10 & ,dtshal=-1.0,trel=2400. ! real(kind=kfpt),parameter:: & dttrigr=-0.0 & ,dtptrigr=dttrigr*pone !<-- average parcel virtual temperature deficit over depth pone !<-- note: capetrigr is scaled by the cloud base temperature ! (see below) ! real(kind=kfpt),parameter:: & dspbfl_base=-3875. & ,dsp0fl_base=-5875. & ,dsptfl_base=-1875. & ,dspbfs_base=-3875. & ,dsp0fs_base=-5875. & ,dsptfs_base=-1875. ! real(kind=kfpt),parameter:: & pl=2500.,plq=70000.,ph=105000. & ,thl=210.,thh=365.,thhq=325. ! integer(kind=kint),parameter:: & itb=76,jtb=134,itbq=152,jtbq=440 ! integer(kind=kint),parameter:: & itrefi_max=3 ! !*** arrays for lookup tables ! real(kind=kfpt),dimension(itb),private,save:: & sthe,the0 real(kind=kfpt),dimension(jtb),private,save:: & qs0,sqs real(kind=kfpt),dimension(itbq),private,save:: & stheq,the0q real(kind=kfpt),dimension(itb,jtb),private,save:: & ptbl real(kind=kfpt),dimension(jtb,itb),private,save:: & ttbl real(kind=kfpt),dimension(jtbq,itbq),private,save:: & ttblq !*** share copies for module_bl_myjpbl ! real(kind=kfpt),dimension(jtb):: & qs0_exp,sqs_exp real(kind=kfpt),dimension(itb,jtb):: & ptbl_exp ! real(kind=kfpt),parameter:: & rdp=(itb-1.)/(ph-pl),rdpq=(itbq-1.)/(ph-plq) & ,rdq=itb-1,rdth=(jtb-1.)/(thh-thl) & ,rdthe=jtb-1.,rdtheq=jtbq-1. & ,rsfcp=1./101300. ! real(kind=kfpt),parameter:: & avgefi=(efimn+1.)*0.5 & ,stefi=1. ! !----------------------------------------------------------------------- ! contains ! !----------------------------------------------------------------------- !&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&& !----------------------------------------------------------------------- ! subroutine bmjdrv( & ids,ide,jds,jde & ,ims,ime,jms,jme & ,its,ite,jts,jte,lm & ,its_b1,ite_b1,jts_b1,jte_b1 & ,entrain,newall,newswap,newupup,nodeep & ,a2,a3,a4,cappa,cp,eliv,elwv,epsq,g & ,p608,pq0,r_d,tiw & ,fres,fr,fsl,fss & ,dt,dyh,ntsd,ncnvc & ,raincv,cutop,cubot,dxh,kpbl & ,th,t,q,u_phy,v_phy,dudt_phy,dvdt_phy & ,phint,phmid,exner & ,cldefi,xland,cu_act_flag & ! optional ,rthcuten,rqcuten & ) !----------------------------------------------------------------------- implicit none !----------------------------------------------------------------------- ! logical(kind=klog),intent(in):: & entrain,newall,newswap,newupup,nodeep ! logical(kind=klog),dimension(ims:ime,jms:jme),intent(inout):: & cu_act_flag ! integer(kind=kint),intent(in):: & ids,ide,jds,jde & ,ims,ime,jms,jme & ,its,ite,jts,jte,lm & ,its_b1,ite_b1,jts_b1,jte_b1 & ,ntsd,ncnvc ! integer(kind=kint),dimension(ims:ime,jms:jme),intent(in):: & kpbl ! real(kind=kfpt),intent(in):: & a2,a3,a4,cappa,cp,dt,dyh,eliv,elwv,epsq & ,fres,fr,fsl,fss,g,p608,pq0,r_d,tiw ! real(kind=kfpt),dimension(jds:jde),intent(in):: & dxh ! real(kind=kfpt),dimension(ims:ime,jms:jme),intent(in):: & xland ! real(kind=kfpt),dimension(ims:ime,jms:jme,1:lm),intent(in):: & q & ,exner,phmid,t,th,u_phy,v_phy ! real(kind=kfpt),dimension(ims:ime,jms:jme,1:lm+1),intent(in):: & phint ! real(kind=kfpt),dimension(ims:ime,jms:jme,1:lm) & ,optional & ,intent(inout):: & rqcuten,rthcuten ! real(kind=kfpt),dimension(ims:ime,jms:jme),intent(inout):: & cldefi,raincv ! real(kind=kfpt),dimension(ims:ime,jms:jme),intent(out):: & cubot,cutop real(kind=kfpt),dimension(ims:ime,jms:jme,1:lm),intent(out):: & dudt_phy,dvdt_phy ! !----------------------------------------------------------------------- !*** !*** local variables !*** !----------------------------------------------------------------------- integer(kind=kint):: & i,icldck,ierr,j,k,lbot,lmh,lpbl,ltop ! real(kind=kfpt):: & dspbfl,dsp0fl,dsptfl,dspbfs,dsp0fs,dsptfs & ,dspbsl,dsp0sl,dsptsl,dspbss,dsp0ss,dsptss & ,rdxeq,fresp & ,slopbl,slop0l,sloptl,slopbs,slop0s,slopts & ,dtcnvc,seamask,pcpcol,psfc,ptop,qnew,qold ! real(kind=kfpt),dimension(1:lm):: & dpcol,dqdt,dtdt,dudt,dvdt,pcol,qcol,tcol,exnercol,ucol,vcol ! !*** begin debugging convection real(kind=kfpt) :: delq,delt,plyr integer(kind=kint) :: imd,jmd logical(kind=klog) :: print_diag !*** end debugging convection ! !----------------------------------------------------------------------- !*********************************************************************** !----------------------------------------------------------------------- ! !*** prepare to call bmj convection scheme ! !----------------------------------------------------------------------- ! !*** check to see if this is a convection timestep ! icldck=mod(ntsd,ncnvc) !----------------------------------------------------------------------- ! !*** compute convection every ncnvc*dt/60.0 minutes ! !*** begin debugging convection imd=(ims+ime)/2 jmd=(jms+jme)/2 print_diag=.false. !*** end debugging convection if(icldck==0.or.ntsd==0)then ! do j=jts,jte do i=its,ite cu_act_flag(i,j)=.true. enddo enddo ! dtcnvc=dt*ncnvc rdxeq=1./dxh((jds+jde)/2+1) !....................................................................... !zj$omp parallel do & !zj$omp private (j,i,dqdt,dtdt,dudt,dvdt & !zj$omp ,dxh,fresp,pcpcol,psfc,ptop,seamask,k & !zj$omp ,tcol,pcol,dpcol,qcol,ucol,vcol & !zj$omp ,lmh,lpbl,delt,delq,plyr,lbot,ltop) !....................................................................... do j=jts_b1,jte_b1 fresp=fres !*(rdxeq*dxh(j))**0.125 do i=its_b1,ite_b1 !--set up deficit saturation pressures depending on resolution---------- dspbfl=dspbfl_base*fresp*fr dsp0fl=dsp0fl_base*fresp*fr dsptfl=dsptfl_base*fresp*fr dspbfs=dspbfs_base*fresp dsp0fs=dsp0fs_base*fresp dsptfs=dsptfs_base*fresp ! dspbsl=dspbfl*fsl dsp0sl=dsp0fl*fsl dsptsl=dsptfl*fsl dspbss=dspbfs*fss dsp0ss=dsp0fs*fss dsptss=dsptfs*fss ! slopbl=(dspbfl-dspbsl)/(1.-efimn) slop0l=(dsp0fl-dsp0sl)/(1.-efimn) sloptl=(dsptfl-dsptsl)/(1.-efimn) slopbs=(dspbfs-dspbss)/(1.-efimn) slop0s=(dsp0fs-dsp0ss)/(1.-efimn) slopts=(dsptfs-dsptss)/(1.-efimn) !------------------------------------------------------------------------ do k=1,lm dqdt(k)=0. dtdt(k)=0. dudt(k)=0. dvdt(k)=0. enddo ! raincv(i,j)=0. pcpcol=0. psfc=phint(i,j,lm+1) ptop=phint(i,j,1) ! !*** convert to bmj land mask (1.0 for sea; 0.0 for land) ! seamask=xland(i,j)-1. ! !*** fill 1-d vertical arrays ! do k=1,lm ! ucol (k)=u_phy(i,j,k) vcol (k)=v_phy(i,j,k) qcol (k)=max(epsq,q(i,j,k)) tcol (k)=t(i,j,k) pcol (k)=phmid(i,j,k) exnercol(k)=exner(i,j,k) dpcol (k)=phint(i,j,k+1)-phint(i,j,k) enddo ! !*** lowest layer above ground must also be flipped ! lmh=lm lpbl=kpbl(i,j) !----------------------------------------------------------------------- !*** !*** call convection !*** !----------------------------------------------------------------------- !*** begin debugging convection ! print_diag=.false. ! if(i==imd.and.j==jmd)print_diag=.true. !*** end debugging convection !----------------------------------------------------------------------- call bmj( & entrain,newall,newswap,newupup,nodeep,print_diag & ,i,j,lm,lmh,lpbl,ntsd & ,lbot,ltop & ,a2,a3,a4,cappa,cldefi(i,j),cp,dtcnvc,eliv,elwv & ,g,p608,pq0,r_d,tiw & ,psfc,ptop,seamask & ,dspbfl,dsp0fl,dsptfl,dspbfs,dsp0fs,dsptfs & ,dspbsl,dsp0sl,dsptsl,dspbss,dsp0ss,dsptss & ,slopbl,slop0l,sloptl,slopbs,slop0s,slopts & ,dpcol,pcol,qcol,tcol,exnercol,ucol,vcol & ,dqdt,dtdt,dudt,dvdt,pcpcol & ) !----------------------------------------------------------------------- !*** compute momentum tendencies ! do k=1,lm dudt_phy(i,j,k)=dudt(k) dvdt_phy(i,j,k)=dvdt(k) enddo ! !*** compute heating and moistening tendencies ! if(present(rthcuten).and.present(rqcuten))then do k=1,lm rthcuten(i,j,k)=dtdt(k)/exner(i,j,k) rqcuten(i,j,k)=dqdt(k) enddo endif ! !*** all units in bmj scheme are mks, thus convert precip from meters !*** to millimeters per step for output. ! raincv(i,j)=pcpcol*1.e3/ncnvc ! !*** convective cloud top and bottom from this call ! cutop(i,j)=float(lm+1-ltop) cubot(i,j)=float(lm+1-lbot) ! !----------------------------------------------------------------------- !*** begin debugging convection if(print_diag)then delt=0. delq=0. plyr=0. if(lbot>0.and.ltop cape) then ltp1=l cape=cpe(l) endif enddo !-- end do l=kb,1,-1 ! ltop=max(min(ltp1,lbot),1) ! !----------------------------------------------------------------------- !--------------- check for maximum instability ------------------------ !----------------------------------------------------------------------- if(cape > capecnv) then capecnv=cape pspcnv=psp thbtcnv=thbt lbotcnv=lbot ltopcnv=ltop do l=lmh,1,-1 cpecnv(l)=cpe(l) dtvcnv(l)=dtv(l) thescnv(l)=thes(l) enddo endif ! end if(cape > capecnv) then ! !----------------------------------------------------------------------- ! enddo max_buoy_loop ! !----------------------------------------------------------------------- !------------------------ maximum instability ------------------------ !----------------------------------------------------------------------- ! if(capecnv > 0.) then psp=pspcnv thbt=thbtcnv lbot=lbotcnv ltop=ltopcnv pbot=prsmid(lbot) ptop=prsmid(ltop) ! do l=lmh,1,-1 cpe(l)=cpecnv(l) dtv(l)=dtvcnv(l) thes(l)=thescnv(l) enddo ! endif ! !----------------------------------------------------------------------- !----- quick exit if cloud is too thin or no cape is present --------- !----------------------------------------------------------------------- ! if(ptop>pbot-pno.or.ltop>lbot-2.or.capecnv<=0.)then lbot=0 ltop=lm pbot=prsmid(lmh) ptop=pbot cldefi=avgefi*sm+stefi*(1.-sm) return endif ! !*** depth of cloud required to make the point a deep convection point !*** is a scaled value of psfc. ! depth=pbot-ptop ! !zj if(depth.ge.depmin*0.50) then !zj if(depth.ge.depmin*0.25) then !zj if(depth.ge.depmin*0.625) then if(depth.ge.depmin*0.75) then !zj if(depth.ge.depmin*0.85) then plume=.true. endif ! if(depth>=depmin) then deep=.true. else shallow=.true. go to 600 endif ! !----------------------------------------------------------------------- !dcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcd !dcdcdcdcdcdcdcdcdcdcdc deep convection dcdcdcdcdcdcdcdcdcdcdcdcdcd !dcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcd !----------------------------------------------------------------------- ! 300 continue ! lb =lbot efi=cldefi !----------------------------------------------------------------------- !--------------initialize variables in the convective column------------ !----------------------------------------------------------------------- !*** !*** one should note that the values assigned to the array trefk !*** in the following loop are really only relevant in anchoring the !*** reference temperature profile at level lb. when building the !*** reference profile from cloud base, then assigning the !*** ambient temperature to trefk is acceptable. however, when !*** building the reference profile from some other level (such as !*** one level above the ground), then trefk should be filled with !*** the temperatures in tref(l) which are the temperatures of !*** the moist adiabat through cloud base. by the time the line !*** numbered 450 has been reached, trefk actually does hold the !*** reference temperature profile. !*** do l=1,lmh dift(l)=0. difq(l)=0. tkl=t(l) tk(l)=tkl trefk(l)=tkl qkl=q(l) qk(l)=qkl qrefk(l)=qkl pkl=prsmid(l) pk(l)=pkl psk(l)=pkl rxnerp=rxner(l) rxnerk(l)=rxnerp ! !--- calculate temperature along moist adiabat (tref) ! if(pkl=depth)then if(lpqm)then psk(l)=pk(l)+dsp rxnersk(l)=(1.e5/psk(l))**cappa thsk(l)=trefk(l)*rxnerk(l) qrefk(l)=pq0/psk(l)*exp(a2*(thsk(l)-a3*rxnersk(l)) & /(thsk(l)-a4*rxnersk(l))) else qrefk(l)=qk(l) endif ! enddo !----------------------------------------------------------------------- !*** !*** enthalpy conservation integral !*** !----------------------------------------------------------------------- enthalpy_conservation : do iter=1,2 ! sumde=0. sumdp=0. ! do l=ltop,lb sumde=((tk(l)-trefk(l))*cp+(qk(l)-qrefk(l))*el(l))*dprs(l) & & +sumde sumdp=sumdp+dprs(l) enddo ! hcorr=sumde/(sumdp-dprs(ltop)) lcor=ltop+1 !*** !*** find lqm !*** lqm=1 do l=1,lb if(pk(l)<=pqm)lqm=l enddo !*** !*** above lqm correct temperature only !*** if(lcor<=lqm)then do l=lcor,lqm trefk(l)=trefk(l)+hcorr*rcp enddo lcor=lqm+1 endif !*** !*** below lqm correct both temperature and moisture !*** do l=lcor,lb tskl=trefk(l)*rxnerk(l)/rxnersk(l) dhdt=qrefk(l)*a23m4l/(tskl-a4)**2+cp trefk(l)=hcorr/dhdt+trefk(l) thskl=trefk(l)*rxnerk(l) qrefk(l)=pq0/psk(l)*exp(a2*(thskl-a3*rxnersk(l)) & /(thskl-a4*rxnersk(l))) enddo ! enddo enthalpy_conservation !----------------------------------------------------------------------- ! !*** heating, moistening, precipitation ! !----------------------------------------------------------------------- avrgt=0. preck=0. dsq=0. dst=0. ! do l=ltop,lb tkl=tk(l) diftl=(trefk(l)-tkl )*tauk difql=(qrefk(l)-qk(l))*tauk avrgtl=(tkl+tkl+diftl) dpot=dprs(l)/avrgtl dst=diftl*dpot+dst dsq=difql*el(l)*dpot+dsq avrgt=avrgtl*dprs(l)+avrgt preck=diftl*dprs(l)+preck dift(l)=diftl difq(l)=difql enddo ! dst=(dst+dst)*cp dsq=dsq+dsq dentpy=dst+dsq avrgt=avrgt/(sumdp+sumdp) ! ! drheat=preck*cp/avrgt drheat=(preck*sm+max(1.e-7,preck)*(1.-sm))*cp/avrgt drheat=max(drheat,1.e-20) efi=efifc*dentpy/drheat !----------------------------------------------------------------------- efi=min(efi,1.) efi=max(efi,efimn) !----------------------------------------------------------------------- ! enddo cloud_efficiency ! !----------------------------------------------------------------------- ! !----------------------------------------------------------------------- !---------------------- deep convection -------------------------------- !----------------------------------------------------------------------- if(dentpy>=epsntp.and.preck>epspr.and..not.nodeep) then ! iswap=0 ! deep convection, no swap cldefi=efi ! if(sm.gt.0.5) then fefi=(cldefi-efimn)*slopes+efmnts else fefi=(cldefi-efimn)*slopel+efmntl endif ! fefi=(dentpy-epsntp)*fefi/dentpy preck=preck*fefi ! !*** update precipitation and tendencies of temperature and moisture ! cup=preck*cprlg pcpcol=cup ! do l=ltop,lb dtdt(l)=dift(l)*fefi*rdtcnvc dqdt(l)=difq(l)*fefi*rdtcnvc enddo !----------------------------------------------------------------------- if(mmntdeep) then facuv=fefi*rdtcnvc*uvscald if(l0.gt.ltop.and.l0.lt.lb) then ubar=0. vbar=0. sumdp=0. ! do l=l0,lb ubar=u(l)*dprs(l)+ubar vbar=v(l)*dprs(l)+vbar sumdp=dprs(l)+sumdp enddo ! rdpsum=1./sumdp ubar=ubar*rdpsum vbar=vbar*rdpsum ! do l=l0,lb dudt(l)=(ubar-u(l))*facuv dvdt(l)=(vbar-v(l))*facuv enddo ! dum=ubar-u(l0) dvm=vbar-v(l0) ! do l=ltop,l0-1 dudt(l)=(pk(l)-pkt)*dum*rdp0t*facuv dvdt(l)=(pk(l)-pkt)*dvm*rdp0t*facuv enddo else ubar=0. vbar=0. sumdp=0. ! do l=ltop,lb ubar=u(l)*dprs(l)+ubar vbar=v(l)*dprs(l)+vbar sumdp=dprs(l)+sumdp enddo ! rdpsum=1./sumdp ubar=ubar*rdpsum vbar=vbar*rdpsum ! do l=ltop,lb dudt(l)=(ubar-u(l))*facuv dvdt(l)=(vbar-v(l))*facuv enddo endif endif !----------------------------------------------------------------------- else !----------------------------------------------------------------------- !*** reduce the cloud top !----------------------------------------------------------------------- ! ! ltop=ltop+3 !iterate cloud top ! ptop=prsmid(ltop) !iterate cloud top ! depmin=psh*psfc*rsfcp !iterate cloud top ! depth=pbot-ptop !iterate cloud top !*** !*** iterate deep convection procedure if needed !*** ! if(depth>=depmin)then !iterate cloud top ! go to 300 !iterate cloud top ! endif !iterate cloud top ! ! cldefi=avgefi cldefi=efimn*sm+stefi*(1.-sm) !*** !*** search for shallow cloud top !*** ! ltsh=lbot ! lbm1=lbot-1 ! pbtk=pk(lbot) ! depmin=psh*psfc*rsfcp ! ptpk=pbtk-depmin ptpk=max(pshu, pk(lbot)-depmin) !*** !*** cloud top is the level just below pbtk-psh or just below pshu !*** do l=1,lmh if(pk(l)<=ptpk)ltop=l+1 enddo ! ! ptpk=pk(ltop) !!*** !!*** highest level allowed is level just below pshu !!*** ! if(ptpk<=pshu)then !! ! do l=1,lmh ! if(pk(l)<=pshu)lshu=l+1 ! enddo !! ! ltop=lshu ! ptpk=pk(ltop) ! endif ! ! if(ltop>=lbot)then !!!!!! lbot=0 ! ltop=lmh !!!!!! pbot=pk(lbot) ! ptop=pk(ltop) ! pbot=ptop ! go to 600 ! endif ! ! ltp1=ltop+1 ! ltp2=ltop+2 !! ! do l=ltop,lbot ! qsatk(l)=pq0/pk(l)*exp(a2*(tk(l)-a3)/(tk(l)-a4)) ! enddo !! ! rhh=qk(ltop)/qsatk(ltop) ! rhmax=0. ! ltsh=ltop !! ! do l=ltp1,lbm1 ! rhl=qk(l)/qsatk(l) ! drhdp=(rhh-rhl)/(pk(l-1)-pk(l)) !! ! if(drhdp>rhmax)then ! ltsh=l-1 ! rhmax=drhdp ! endif !! ! rhh=rhl ! enddo ! !----------------------------------------------------------------------- !-- make shallow cloud top a function of virtual temperature excess (dtv) !----------------------------------------------------------------------- ! ltp1=lbot do l=lbot-1,ltop,-1 if (dtv(l) > 0.) then ltp1=l else exit endif enddo ltop=min(ltp1,lbot) !*** !*** cloud must be at least two layers thick !*** ! if(lbot-ltop<2)ltop=lbot-2 (eliminate this criterion) ! !-- end: buoyancy check (24 aug 2006) ! iswap=1 ! failed deep convection, shallow swap point ptop=pk(ltop) shallow=.true. deep=.false. ! endif !dcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcd !dcdcdcdcdcdcdc end of deep convection dcdcdcdcdcdcd !dcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcdcd ! !----------------------------------------------------------------------- !----------------------------------------------------------------------- 600 continue !----------------------------------------------------------------------- !----------------------------------------------------------------------- ! !----------------gather shallow convection points----------------------- ! ! if(ptop<=pbot-pno.and.ltop<=lbot-2)then ! depmin=psh*psfc*rsfcp !! !! if(lpbllmh-1)lbot=lmh-1 !! pbot=prsmid(lbot) !! ! if(ptop+1.>=pbot-depmin)shallow=.true. ! else ! lbot=0 ! ltop=lm ! endif ! !*********************************************************************** !----------------------------------------------------------------------- !*** begin debugging convection if(print_diag)then write(6,"(a,2i3,l2,3e12.4)") & '{cu2a lbot,ltop,shallow,pbot,ptop,depmin = ' & ,lbot,ltop,shallow,pbot,ptop,depmin endif !*** end debugging convection !----------------------------------------------------------------------- ! if(.not.shallow)return ! !----------------------------------------------------------------------- !*********************************************************************** !----------------------------------------------------------------------- !scscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscs !scscscscscscsc shallow convection cscscscscscscscscscs !scscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscscs !----------------------------------------------------------------------- do l=1,lmh pk(l)=prsmid(l) tk(l)=t(l) qk(l)=q(l) trefk(l)=t(l) qrefk(l)=q(l) qsatk(l)=pq0/pk(l)*exp(a2*(tk(l)-a3)/(tk(l)-a4)) rxnerk(l)=rxner(l) thvref(l)=tk(l)*rxnerk(l)*(qk(l)*p608+1.) ! if(tk(l)>=tiw)then el(l)=elwv else el(l)=eliv endif enddo ! !----------------------------------------------------------------------- !-- begin: raise cloud top if avg rh>rhshmax and cape>0 ! rhshmax=rh at cloud base associated with a dsp of pone !----------------------------------------------------------------------- ! tlev2=t(lbot)*((pk(lbot)-pone)/pk(lbot))**cappa qsat1=pq0/pk(lbot)*exp(a2*(t(lbot)-a3)/(tk(lbot)-a4)) qsat2=pq0/(pk(lbot)-pone)*exp(a2*(tlev2-a3)/(tlev2-a4)) rhshmax=qsat2/qsat1 sumdp=0. rhavg=0. ! do l=lbot,ltop,-1 rhavg=rhavg+dprs(l)*qk(l)/qsatk(l) sumdp=sumdp+dprs(l) enddo ! if (rhavg/sumdp > rhshmax) then ltsh=ltop do l=ltop-1,1,-1 rhavg=rhavg+dprs(l)*qk(l)/qsatk(l) sumdp=sumdp+dprs(l) if (cpe(l) > 0.) then ltsh=l else exit endif if (rhavg/sumdp <= rhshmax) exit if (pk(l) <= pshu) exit enddo ltop=ltsh !swapnomoist iswap=0 ! old cloud for moist clouds endif ! !-- end: raise cloud top if avg rh>rhshmax and cape>0 ! !---------------------------shallow cloud top--------------------------- lbm1=lbot-1 ptpk=ptop ltp1=ltop-1 depth=pbot-ptop !----------------------------------------------------------------------- !zj if(depth.ge.depmin*0.50) then !zj if(depth.ge.depmin*0.25) then !zj if(depth.ge.depmin*0.625) then ! if(depth.ge.depmin*0.75) then !zj if(depth.ge.depmin*0.85) then ! plume=.true. ! endif !----------------------------------------------------------------------- !*** begin debugging convection if(print_diag)then write(6,"(a,4e12.4)") '{cu2b pbot,ptop,depth,depmin= ' & ,pbot,ptop,depth,depmin endif !*** end debugging convection !----------------------------------------------------------------------- ! !bsf if(depthpbot-pno.or.ltop>lbot-2)then lbot=0 ltop=lm ptop=pbot return endif !----------------------------------------------------------------------- !*** new cloud at all shallow points !----------------------------------------------------------------------- !zj if(newall) go to 810 ! new cloud at all shallow points !zj if(newall.and.sm.lt.0.5) go to 810 ! new cloud at land points if(newall.and.plume) go to 810 ! new cloud at plume points !zj if(newall.and.plume.and.sm.lt.0.5) go to 810 ! new cloud at plume land points !----------------------------------------------------------------------- !*** new cloud at swap shallow points !----------------------------------------------------------------------- !zj if(newswap.and.iswap.gt.0) go to 810 ! new cloud only at swap pts. !zj if(newswap.and.iswap.gt.0.and.sm.lt.0.5) go to 810 ! new cloud only at swap pts. !zj if(newswap.and.iswap.gt.0.and.plume) go to 810 ! new cloud if plume at swap pts. !zj if(newswap.and.iswap.gt.0.and.plume.and.sm.lt.0.5) go to 810 ! new cloud only at swap pts. !----------------------------------------------------------------------- ! !--------------scaling potential temperature & table index at top------- ! thtpk=t(ltp1)*rxner(ltp1) ! tthk=(thtpk-thl)*rdth qqk =tthk-aint(tthk) it =int(tthk)+1 ! if(it<1)then it=1 qqk=0. endif ! if(it>=jtb)then it=jtb-1 qqk=0. endif ! !--------------base and scaling factor for spec. humidity at top-------- ! bqs00k=qs0(it) sqs00k=sqs(it) bqs10k=qs0(it+1) sqs10k=sqs(it+1) ! !--------------scaling spec. humidity & table index at top-------------- ! bqk=(bqs10k-bqs00k)*qqk+bqs00k sqk=(sqs10k-sqs00k)*qqk+sqs00k ! ! tqk=(q(ltop)-bqk)/sqk*rdq tqk=(q(ltp1)-bqk)/sqk*rdq ! ppk=tqk-aint(tqk) iq =int(tqk)+1 ! if(iq<1)then iq=1 ppk=0. endif ! if(iq>=itb)then iq=itb-1 ppk=0. endif ! !----------------cloud top saturation point pressure-------------------- part1=(ptbl(iq+1,it)-ptbl(iq,it))*ppk part2=(ptbl(iq,it+1)-ptbl(iq,it))*qqk part3=(ptbl(iq ,it )-ptbl(iq+1,it ) & -ptbl(iq ,it+1)+ptbl(iq+1,it+1))*ppk*qqk ptpk=ptbl(iq,it)+part1+part2+part3 !----------------------------------------------------------------------- dpmix=ptpk-psp if(abs(dpmix).lt.3000.)dpmix=-3000. !----------------temperature profile slope------------------------------ smix=(thtpk-thbt)/dpmix*stabs ! treflo=trefk(lbot+1) pklo=pk(lbot+1) pkhi=pk(lbot) rxnerlo=rxnerk(lbot+1) rxnerhi=rxnerk(lbot) ! lmid=.5*(lbot+ltop) ! do l=lbot,ltop,-1 treflo=((pkhi-pklo)*smix+treflo*rxnerlo)/rxnerhi trefk(l)=treflo if(l<=lmid) trefk(l)=max(trefk(l),tk(l)+dtshal) rxnerlo=rxnerhi pklo=pkhi rxnerhi=rxnerk(l-1) pkhi=pk(l-1) enddo !----------------temperature reference profile correction--------------- sumdt=0. sumdp=0. ! do l=ltop,lbot sumdt=(tk(l)-trefk(l))*dprs(l)+sumdt sumdp=sumdp+dprs(l) enddo ! rdpsum=1./sumdp fpk(lbot)=trefk(lbot) ! tcorr=sumdt*rdpsum ! do l=ltop,lbot trfkl =trefk(l)+tcorr trefk(l)=trfkl fpk (l)=trfkl enddo !----------------humidity profile equations----------------------------- psum =0. qsum =0. potsum=0. qotsum=0. otsum =0. dst =0. fptk =fpk(ltop) ! do l=ltop,lbot dpkl =fpk(l)-fptk psum =dpkl *dprs(l)+psum qsum =qk(l)*dprs(l)+qsum rtbar =2./(trefk(l)+tk(l)) otsum =dprs(l)*rtbar+otsum potsum=dpkl *rtbar*dprs(l)+potsum qotsum=qk(l) *rtbar*dprs(l)+qotsum dst =(trefk(l)-tk(l))*rtbar*dprs(l)/el(l)+dst enddo ! psum =psum*rdpsum qsum =qsum*rdpsum rotsum=1./otsum potsum=potsum*rotsum qotsum=qotsum*rotsum dst =dst*rotsum*cp ! !----------------------------------------------------------------------- !*** begin debugging convection if(print_diag)then write(6,"(a,5e12.4)") '{cu2c dst,psum,qsum,potsum,qotsum = ' & ,dst,psum,qsum,potsum,qotsum endif !*** end debugging convection !----------------------------------------------------------------------- !*** if upward transport of temperature go to new cloud !----------------------------------------------------------------------- !zj if(newupup.and.dst.gt.0.) go to 810 ! new shallow cloud for both heat and moisture up !zj if(newupup.and.dst.gt.0..and.sm.lt.0.5) go to 810 ! new shallow cloud for both heat and moisture up if(newupup.and.dst.gt.0..and.plume) go to 810 ! new shallow cloud if plume for both heat and moisture up !zj if(newupup.and.dst.gt.0..and.plume.and.sm.lt.0.5) go to 810 ! new shallow cloud for both heat and moisture up !----------------------------------------------------------------------- !*** otherwise old cloud !----------------------------------------------------------------------- if(dst.gt.0.) then lbot=0 ltop=lm ptop=pbot return endif !----------------------------------------------------------------------- !*** otherwise continue with old cloud !----------------ensure positive entropy change------------------------- dstq=dst*epsdn !----------------check for isothermal atmosphere------------------------ den=potsum-psum ! if(-den/psum<5.e-5)then lbot=0 ltop=lm ptop=pbot return !----------------slope of the reference humidity profile---------------- ! else dqref=(qotsum-dstq-qsum)/den endif ! !-------------- humidity does not increase with height------------------ ! if(dqref<0.)then lbot=0 ltop=lm ptop=pbot return endif ! !----------------humidity at the cloud top------------------------------ ! qrftp=qsum-dqref*psum ! !----------------humidity profile--------------------------------------- ! do l=ltop,lbot qrfkl=(fpk(l)-fptk)*dqref+qrftp ! !*** too dry clouds not allowed ! tnew=(trefk(l)-tk(l))*tauksc+tk(l) qsatk(l)=pq0/pk(l)*exp(a2*(tnew-a3)/(tnew-a4)) qnew=(qrfkl-qk(l))*tauksc+qk(l) ! if(qnewqsatk(l)*rhhsc)then lbot=0 ltop=lm ptop=pbot return endif ! thvref(l)=trefk(l)*rxnerk(l)*(qrfkl*p608+1.) qrefk(l)=qrfkl enddo ! !------------------ eliminate clouds with bottoms too dry -------------- !! ! qnew=(qrefk(lbot)-qk(lbot))*tauksc+qk(lbot) !! ! if(qnew