! MODULE surf_landice_mod IMPLICIT NONE CONTAINS ! !**************************************************************************************** ! SUBROUTINE surf_landice(itime, dtime, knon, knindex, & rlon, rlat, debut, lafin, & rmu0, lwdownm, albedo, pphi1, & swnet, lwnet, tsurf, p1lay, & cdragh, cdragm, precip_rain, precip_snow, precip_bs, temp_air, spechum, & AcoefH, AcoefQ, BcoefH, BcoefQ, & AcoefU, AcoefV, BcoefU, BcoefV, & AcoefQBS, BcoefQBS, & ps, u1, v1, gustiness, rugoro, pctsrf, & snow, qsurf, qsol, qbs1, agesno, & tsoil, z0m, z0h, SFRWL, alb_dir, alb_dif, evap, fluxsens, fluxlat, fluxbs, & tsurf_new, dflux_s, dflux_l, & alt, slope, cloudf, & snowhgt, qsnow, to_ice, sissnow, & alb3, runoff, & flux_u1, flux_v1) USE dimphy USE geometry_mod, ONLY : longitude,latitude USE surface_data, ONLY : type_ocean, calice, calsno, landice_opt, iflag_albcalc USE fonte_neige_mod, ONLY : fonte_neige,run_off_lic,fqcalving_global,ffonte_global,fqfonte_global,runofflic_global USE cpl_mod, ONLY : cpl_send_landice_fields USE calcul_fluxs_mod USE phys_local_var_mod, ONLY : zxrhoslic, zxustartlic USE phys_output_var_mod, ONLY : snow_o,zfra_o !FC USE ioipsl_getin_p_mod, ONLY : getin_p USE blowing_snow_ini_mod, ONLY : zeta_bs, pbst_bs, prt_bs, iflag_saltation_bs #ifdef CPP_INLANDSIS USE surf_inlandsis_mod, ONLY : surf_inlandsis #endif USE indice_sol_mod ! INCLUDE "indicesol.h" INCLUDE "dimsoil.h" INCLUDE "YOMCST.h" INCLUDE "clesphys.h" ! Input variables !**************************************************************************************** INTEGER, INTENT(IN) :: itime, knon INTEGER, DIMENSION(klon), INTENT(in) :: knindex REAL, INTENT(in) :: dtime REAL, DIMENSION(klon), INTENT(IN) :: swnet ! net shortwave radiance REAL, DIMENSION(klon), INTENT(IN) :: lwnet ! net longwave radiance REAL, DIMENSION(klon), INTENT(IN) :: tsurf REAL, DIMENSION(klon), INTENT(IN) :: p1lay REAL, DIMENSION(klon), INTENT(IN) :: cdragh, cdragm REAL, DIMENSION(klon), INTENT(IN) :: precip_rain, precip_snow, precip_bs REAL, DIMENSION(klon), INTENT(IN) :: temp_air, spechum REAL, DIMENSION(klon), INTENT(IN) :: AcoefH, AcoefQ REAL, DIMENSION(klon), INTENT(IN) :: BcoefH, BcoefQ REAL, DIMENSION(klon), INTENT(IN) :: AcoefU, AcoefV, BcoefU, BcoefV REAL, DIMENSION(klon), INTENT(IN) :: AcoefQBS, BcoefQBS REAL, DIMENSION(klon), INTENT(IN) :: ps REAL, DIMENSION(klon), INTENT(IN) :: u1, v1, gustiness, qbs1 REAL, DIMENSION(klon), INTENT(IN) :: rugoro REAL, DIMENSION(klon,nbsrf), INTENT(IN) :: pctsrf LOGICAL, INTENT(IN) :: debut !true if first step LOGICAL, INTENT(IN) :: lafin !true if last step REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat REAL, DIMENSION(klon), INTENT(IN) :: rmu0 REAL, DIMENSION(klon), INTENT(IN) :: lwdownm !ylwdown REAL, DIMENSION(klon), INTENT(IN) :: albedo !mean albedo REAL, DIMENSION(klon), INTENT(IN) :: pphi1 REAL, DIMENSION(klon), INTENT(IN) :: alt !mean altitude of the grid box REAL, DIMENSION(klon), INTENT(IN) :: slope !mean slope in grid box REAL, DIMENSION(klon), INTENT(IN) :: cloudf !total cloud fraction ! In/Output variables !**************************************************************************************** REAL, DIMENSION(klon), INTENT(INOUT) :: snow, qsol REAL, DIMENSION(klon), INTENT(INOUT) :: agesno REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil ! Output variables !**************************************************************************************** REAL, DIMENSION(klon), INTENT(OUT) :: qsurf REAL, DIMENSION(klon), INTENT(OUT) :: z0m, z0h !albedo SB >>> ! REAL, DIMENSION(klon), INTENT(OUT) :: alb1 ! new albedo in visible SW interval ! REAL, DIMENSION(klon), INTENT(OUT) :: alb2 ! new albedo in near IR interval REAL, DIMENSION(6), INTENT(IN) :: SFRWL REAL, DIMENSION(klon,nsw), INTENT(OUT) :: alb_dir,alb_dif !albedo SB <<< REAL, DIMENSION(klon), INTENT(OUT) :: evap, fluxsens, fluxlat REAL, DIMENSION(klon), INTENT(OUT) :: fluxbs REAL, DIMENSION(klon), INTENT(OUT) :: tsurf_new REAL, DIMENSION(klon), INTENT(OUT) :: dflux_s, dflux_l REAL, DIMENSION(klon), INTENT(OUT) :: flux_u1, flux_v1 REAL, DIMENSION(klon), INTENT(OUT) :: alb3 REAL, DIMENSION(klon), INTENT(OUT) :: qsnow !column water in snow [kg/m2] REAL, DIMENSION(klon), INTENT(OUT) :: snowhgt !Snow height (m) REAL, DIMENSION(klon), INTENT(OUT) :: to_ice REAL, DIMENSION(klon), INTENT(OUT) :: sissnow REAL, DIMENSION(klon), INTENT(OUT) :: runoff !Land ice runoff ! Local variables !**************************************************************************************** REAL, DIMENSION(klon) :: soilcap, soilflux REAL, DIMENSION(klon) :: cal, beta, dif_grnd REAL, DIMENSION(klon) :: zfra, alb_neig REAL, DIMENSION(klon) :: radsol REAL, DIMENSION(klon) :: u0, v0, u1_lay, v1_lay, ustar INTEGER :: i,j,nt REAL, DIMENSION(klon) :: fqfonte,ffonte REAL, DIMENSION(klon) :: run_off_lic_frac REAL, DIMENSION(klon) :: emis_new !Emissivity REAL, DIMENSION(klon) :: swdown,lwdown REAL, DIMENSION(klon) :: precip_snow_adv, snow_adv !Snow Drift precip./advection (not used in inlandsis) REAL, DIMENSION(klon) :: erod !erosion of surface snow (flux, kg/m2/s like evap) REAL, DIMENSION(klon) :: zsl_height, wind_velo !surface layer height, wind spd REAL, DIMENSION(klon) :: dens_air, snow_cont_air !air density; snow content air REAL, DIMENSION(klon) :: alb_soil !albedo of underlying ice REAL, DIMENSION(klon) :: pexner !Exner potential REAL :: pref REAL, DIMENSION(klon,nsoilmx) :: tsoil0 !modif REAL :: dtis ! subtimestep LOGICAL :: debut_is, lafin_is ! debut and lafin for inlandsis CHARACTER (len = 20) :: modname = 'surf_landice' CHARACTER (len = 80) :: abort_message REAL,DIMENSION(klon) :: alb1,alb2 REAL,DIMENSION(klon) :: precip_totsnow, evap_totsnow REAL, DIMENSION (klon,6) :: alb6 REAL :: rho0, rhoice, ustart0, hsalt, esalt, rhod REAL :: lambdasalt,fluxsalt, csalt, nunu, aa, bb, cc REAL :: tau_dens, tau_dens0, tau_densmin, rhomax, rhohard REAL, DIMENSION(klon) :: ws1, rhos, ustart, qsalt ! End definition !**************************************************************************************** !FC !FC REAL,SAVE :: alb_vis_sno_lic !$OMP THREADPRIVATE(alb_vis_sno_lic) REAL,SAVE :: alb_nir_sno_lic !$OMP THREADPRIVATE(alb_nir_sno_lic) LOGICAL, SAVE :: firstcall = .TRUE. !$OMP THREADPRIVATE(firstcall) !FC firtscall initializations !****************************************************************************************** IF (firstcall) THEN alb_vis_sno_lic=0.77 CALL getin_p('alb_vis_sno_lic',alb_vis_sno_lic) PRINT*, 'alb_vis_sno_lic',alb_vis_sno_lic alb_nir_sno_lic=0.77 CALL getin_p('alb_nir_sno_lic',alb_nir_sno_lic) PRINT*, 'alb_nir_sno_lic',alb_nir_sno_lic firstcall=.false. ENDIF !****************************************************************************************** ! Initialize output variables alb3(:) = 999999. alb2(:) = 999999. alb1(:) = 999999. fluxbs(:)=0. runoff(:) = 0. !**************************************************************************************** ! Calculate total absorbed radiance at surface ! !**************************************************************************************** radsol(:) = 0.0 radsol(1:knon) = swnet(1:knon) + lwnet(1:knon) !**************************************************************************************** !**************************************************************************************** ! landice_opt = 0 : soil_model, calcul_flux, fonte_neige, ... ! landice_opt = 1 : prepare and call INterace Lmdz SISvat (INLANDSIS) !**************************************************************************************** IF (landice_opt .EQ. 1) THEN !**************************************************************************************** ! CALL to INLANDSIS interface !**************************************************************************************** #ifdef CPP_INLANDSIS debut_is=debut lafin_is=.false. ! Suppose zero surface speed u0(:) = 0.0 v0(:) = 0.0 CALL calcul_flux_wind(knon, dtime, & u0, v0, u1, v1, gustiness, cdragm, & AcoefU, AcoefV, BcoefU, BcoefV, & p1lay, temp_air, & flux_u1, flux_v1) ! Set constants and compute some input for SISVAT ! = 1000 hPa ! and calculate incoming flux for SW and LW interval: swdown, lwdown swdown(:) = 0.0 lwdown(:) = 0.0 snow_cont_air(:) = 0. ! the snow content in air is not a prognostic variable of the model alb_soil(:) = 0.4 ! before albedo(:) but here it is the ice albedo that we have to set ustar(:) = 0. pref = 100000. DO i = 1, knon swdown(i) = swnet(i)/(1-albedo(i)) lwdown(i) = lwdownm(i) wind_velo(i) = u1(i)**2 + v1(i)**2 wind_velo(i) = wind_velo(i)**0.5 pexner(i) = (p1lay(i)/pref)**(RD/RCPD) dens_air(i) = p1lay(i)/RD/temp_air(i) ! dry air density zsl_height(i) = pphi1(i)/RG tsoil0(i,:) = tsoil(i,:) ustar(i)= (cdragm(i)*(wind_velo(i)**2))**0.5 END DO dtis=dtime IF (lafin) THEN lafin_is=.true. END IF CALL surf_inlandsis(knon, rlon, rlat, knindex, itime, dtis, debut_is, lafin_is,& rmu0, swdown, lwdown, albedo, pexner, ps, p1lay, precip_rain, precip_snow, & zsl_height, wind_velo, ustar, temp_air, dens_air, spechum, tsurf,& rugoro, snow_cont_air, alb_soil, alt, slope, cloudf, & radsol, qsol, tsoil0, snow, zfra, snowhgt, qsnow, to_ice, sissnow,agesno, & AcoefH, AcoefQ, BcoefH, BcoefQ, cdragm, cdragh, & run_off_lic, fqfonte, ffonte, evap, erod, fluxsens, fluxlat,dflux_s, dflux_l, & tsurf_new, alb1, alb2, alb3, alb6, & emis_new, z0m, z0h, qsurf) debut_is=.false. ! Treatment of snow melting and calving ! for consistency with standard LMDZ, add calving to run_off_lic run_off_lic(:)=run_off_lic(:) + to_ice(:) DO i = 1, knon ffonte_global(knindex(i),is_lic) = ffonte(i) fqfonte_global(knindex(i),is_lic) = fqfonte(i)! net melting= melting - refreezing fqcalving_global(knindex(i),is_lic) = to_ice(i) ! flux runofflic_global(knindex(i)) = run_off_lic(i) ENDDO ! Here, we assume that the calving term is equal to the to_ice term ! (no ice accumulation) #else abort_message='Pb de coherence: landice_opt = 1 mais CPP_INLANDSIS = .false.' CALL abort_physic(modname,abort_message,1) #endif ELSE !**************************************************************************************** ! Soil calculations ! !**************************************************************************************** ! EV: use calbeta CALL calbeta(dtime, is_lic, knon, snow, qsol, beta, cal, dif_grnd) ! use soil model and recalculate properly cal IF (soil_model) THEN CALL soil(dtime, is_lic, knon, snow, tsurf, qsol, & & longitude(knindex(1:knon)), latitude(knindex(1:knon)), tsoil, soilcap, soilflux) cal(1:knon) = RCPD / soilcap(1:knon) radsol(1:knon) = radsol(1:knon) + soilflux(1:knon) ELSE cal = RCPD * calice WHERE (snow > 0.0) cal = RCPD * calsno ENDIF !**************************************************************************************** ! Calulate fluxes ! !**************************************************************************************** ! beta(:) = 1.0 ! dif_grnd(:) = 0.0 ! Suppose zero surface speed u0(:)=0.0 v0(:)=0.0 u1_lay(:) = u1(:) - u0(:) v1_lay(:) = v1(:) - v0(:) CALL calcul_fluxs(knon, is_lic, dtime, & tsurf, p1lay, cal, beta, cdragh, cdragh, ps, & precip_rain, precip_snow, snow, qsurf, & radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, & 1.,AcoefH, AcoefQ, BcoefH, BcoefQ, & tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l) CALL calcul_flux_wind(knon, dtime, & u0, v0, u1, v1, gustiness, cdragm, & AcoefU, AcoefV, BcoefU, BcoefV, & p1lay, temp_air, & flux_u1, flux_v1) !**************************************************************************************** ! Calculate albedo ! !**************************************************************************************** CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)) ! EV: following lines are obsolete since we set alb1 and alb2 to constant values ! I therefore comment them ! alb1(1:knon) = alb_neig(1:knon)*zfra(1:knon) + & ! 0.6 * (1.0-zfra(1:knon)) ! !IM: plusieurs choix/tests sur l'albedo des "glaciers continentaux" ! alb1(1 : knon) = 0.6 !IM cf FH/GK ! alb1(1 : knon) = 0.82 ! alb1(1 : knon) = 0.77 !211003 Ksta0.77 ! alb1(1 : knon) = 0.8 !KstaTER0.8 & LMD_ARMIP5 !IM: KstaTER0.77 & LMD_ARMIP6 ! Attantion: alb1 and alb2 are not the same! alb1(1:knon) = alb_vis_sno_lic alb2(1:knon) = alb_nir_sno_lic !**************************************************************************************** ! Rugosity ! !**************************************************************************************** z0m = z0m_landice z0h = z0h_landice !z0m = SQRT(z0m**2+rugoro**2) ! Simple blowing snow param if (ok_bs) then ustart0 = 0.211 rhoice = 920.0 rho0 = 200.0 rhomax=450.0 rhohard=400.0 tau_dens0=86400.0*10. ! 10 days by default, in s tau_densmin=86400.0 ! 1 days according to in situ obs by C. Amory ! computation of threshold friction velocity ! which depends on surface snow density do i = 1, knon ! estimation of snow density ! snow density increases with snow age and ! increases even faster in case of sedimentation of blowing snow tau_dens=max(tau_densmin, tau_dens0*exp(-abs(precip_bs(i))/pbst_bs-abs(precip_rain(i))/prt_bs)) rhos(i)=rho0+(rhohard-rho0)*(1.-exp(-agesno(i)*86400.0/tau_dens)) ! blowing snow flux formula used in MAR ws1(i)=(u1(i)**2+v1(i)**2)**0.5 ustar(i)=(cdragm(i)*(u1(i)**2+v1(i)**2))**0.5 ustart(i)=ustart0*exp(max(rhoice/rho0-rhoice/rhos(i),0.))*exp(max(0.,rhos(i)-rhomax)) ! we have multiplied by exp to prevent erosion when rhos>rhomax (usefull till ! rhohard<450) enddo ! computation of qbs at the top of the saltation layer ! two formulations possible ! pay attention that qbs is a mixing ratio and has to be converted ! to specific content if (iflag_saltation_bs .eq. 1) then ! expression from CRYOWRF (Sharma et al. 2022) aa=2.6 bb=2.5 cc=2.0 lambdasalt=0.45 do i =1, knon rhod=p1lay(i)/RD/temp_air(i) nunu=max(ustar(i)/ustart(i),1.e-3) fluxsalt=rhod/RG*(ustar(i)**3)*(1.-nunu**(-2)) * & (aa+bb*nunu**(-2)+cc*nunu**(-1)) csalt=fluxsalt/(2.8*ustart(i)) hsalt=0.08436*ustar(i)**1.27 qsalt(i)=1./rhod*csalt*lambdasalt*RG/(max(ustar(i)**2,1E-6)) & * exp(-lambdasalt*RG*hsalt/max(ustar(i)**2,1E-6)) qsalt(i)=max(qsalt(i),0.) enddo else ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001) do i=1, knon esalt=1./(3.25*max(ustar(i),0.001)) hsalt=0.08436*ustar(i)**1.27 qsalt(i)=(max(ustar(i)**2-ustart(i)**2,0.))/(RG*hsalt)*esalt !ep=qsalt*cdragm(i)*sqrt(u1(i)**2+v1(i)**2) enddo endif ! calculation of erosion (emission flux towards the first atmospheric level) ! consistent with implicit resolution of turbulent mixing equation do i=1, knon rhod=p1lay(i)/RD/temp_air(i) fluxbs(i)=rhod*ws1(i)*cdragm(i)*zeta_bs*(AcoefQBS(i)-qsalt(i)) & / (1.-rhod*ws1(i)*zeta_bs*BcoefQBS(i)*dtime) !fluxbs(i)= zeta_bs*rhod*ws1(i)*cdragm(i)*(qbs1(i)-qsalt(i)) enddo ! for outputs do j = 1, knon i = knindex(j) zxustartlic(i) = ustart(j) zxrhoslic(i) = rhos(j) enddo endif !**************************************************************************************** ! Calculate surface snow amount ! !**************************************************************************************** IF (ok_bs) THEN precip_totsnow(:)=precip_snow(:)+precip_bs(:) evap_totsnow(:)=evap(:)-fluxbs(:) ! flux bs is positive towards the surface (snow erosion) ELSE precip_totsnow(:)=precip_snow(:) evap_totsnow(:)=evap(:) ENDIF CALL fonte_neige(knon, is_lic, knindex, dtime, & tsurf, precip_rain, precip_totsnow, & snow, qsol, tsurf_new, evap_totsnow) WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0. zfra(1:knon) = MAX(0.0,MIN(1.0,snow(1:knon)/(snow(1:knon)+10.0))) END IF ! landice_opt !**************************************************************************************** ! Send run-off on land-ice to coupler if coupled ocean. ! run_off_lic has been calculated in fonte_neige or surf_inlandsis ! If landice_opt>=2, corresponding call is done from surf_land_orchidee !**************************************************************************************** IF (type_ocean=='couple' .AND. landice_opt .LT. 2) THEN ! Compress fraction where run_off_lic is active (here all pctsrf(is_lic)) run_off_lic_frac(:)=0.0 DO j = 1, knon i = knindex(j) run_off_lic_frac(j) = pctsrf(i,is_lic) ENDDO CALL cpl_send_landice_fields(itime, knon, knindex, run_off_lic, run_off_lic_frac) ENDIF ! transfer runoff rate [kg/m2/s](!) to physiq for output runoff(1:knon)=run_off_lic(1:knon)/dtime snow_o=0. zfra_o = 0. DO j = 1, knon i = knindex(j) snow_o(i) = snow(j) zfra_o(i) = zfra(j) ENDDO !albedo SB >>> select case(NSW) case(2) alb_dir(1:knon,1)=alb1(1:knon) alb_dir(1:knon,2)=alb2(1:knon) case(4) alb_dir(1:knon,1)=alb1(1:knon) alb_dir(1:knon,2)=alb2(1:knon) alb_dir(1:knon,3)=alb2(1:knon) alb_dir(1:knon,4)=alb2(1:knon) case(6) alb_dir(1:knon,1)=alb1(1:knon) alb_dir(1:knon,2)=alb1(1:knon) alb_dir(1:knon,3)=alb1(1:knon) alb_dir(1:knon,4)=alb2(1:knon) alb_dir(1:knon,5)=alb2(1:knon) alb_dir(1:knon,6)=alb2(1:knon) IF ((landice_opt .EQ. 1) .AND. (iflag_albcalc .EQ. 2)) THEN alb_dir(1:knon,1)=alb6(1:knon,1) alb_dir(1:knon,2)=alb6(1:knon,2) alb_dir(1:knon,3)=alb6(1:knon,3) alb_dir(1:knon,4)=alb6(1:knon,4) alb_dir(1:knon,5)=alb6(1:knon,5) alb_dir(1:knon,6)=alb6(1:knon,6) ENDIF end select alb_dif=alb_dir !albedo SB <<< END SUBROUTINE surf_landice ! !**************************************************************************************** ! END MODULE surf_landice_mod