module lmdz_blowing_snow_sublim_sedim contains subroutine blowing_snow_sublim_sedim(ngrid,nlay,dtime,temp,qv,qbs,pplay,paprs,dtemp_bs,dq_bs,dqbs_bs,bsfl,precip_bs) !============================================================================== ! Routine that calculates the evaporation and sedimentation of blowing snow ! inspired by what is done in lscp_mod ! Etienne Vignon, October 2022 !============================================================================== use lmdz_blowing_snow_ini, only : coef_eva_bs,RTT,RD,RG,expo_eva_bs, fallv_bs, qbsmin use lmdz_blowing_snow_ini, only : RCPD, RLSTT, RLMLT, RLVTT, RVTMP2, tbsmelt, taumeltbs0 USE lmdz_lscp_tools, only : calc_qsat_ecmwf implicit none !++++++++++++++++++++++++++++++++++++++++++++++++++++ ! Declarations !++++++++++++++++++++++++++++++++++++++++++++++++++++ !INPUT !===== integer, intent(in) :: ngrid,nlay real, intent(in) :: dtime real, intent(in), dimension(ngrid,nlay) :: temp real, intent(in), dimension(ngrid,nlay) :: qv real, intent(in), dimension(ngrid,nlay) :: qbs real, intent(in), dimension(ngrid,nlay) :: pplay real, intent(in), dimension(ngrid,nlay+1) :: paprs ! OUTPUT !======== real, intent(out), dimension(ngrid,nlay) :: dtemp_bs real, intent(out), dimension(ngrid,nlay) :: dq_bs real, intent(out), dimension(ngrid,nlay) :: dqbs_bs real, intent(out), dimension(ngrid,nlay+1) :: bsfl real, intent(out), dimension(ngrid) :: precip_bs ! LOCAL !====== integer :: k,i,n real :: zqev0, zqevi, zqevti, zcpair, zcpeau, dqbsmelt real :: dqsedim,precbs, deltaqchaud, zmelt, taumeltbs real :: maxdeltaqchaud real, dimension(ngrid) :: zt,zq,zqbs,qsi,dqsi,qsl, dqsl,qzero,sedim,sedimn real, dimension(ngrid) :: zqbsi,zmqc, zmair, zdz real, dimension(ngrid,nlay) :: velo, zrho !++++++++++++++++++++++++++++++++++++++++++++++++++ ! Initialisation !++++++++++++++++++++++++++++++++++++++++++++++++++ qzero(:)=0. dtemp_bs(:,:)=0. dq_bs(:,:)=0. dqbs_bs(:,:)=0. velo(:,:)=0. zt(:)=0. zq(:)=0. zqbs(:)=0. sedim(:)=0. ! begin of top-down loop DO k = nlay, 1, -1 DO i=1,ngrid zt(i)=temp(i,k) zq(i)=qv(i,k) zqbs(i)=qbs(i,k) ENDDO IF (k.LE.nlay-1) THEN ! thermalization of blowing snow precip coming from above DO i = 1, ngrid zmair(i)=(paprs(i,k)-paprs(i,k+1))/RG ! RVTMP2=rcpv/rcpd-1 zcpair=RCPD*(1.0+RVTMP2*zq(i)) zcpeau=RCPD*RVTMP2 ! zmqc: precipitation mass that has to be thermalized with ! layer's air so that precipitation at the ground has the ! same temperature as the lowermost layer zmqc(i) = (sedim(i))*dtime/zmair(i) ! t(i,k+1)+d_t(i,k+1): new temperature of the overlying layer zt(i) = ( (temp(i,k+1)+dtemp_bs(i,k+1))*zmqc(i)*zcpeau + zcpair*zt(i) ) & / (zcpair + zmqc(i)*zcpeau) ENDDO ELSE DO i = 1, ngrid zmair(i)=(paprs(i,k)-paprs(i,k+1))/RG zmqc(i) = 0. ENDDO ENDIF ! calulation saturation specific humidity CALL CALC_QSAT_ECMWF(ngrid,zt(:),qzero(:),pplay(:,k),RTT,2,.false.,qsi(:),dqsi(:)) CALL CALC_QSAT_ECMWF(ngrid,zt(:),qzero(:),pplay(:,k),RTT,1,.false.,qsl(:),dqsl(:)) ! sublimation calculation ! SUndqvist formula dP/dz=beta*(1-q/qsat)*sqrt(P) DO i = 1, ngrid zrho(i,k) = pplay(i,k) / zt(i) / RD zdz(i) = (paprs(i,k)-paprs(i,k+1)) / (zrho(i,k)*RG) ! BS fall velocity velo(i,k) = fallv_bs IF (zt(i) .GT. RTT) THEN ! if positive celcius temperature, we assume ! that part of the the blowing snow flux melts and evaporates ! vapor, bs, temperature, precip fluxes update zmelt = ((zt(i)-RTT)/(tbsmelt-RTT)) zmelt = MIN(MAX(zmelt,0.),1.) sedimn(i)=sedim(i)*zmelt deltaqchaud=-(sedimn(i)-sedim(i))*(RG/(paprs(i,k)-paprs(i,k+1)))*dtime ! max evap such as celcius temperature cannot become negative maxdeltaqchaud= max(RCPD*(1.0+RVTMP2*(zq(i)+zqbs(i)))*(zt(i)-RTT)/(RLMLT+RLVTT),0.) deltaqchaud=min(max(deltaqchaud,0.),maxdeltaqchaud) zq(i) = zq(i) + deltaqchaud ! melting + evaporation zt(i) = zt(i) - deltaqchaud * (RLMLT+RLVTT)/RCPD/(1.0+RVTMP2*(zq(i)+zqbs(i))) sedim(i)=sedimn(i) ! if temperature still positive, we assume that part of the blowing snow ! already present in the mesh melts and evaporates with a typical time ! constant between taumeltbs0 and 0 IF ( zt(i) .GT. RTT ) THEN taumeltbs=taumeltbs0*exp(-max(0.,(zt(i)-RTT)/(tbsmelt-RTT))) deltaqchaud=min(zqbs(i),zqbs(i)/taumeltbs*dtime) maxdeltaqchaud= max(RCPD*(1.0+RVTMP2*(zq(i)+zqbs(i)))*(zt(i)-RTT)/(RLMLT+RLVTT),0.) deltaqchaud=min(max(deltaqchaud,0.),maxdeltaqchaud) zq(i) = zq(i) + deltaqchaud ! melting + evaporation zt(i) = zt(i) - deltaqchaud * (RLMLT+RLVTT)/RCPD/(1.0+RVTMP2*(zq(i)+zqbs(i))) ! qbs update zqbs(i)=max(zqbs(i)-deltaqchaud,0.) ENDIF ! now sedimentation scheme with an exact numerical resolution ! (correct if fall velocity is constant) zqbsi(i)=zqbs(i) zqbs(i) = zqbsi(i)*exp(-velo(i,k)/zdz(i)*dtime)+sedim(i)/zrho(i,k)/velo(i,k) zqbs(i) = max(zqbs(i),0.) ! flux update sedim(i) = sedim(i) + zrho(i,k)*zdz(i)/dtime*(zqbsi(i)-zqbs(i)) sedim(i) = max(0.,sedim(i)) ELSE ! negative celcius temperature ! Sublimation scheme zqevti = coef_eva_bs*(1.0-zq(i)/qsi(i))*(sedim(i)**expo_eva_bs) & *(paprs(i,k)-paprs(i,k+1))/pplay(i,k)*zt(i)*RD/RG zqevti = MAX(0.0,MIN(zqevti,sedim(i)))*RG*dtime/(paprs(i,k)-paprs(i,k+1)) ! Sublimation limit: we ensure that the whole mesh does not exceed saturation wrt ice zqev0 = MAX(0.0, qsi(i)-zq(i)) zqevi = MIN(zqev0,zqevti) ! New solid precipitation fluxes sedimn(i) = Max(0.,sedim(i) - zqevi*(paprs(i,k)-paprs(i,k+1))/RG/dtime) ! vapor, temperature, precip fluxes update following sublimation zq(i) = zq(i) - (sedimn(i)-sedim(i))*(RG/(paprs(i,k)-paprs(i,k+1)))*dtime zq(i) = max(0., zq(i)) zt(i) = zt(i) & + (sedimn(i)-sedim(i)) & * (RG/(paprs(i,k)-paprs(i,k+1)))*dtime & * RLSTT/RCPD/(1.0+RVTMP2*(zq(i)+zqbs(i))) sedim(i)=sedimn(i) zqbsi(i)=zqbs(i) ! now sedimentation scheme with an exact numerical resolution ! (correct if fall velocity is constant) zqbs(i) = zqbsi(i)*exp(-velo(i,k)/zdz(i)*dtime)+sedim(i)/zrho(i,k)/velo(i,k) zqbs(i) = max(zqbs(i),0.) ! flux update sedim(i) = sedim(i) + zrho(i,k)*zdz(i)/dtime*(zqbsi(i)-zqbs(i)) sedim(i) = max(0.,sedim(i)) ENDIF ! if qbs