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 & #ifdef ISO ,xtprecip_rain, xtprecip_snow,xtspechum,Rland_ice & ,xtsnow,xtsol,xtevap & #endif ) USE dimphy USE lmdz_geometry, 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, zxqsaltlic, tempsmoothlic USE phys_output_var_mod, ONLY: snow_o, zfra_o #ifdef ISO USE fonte_neige_mod, ONLY: xtrun_off_lic USE infotrac_phy, ONLY: ntiso,niso USE isotopes_routines_mod, ONLY: calcul_iso_surf_lic_vectall #ifdef ISOVERIF USE isotopes_mod, ONLY: iso_eau,ridicule USE isotopes_verif_mod #endif #endif !FC USE lmdz_ioipsl_getin_p, ONLY: getin_p USE lmdz_blowing_snow_ini, ONLY: c_esalt_bs, zeta_bs, pbst_bs, prt_bs, rhoice_bs, rhohard_bs USE lmdz_blowing_snow_ini, ONLY: rhofresh_bs, tau_eqsalt_bs, tau_dens0_bs, tau_densmin_bs USE surf_inlandsis_mod, ONLY: surf_inlandsis USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_INLANDSIS USE lmdz_abort_physic, ONLY: abort_physic USE indice_sol_mod USE lmdz_clesphys USE lmdz_yomcst ! INCLUDE "indicesol.h" INCLUDE "dimsoil.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 #ifdef ISO REAL, DIMENSION(ntiso,klon), INTENT(IN) :: xtprecip_rain, xtprecip_snow REAL, DIMENSION(ntiso,klon), INTENT(IN) :: xtspechum #endif 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 #ifdef ISO REAL, DIMENSION(niso,klon), INTENT(INOUT) :: xtsnow, xtsol REAL, DIMENSION(niso,klon), INTENT(INOUT) :: Rland_ice #endif ! 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 #ifdef ISO REAL, DIMENSION(ntiso,klon), INTENT(OUT) :: xtevap ! REAL, DIMENSION(niso,klon) :: xtrun_off_lic_0_diag ! est une variable globale de ! fonte_neige #endif ! 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 #ifdef ISO REAL, PARAMETER :: t_coup = 273.15 REAL, DIMENSION(klon) :: fqfonte_diag REAL, DIMENSION(klon) :: fq_fonte_diag REAL, DIMENSION(klon) :: snow_evap_diag REAL, DIMENSION(klon) :: fqcalving_diag REAL max_eau_sol_diag REAL, DIMENSION(klon) :: runoff_diag REAL, DIMENSION(klon) :: run_off_lic_diag REAL :: coeff_rel_diag INTEGER :: ixt REAL, DIMENSION(niso,klon) :: xtsnow_prec,xtsol_prec REAL, DIMENSION(klon) :: snow_prec,qsol_prec ! REAL, DIMENSION(klon) :: run_off_lic_0_diag #endif 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 :: time_tempsmooth, coef_tempsmooth REAL, DIMENSION(klon) :: precip_totsnow, evap_totsnow REAL, DIMENSION (klon, 6) :: alb6 REAL :: esalt REAL :: lambdasalt, fluxsalt, csalt, nunu, aa, bb, cc REAL :: tau_dens, maxerosion REAL, DIMENSION(klon) :: ws1, rhod, rhos, ustart0, ustart, qsalt, hsalt REAL, DIMENSION(klon) :: fluxbs_1, fluxbs_2, bsweight_fresh LOGICAL, DIMENSION(klon) :: ok_remaining_freshsnow REAL :: ta1, ta2, ta3, z01, z02, z03, coefa, coefb, coefc, coefd ! 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 !****************************************************************************************** #ifdef ISO #ifdef ISOVERIF ! WRITE(*,*) 'surf_land_ice 1499' DO i=1,knon IF (iso_eau > 0) THEN CALL iso_verif_egalite_choix(xtsnow(iso_eau,i),snow(i), & 'surf_land_ice 126',errmax,errmaxrel) ENDIF !IF (iso_eau > 0) THEN ENDDO !DO i=1,knon #endif #endif 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 DO j=1,knon i = knindex(j) tempsmoothlic(i) = temp_air(j) END DO 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 == 1) THEN !**************************************************************************************** ! CALL to INLANDSIS interface !**************************************************************************************** IF (CPPKEY_INLANDSIS) THEN #ifdef ISO CALL abort_gcm('surf_landice 235','isotopes pas dans INLANDSIS',1) #endif 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) END IF 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) #ifdef ISO #ifdef ISOVERIF !WRITE(*,*) 'surf_land_ice 1499' DO i=1,knon IF (iso_eau > 0) THEN IF (snow(i) > ridicule) THEN CALL iso_verif_egalite_choix(xtsnow(iso_eau,i),snow(i), & 'surf_land_ice 1151',errmax,errmaxrel) ENDIF !IF ((snow(i) > ridicule)) THEN ENDIF !IF (iso_eau > 0) THEN ENDDO !DO i=1,knon #endif DO i=1,knon snow_prec(i)=snow(i) DO ixt=1,niso xtsnow_prec(ixt,i)=xtsnow(ixt,i) ENDDO !DO ixt=1,niso ! initialisation: fq_fonte_diag(i)=0.0 fqfonte_diag(i)=0.0 snow_evap_diag(i)=0.0 ENDDO !DO i=1,knon #endif 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 !**************************************************************************************** !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 !**************************************************************************************** IF (z0m_landice > 0.) THEN z0m(1:knon) = z0m_landice z0h(1:knon) = z0h_landice ELSE ! parameterization of z0=f(T) following measurements in Adelie Land by Amory et al 2018 coefa = 0.1658 !0.1862 !Ant coefb = -50.3869 !-55.7718 !Ant ta1 = 253.15 !255. Ant ta2 = 273.15 ta3 = 273.15 + 3 z01 = exp(coefa * ta1 + coefb) !~0.2 ! ~0.25 mm z02 = exp(coefa * ta2 + coefb) !~6 !~7 mm z03 = z01 coefc = log(z03 / z02) / (ta3 - ta2) coefd = log(z03) - coefc * ta3 time_tempsmooth = 2. * 86400. coef_tempsmooth = min(1., dtime / time_tempsmooth) !coef_tempsmooth=0. DO j = 1, knon i=knindex(j) print*, "tempsmoothlic", tempsmoothlic(i) tempsmoothlic(i) = temp_air(j) * coef_tempsmooth + tempsmoothlic(i) * (1. - coef_tempsmooth) IF (tempsmoothlic(i) < ta1) THEN z0m(j) = z01 ELSE IF (tempsmoothlic(i) >= ta1 .and. tempsmoothlic(i) < ta2) THEN z0m(j) = exp(coefa * tempsmoothlic(i) + coefb) ELSE IF (tempsmoothlic(i) >= ta2 .and. tempsmoothlic(i) < ta3) THEN ! if st > 0, melting induce smooth surface z0m(j) = exp(coefc * tempsmoothlic(i) + coefd) ELSE z0m(j) = z03 END IF z0h(j) = z0m(j) END DO END IF !**************************************************************************************** ! Simple blowing snow param !**************************************************************************************** ! we proceed in 2 steps: ! first we erode - if possible -the accumulated snow during the time step ! then we update the density of the underlying layer and see if we can also erode ! this layer IF (ok_bs) THEN fluxbs(:) = 0. DO j = 1, knon ws1(j) = (u1(j)**2 + v1(j)**2)**0.5 ustar(j) = (cdragm(j) * (u1(j)**2 + v1(j)**2))**0.5 rhod(j) = p1lay(j) / RD / temp_air(j) ustart0(j) = (log(2.868) - log(1.625)) / 0.085 * sqrt(cdragm(j)) enddo ! 1st step: erosion of fresh snow accumulated during the time step DO j = 1, knon IF (precip_snow(j) > 0.) THEN rhos(j) = rhofresh_bs ! blowing snow flux formula used in MAR ustart(j) = ustart0(j) * exp(max(rhoice_bs / rhofresh_bs - rhoice_bs / rhos(j), 0.)) * exp(max(0., rhos(j) - rhohard_bs)) ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs ! computation of qbs at the top of the saltation layer ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001) esalt = 1. / (c_esalt_bs * max(1.e-6, ustar(j))) hsalt(j) = 0.08436 * (max(1.e-6, ustar(j))**1.27) qsalt(j) = (max(ustar(j)**2 - ustart(j)**2, 0.)) / (RG * hsalt(j)) * esalt ! calculation of erosion (flux positive towards the surface here) ! consistent with implicit resolution of turbulent mixing equation ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level) ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer ! (rho*qsalt*hsalt) ! during this first step we also lower bound the erosion to the amount of fresh snow accumulated during the time step maxerosion = min(precip_snow(j), hsalt(j) * qsalt(j) * rhod(j) / tau_eqsalt_bs) fluxbs_1(j) = rhod(j) * ws1(j) * cdragh(j) * zeta_bs * (AcoefQBS(j) - qsalt(j)) & / (1. - rhod(j) * ws1(j) * cdragh(j) * zeta_bs * BcoefQBS(j) * dtime) fluxbs_1(j) = max(-maxerosion, fluxbs_1(j)) IF (precip_snow(j) > abs(fluxbs_1(j))) THEN ok_remaining_freshsnow(j) = .TRUE. bsweight_fresh(j) = 1. else ok_remaining_freshsnow(j) = .FALSE. bsweight_fresh(j) = exp(-(abs(fluxbs_1(j)) - precip_snow(j)) / precip_snow(j)) endif else ok_remaining_freshsnow(j) = .FALSE. fluxbs_1(j) = 0. bsweight_fresh(j) = 0. endif enddo ! we now compute the snow age of the overlying layer (snow surface after erosion of the fresh snow accumulated during the time step) ! this is done through the routine albsno CALL albsno(klon, knon, dtime, agesno(:), alb_neig(:), precip_snow(:) + fluxbs_1(:)) ! 2nd step: ! computation of threshold friction velocity ! which depends on surface snow density DO j = 1, knon IF (ok_remaining_freshsnow(j)) THEN fluxbs_2(j) = 0. else ! we start eroding the underlying layer ! estimation of snow density ! snow density increases with snow age and ! increases even faster in case of sedimentation of blowing snow or rain tau_dens = max(tau_densmin_bs, tau_dens0_bs * exp(-abs(precip_bs(j)) / pbst_bs - & abs(precip_rain(j)) / prt_bs) * exp(-max(tsurf(j) - RTT, 0.))) rhos(j) = rhofresh_bs + (rhohard_bs - rhofresh_bs) * (1. - exp(-agesno(j) * 86400.0 / tau_dens)) ! blowing snow flux formula used in MAR ustart(j) = ustart0(j) * exp(max(rhoice_bs / rhofresh_bs - rhoice_bs / rhos(j), 0.)) * exp(max(0., rhos(j) - rhohard_bs)) ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs ! computation of qbs at the top of the saltation layer ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001) esalt = 1. / (c_esalt_bs * max(1.e-6, ustar(j))) hsalt(j) = 0.08436 * (max(1.e-6, ustar(j))**1.27) qsalt(j) = (max(ustar(j)**2 - ustart(j)**2, 0.)) / (RG * hsalt(j)) * esalt ! calculation of erosion (flux positive towards the surface here) ! consistent with implicit resolution of turbulent mixing equation ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level) ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer ! (rho*qsalt*hsalt) maxerosion = hsalt(j) * qsalt(j) * rhod(j) / tau_eqsalt_bs fluxbs_2(j) = rhod(j) * ws1(j) * cdragh(j) * zeta_bs * (AcoefQBS(j) - qsalt(j)) & / (1. - rhod(j) * ws1(j) * cdragh(j) * zeta_bs * BcoefQBS(j) * dtime) fluxbs_2(j) = max(-maxerosion, fluxbs_2(j)) endif enddo ! final flux and outputs DO j = 1, knon ! total flux is the erosion of fresh snow + ! a fraction of the underlying snow (if all the fresh snow has been eroded) ! the calculation of the fraction is quite delicate since we do not know ! how much time was needed to erode the fresh snow. We assume that this time ! is dt*exp(-(abs(fluxbs1)-precipsnow)/precipsnow)=dt*bsweight_fresh fluxbs(j) = fluxbs_1(j) + fluxbs_2(j) * (1. - bsweight_fresh(j)) i = knindex(j) zxustartlic(i) = ustart(j) zxrhoslic(i) = rhos(j) zxqsaltlic(i) = qsalt(j) enddo else ! not ok_bs ! those lines are useful to calculate the snow age CALL albsno(klon, knon, dtime, agesno(:), alb_neig(:), precip_snow(:)) endif ! if ok_bs !**************************************************************************************** ! Calculate 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 & #ifdef ISO ,fq_fonte_diag,fqfonte_diag,snow_evap_diag,fqcalving_diag & ,max_eau_sol_diag,runoff_diag,run_off_lic_diag,coeff_rel_diag & #endif ) #ifdef ISO #ifdef ISOVERIF DO i=1,knon IF (iso_eau > 0) THEN CALL iso_verif_egalite_choix(Rland_ice(iso_eau,i),1.0, & 'surf_landice_mod 217',errmax,errmaxrel) ENDIF !IF (iso_eau > 0) THEN ENDDO !DO i=1,knon #endif CALL calcul_iso_surf_lic_vectall(klon,knon, & evap,snow_evap_diag,Tsurf_new,snow, & fq_fonte_diag,fqfonte_diag,dtime,t_coup, & precip_snow,xtprecip_snow,precip_rain,xtprecip_rain, snow_prec,xtsnow_prec, & xtspechum,spechum,ps,Rland_ice, & xtevap,xtsnow,fqcalving_diag, & knindex,is_lic,run_off_lic_diag,coeff_rel_diag & ) ! CALL fonte_neige_export_xtrun_off_lic_0(knon,xtrun_off_lic_0_diag) #endif WHERE (snow(1:knon) < 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 < 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 == 1) .AND. (iflag_albcalc == 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