subroutine blowing_snow_sublim_sedim(ngrid,nlay,dtime,temp,q,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 blowing_snow_ini_mod, only : coef_eva_bs,RTT,RD,RG,expo_eva_bs, fallv_bs use blowing_snow_ini_mod, only : RCPD, RLSTT, RLMLT, RLVTT, RVTMP2 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) :: q 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 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)=q(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 ! if sedimentation: IF (sedim(i) .GT. 0.) THEN IF (zt(i) .GT. RTT) THEN ! if positive celcius temperature, we assume ! that all the blowing snow melt and evaporate zqevti=sedim(i)*(paprs(i,k)-paprs(i,k+1))/pplay(i,k)*zt(i)*RD/RG ! we ensure that the whole mesh does not exceed saturation wrt liquid zqev0 = MAX(0.0, qsl(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 zq(i) = zq(i) - (sedimn(i)-sedim(i))*(RG/(paprs(i,k)-paprs(i,k+1)))*dtime zq(i) = max(0., zq(i)) !zqbs(i) = zqbs(i) + (sedimn(i)-sedim(i)) * (RG/(paprs(i,k)-paprs(i,k+1)))*dtime !zqbs(i) = max(0., zqbs(i)) ! melting zt(i) = zt(i) & + (sedimn(i)-sedim(i)) & * (RG/(paprs(i,k)-paprs(i,k+1)))*dtime & * RLMLT/RCPD/(1.0+RVTMP2*(zq(i)+zqbs(i))) ! evaporation zt(i) = zt(i) & + (sedimn(i)-sedim(i)) & * (RG/(paprs(i,k)-paprs(i,k+1)))*dtime & * RLVTT/RCPD/(1.0+RVTMP2*(zq(i)+zqbs(i))) ELSE 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)) !zqbs(i) = zqbs(i) + (sedimn(i)-sedim(i)) * (RG/(paprs(i,k)-paprs(i,k+1)))*dtime !zqbs(i) = max(0., zqbs(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))) ENDIF ! sedim update sedim(i)=sedimn(i) ELSE sedim(i)=0. ENDIF ! if sedim > 0 zqbsi(i)=zqbs(i) ENDDO ! loop on ngrid ! Now sedimention scheme ! exact resolution of the conservation equation for qbs with the updated flux from the top (after evap) ! valid only if the fall velocity is constant 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) velo(i,k) = fallv_bs 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)) ENDDO ! old version with bugs ! DO n = 1, niter_bs ! 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) ! velo(i,k) = fallv_bs ! dqsedim = dtime/REAL(niter_bs)/zdz(i)*zqbs(i)*velo(i,k) ! dqice/dt=1/rho*d(rho*wice*qice)/dz ! precbs = MIN(MAX(dqsedim,0.0),zqbs(i)) ! zqbs(i) = MAX(zqbs(i)-1.*precbs , 0.0) ! ENDDO !loop on ngrid ! ENDDO ! loop on niter_bs ! ! add to non sublimated precip ! DO i=1,ngrid ! sedim(i) = sedim(i)+max(0.,zqbsi(i)-zqbs(i))*(paprs(i,k)-paprs(i,k+1))/(RG*dtime) ! ENDDO ! ! Outputs: DO i = 1, ngrid bsfl(i,k)=sedim(i) dqbs_bs(i,k) = zqbs(i)-qbs(i,k) dq_bs(i,k) = zq(i) - q(i,k) dtemp_bs(i,k) = zt(i) - temp(i,k) ENDDO ENDDO ! vertical loop !surface bs flux DO i = 1, ngrid precip_bs(i) = sedim(i) ENDDO return end subroutine blowing_snow_sublim_sedim