source: LMDZ6/trunk/libf/phylmd/surf_landice_mod.F90 @ 5022

Last change on this file since 5022 was 5022, checked in by Sebastien Nguyen, 2 months ago

include ISO keys in pbl_surface and associated routines in phylmd

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 27.8 KB
Line 
1!
2MODULE surf_landice_mod
3 
4  IMPLICIT NONE
5
6CONTAINS
7!
8!****************************************************************************************
9!
10  SUBROUTINE surf_landice(itime, dtime, knon, knindex, &
11       rlon, rlat, debut, lafin, &
12       rmu0, lwdownm, albedo, pphi1, &
13       swnet, lwnet, tsurf, p1lay, &
14       cdragh, cdragm, precip_rain, precip_snow, precip_bs, temp_air, spechum, &
15       AcoefH, AcoefQ, BcoefH, BcoefQ, &
16       AcoefU, AcoefV, BcoefU, BcoefV, &
17       AcoefQBS, BcoefQBS, &
18       ps, u1, v1, gustiness, rugoro, pctsrf, &
19       snow, qsurf, qsol, qbs1, agesno, &
20       tsoil, z0m, z0h, SFRWL, alb_dir, alb_dif, evap, fluxsens, fluxlat, fluxbs, &
21       tsurf_new, dflux_s, dflux_l, &
22       alt, slope, cloudf, &
23       snowhgt, qsnow, to_ice, sissnow, &
24       alb3, runoff, &
25       flux_u1, flux_v1 &
26#ifdef ISO
27         &      ,xtprecip_rain, xtprecip_snow,xtspechum,Rland_ice &
28         &      ,xtsnow,xtsol,xtevap &
29#endif               
30           &    )
31
32    USE dimphy
33    USE geometry_mod,     ONLY : longitude,latitude
34    USE surface_data,     ONLY : type_ocean, calice, calsno, landice_opt, iflag_albcalc
35    USE fonte_neige_mod,  ONLY : fonte_neige,run_off_lic,fqcalving_global,ffonte_global,fqfonte_global,runofflic_global
36    USE cpl_mod,          ONLY : cpl_send_landice_fields
37    USE calcul_fluxs_mod
38    USE phys_local_var_mod, ONLY : zxrhoslic, zxustartlic, zxqsaltlic
39    USE phys_output_var_mod, ONLY : snow_o,zfra_o
40#ifdef ISO   
41    USE fonte_neige_mod,  ONLY : xtrun_off_lic
42    USE infotrac_phy,     ONLY : ntiso,niso
43    USE isotopes_routines_mod, ONLY: calcul_iso_surf_lic_vectall
44#ifdef ISOVERIF
45    USE isotopes_mod, ONLY: iso_eau,ridicule
46    USE isotopes_verif_mod
47#endif
48#endif
49 
50!FC
51    USE ioipsl_getin_p_mod, ONLY : getin_p
52    USE lmdz_blowing_snow_ini, ONLY : c_esalt_bs, zeta_bs, pbst_bs, prt_bs, rhoice_bs, rhohard_bs
53    USE lmdz_blowing_snow_ini, ONLY : rhofresh_bs, tau_eqsalt_bs, tau_dens0_bs, tau_densmin_bs
54#ifdef CPP_INLANDSIS
55    USE surf_inlandsis_mod,  ONLY : surf_inlandsis
56#endif
57
58    USE indice_sol_mod
59
60!    INCLUDE "indicesol.h"
61    INCLUDE "dimsoil.h"
62    INCLUDE "YOMCST.h"
63    INCLUDE "clesphys.h"
64
65! Input variables
66!****************************************************************************************
67    INTEGER, INTENT(IN)                           :: itime, knon
68    INTEGER, DIMENSION(klon), INTENT(in)          :: knindex
69    REAL, INTENT(in)                              :: dtime
70    REAL, DIMENSION(klon), INTENT(IN)             :: swnet ! net shortwave radiance
71    REAL, DIMENSION(klon), INTENT(IN)             :: lwnet ! net longwave radiance
72    REAL, DIMENSION(klon), INTENT(IN)             :: tsurf
73    REAL, DIMENSION(klon), INTENT(IN)             :: p1lay
74    REAL, DIMENSION(klon), INTENT(IN)             :: cdragh, cdragm
75    REAL, DIMENSION(klon), INTENT(IN)             :: precip_rain, precip_snow, precip_bs
76    REAL, DIMENSION(klon), INTENT(IN)             :: temp_air, spechum
77    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefH, AcoefQ
78    REAL, DIMENSION(klon), INTENT(IN)             :: BcoefH, BcoefQ
79    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefU, AcoefV, BcoefU, BcoefV
80    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefQBS, BcoefQBS
81    REAL, DIMENSION(klon), INTENT(IN)             :: ps
82    REAL, DIMENSION(klon), INTENT(IN)             :: u1, v1, gustiness, qbs1
83    REAL, DIMENSION(klon), INTENT(IN)             :: rugoro
84    REAL, DIMENSION(klon,nbsrf), INTENT(IN)       :: pctsrf
85#ifdef ISO
86    REAL, DIMENSION(ntiso,klon), INTENT(IN)       :: xtprecip_rain, xtprecip_snow
87    REAL, DIMENSION(ntiso,klon), INTENT(IN)       :: xtspechum
88#endif
89
90
91    LOGICAL,  INTENT(IN)                          :: debut   !true if first step
92    LOGICAL,  INTENT(IN)                          :: lafin   !true if last step
93    REAL, DIMENSION(klon), INTENT(IN)             :: rlon, rlat
94    REAL, DIMENSION(klon), INTENT(IN)             :: rmu0
95    REAL, DIMENSION(klon), INTENT(IN)             :: lwdownm !ylwdown
96    REAL, DIMENSION(klon), INTENT(IN)             :: albedo  !mean albedo
97    REAL, DIMENSION(klon), INTENT(IN)             :: pphi1   
98    REAL, DIMENSION(klon), INTENT(IN)             :: alt   !mean altitude of the grid box 
99    REAL, DIMENSION(klon), INTENT(IN)             :: slope   !mean slope in grid box 
100    REAL, DIMENSION(klon), INTENT(IN)             :: cloudf  !total cloud fraction
101
102! In/Output variables
103!****************************************************************************************
104    REAL, DIMENSION(klon), INTENT(INOUT)          :: snow, qsol
105    REAL, DIMENSION(klon), INTENT(INOUT)          :: agesno
106    REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil
107#ifdef ISO
108    REAL, DIMENSION(niso,klon), INTENT(INOUT)     :: xtsnow, xtsol
109    REAL, DIMENSION(niso,klon), INTENT(INOUT)     :: Rland_ice
110#endif
111
112
113! Output variables
114!****************************************************************************************
115    REAL, DIMENSION(klon), INTENT(OUT)            :: qsurf
116    REAL, DIMENSION(klon), INTENT(OUT)            :: z0m, z0h
117!albedo SB >>>
118!    REAL, DIMENSION(klon), INTENT(OUT)            :: alb1  ! new albedo in visible SW interval
119!    REAL, DIMENSION(klon), INTENT(OUT)            :: alb2  ! new albedo in near IR interval
120    REAL, DIMENSION(6), INTENT(IN)                :: SFRWL
121    REAL, DIMENSION(klon,nsw), INTENT(OUT)        :: alb_dir,alb_dif
122!albedo SB <<<
123    REAL, DIMENSION(klon), INTENT(OUT)            :: evap, fluxsens, fluxlat
124    REAL, DIMENSION(klon), INTENT(OUT)            :: fluxbs
125    REAL, DIMENSION(klon), INTENT(OUT)            :: tsurf_new
126    REAL, DIMENSION(klon), INTENT(OUT)            :: dflux_s, dflux_l     
127    REAL, DIMENSION(klon), INTENT(OUT)            :: flux_u1, flux_v1
128
129    REAL, DIMENSION(klon), INTENT(OUT)           :: alb3
130    REAL, DIMENSION(klon), INTENT(OUT)           :: qsnow   !column water in snow [kg/m2]
131    REAL, DIMENSION(klon), INTENT(OUT)           :: snowhgt !Snow height (m)
132    REAL, DIMENSION(klon), INTENT(OUT)           :: to_ice
133    REAL, DIMENSION(klon), INTENT(OUT)           :: sissnow
134    REAL, DIMENSION(klon), INTENT(OUT)           :: runoff  !Land ice runoff
135#ifdef ISO
136    REAL, DIMENSION(ntiso,klon), INTENT(OUT)     :: xtevap     
137!    real, DIMENSION(niso,klon) :: xtrun_off_lic_0_diag ! est une variable globale de
138!    fonte_neige
139#endif
140 
141
142! Local variables
143!****************************************************************************************
144    REAL, DIMENSION(klon)    :: soilcap, soilflux
145    REAL, DIMENSION(klon)    :: cal, beta, dif_grnd
146    REAL, DIMENSION(klon)    :: zfra, alb_neig
147    REAL, DIMENSION(klon)    :: radsol
148    REAL, DIMENSION(klon)    :: u0, v0, u1_lay, v1_lay, ustar
149    INTEGER                  :: i,j,nt
150    REAL, DIMENSION(klon)    :: fqfonte,ffonte
151    REAL, DIMENSION(klon)    :: run_off_lic_frac
152#ifdef ISO       
153    REAL, PARAMETER          :: t_coup = 273.15
154    REAL, DIMENSION(klon)    :: fqfonte_diag
155    REAL, DIMENSION(klon)    :: fq_fonte_diag
156    REAL, DIMENSION(klon)    ::  snow_evap_diag
157    REAL, DIMENSION(klon)    ::  fqcalving_diag
158    REAL max_eau_sol_diag 
159    REAL, DIMENSION(klon)    ::  runoff_diag
160    REAL, DIMENSION(klon)    ::    run_off_lic_diag
161    REAL                     ::  coeff_rel_diag
162    INTEGER                  :: ixt
163    REAL, DIMENSION(niso,klon) :: xtsnow_prec,xtsol_prec
164    REAL, DIMENSION(klon) :: snow_prec,qsol_prec
165!    real, DIMENSION(klon) :: run_off_lic_0_diag
166#endif
167
168
169    REAL, DIMENSION(klon)    :: emis_new                  !Emissivity
170    REAL, DIMENSION(klon)    :: swdown,lwdown
171    REAL, DIMENSION(klon)    :: precip_snow_adv, snow_adv !Snow Drift precip./advection (not used in inlandsis)
172    REAL, DIMENSION(klon)    :: erod                      !erosion of surface snow (flux, kg/m2/s like evap)
173    REAL, DIMENSION(klon)    :: zsl_height, wind_velo     !surface layer height, wind spd
174    REAL, DIMENSION(klon)    :: dens_air,  snow_cont_air  !air density; snow content air
175    REAL, DIMENSION(klon)    :: alb_soil                  !albedo of underlying ice
176    REAL, DIMENSION(klon)    :: pexner                    !Exner potential
177    REAL                     :: pref
178    REAL, DIMENSION(klon,nsoilmx) :: tsoil0               !modif
179    REAL                          :: dtis                ! subtimestep
180    LOGICAL                       :: debut_is, lafin_is  ! debut and lafin for inlandsis
181
182    CHARACTER (len = 20)                      :: modname = 'surf_landice'
183    CHARACTER (len = 80)                      :: abort_message
184
185
186    REAL,DIMENSION(klon) :: alb1,alb2
187    REAL,DIMENSION(klon) :: precip_totsnow, evap_totsnow
188    REAL, DIMENSION (klon,6) :: alb6
189    REAL                   :: esalt
190    REAL                   :: lambdasalt,fluxsalt, csalt, nunu, aa, bb, cc
191    REAL                   :: tau_dens, maxerosion
192    REAL, DIMENSION(klon)  :: ws1, rhod, rhos, ustart0, ustart, qsalt, hsalt
193    REAL, DIMENSION(klon)  :: fluxbs_1, fluxbs_2, bsweight_fresh
194    LOGICAL, DIMENSION(klon) :: ok_remaining_freshsnow
195    REAL  :: ta1, ta2, ta3, z01, z02, z03, coefa, coefb, coefc, coefd
196
197
198! End definition
199!****************************************************************************************
200!FC
201!FC
202   REAL,SAVE :: alb_vis_sno_lic
203  !$OMP THREADPRIVATE(alb_vis_sno_lic)
204   REAL,SAVE :: alb_nir_sno_lic
205  !$OMP THREADPRIVATE(alb_nir_sno_lic)
206  LOGICAL, SAVE :: firstcall = .TRUE.
207  !$OMP THREADPRIVATE(firstcall)
208
209
210!FC firtscall initializations
211!******************************************************************************************
212#ifdef ISO
213#ifdef ISOVERIF
214!     write(*,*) 'surf_land_ice 1499'   
215  DO i=1,knon
216    IF (iso_eau > 0) THEN
217      CALL iso_verif_egalite_choix(xtsnow(iso_eau,i),snow(i), &
218    &                              'surf_land_ice 126',errmax,errmaxrel)
219    ENDIF !IF (iso_eau > 0) THEN     
220  ENDDO !DO i=1,knon 
221#endif
222#endif
223
224  IF (firstcall) THEN
225  alb_vis_sno_lic=0.77
226  CALL getin_p('alb_vis_sno_lic',alb_vis_sno_lic)
227           PRINT*, 'alb_vis_sno_lic',alb_vis_sno_lic
228  alb_nir_sno_lic=0.77
229  CALL getin_p('alb_nir_sno_lic',alb_nir_sno_lic)
230           PRINT*, 'alb_nir_sno_lic',alb_nir_sno_lic
231 
232  firstcall=.false.
233  ENDIF
234!******************************************************************************************
235
236! Initialize output variables
237    alb3(:) = 999999.
238    alb2(:) = 999999.
239    alb1(:) = 999999.
240    fluxbs(:)=0. 
241    runoff(:) = 0.
242!****************************************************************************************
243! Calculate total absorbed radiance at surface
244!
245!****************************************************************************************
246    radsol(:) = 0.0
247    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
248
249!****************************************************************************************
250
251!****************************************************************************************
252!  landice_opt = 0 : soil_model, calcul_flux, fonte_neige, ... 
253!  landice_opt = 1  : prepare and call INterace Lmdz SISvat (INLANDSIS)
254!****************************************************************************************
255
256
257    IF (landice_opt .EQ. 1) THEN
258
259!****************************************************************************************   
260! CALL to INLANDSIS interface
261!****************************************************************************************
262#ifdef CPP_INLANDSIS
263
264#ifdef ISO
265        CALL abort_gcm('surf_landice 235','isotopes pas dans INLANDSIS',1)
266#endif
267
268        debut_is=debut
269        lafin_is=.false.
270        ! Suppose zero surface speed
271        u0(:)            = 0.0
272        v0(:)            = 0.0
273
274
275        CALL calcul_flux_wind(knon, dtime, &
276         u0, v0, u1, v1, gustiness, cdragm, &
277         AcoefU, AcoefV, BcoefU, BcoefV, &
278         p1lay, temp_air, &
279         flux_u1, flux_v1)
280
281       
282       ! Set constants and compute some input for SISVAT
283       ! = 1000 hPa
284       ! and calculate incoming flux for SW and LW interval: swdown, lwdown
285       swdown(:)        = 0.0
286       lwdown(:)        = 0.0
287       snow_cont_air(:) = 0.  ! the snow content in air is not a prognostic variable of the model     
288       alb_soil(:)      = 0.4 ! before albedo(:) but here it is the ice albedo that we have to set
289       ustar(:)         = 0.
290       pref             = 100000.       
291       DO i = 1, knon
292          swdown(i)        = swnet(i)/(1-albedo(i))
293          lwdown(i)        = lwdownm(i)
294          wind_velo(i)     = u1(i)**2 + v1(i)**2
295          wind_velo(i)     = wind_velo(i)**0.5
296          pexner(i)        = (p1lay(i)/pref)**(RD/RCPD)
297          dens_air(i)      = p1lay(i)/RD/temp_air(i)  ! dry air density
298          zsl_height(i)    = pphi1(i)/RG     
299          tsoil0(i,:)      = tsoil(i,:) 
300          ustar(i)= (cdragm(i)*(wind_velo(i)**2))**0.5   
301       END DO
302       
303
304
305        dtis=dtime
306
307          IF (lafin) THEN
308            lafin_is=.true.
309          END IF
310
311          CALL surf_inlandsis(knon, rlon, rlat, knindex, itime, dtis, debut_is, lafin_is,&
312            rmu0, swdown, lwdown, albedo, pexner, ps, p1lay, precip_rain, precip_snow,   &
313            zsl_height, wind_velo, ustar, temp_air, dens_air, spechum, tsurf,&
314            rugoro, snow_cont_air, alb_soil, alt, slope, cloudf, &
315            radsol, qsol, tsoil0, snow, zfra, snowhgt, qsnow, to_ice, sissnow,agesno,   &
316            AcoefH, AcoefQ, BcoefH, BcoefQ, cdragm, cdragh, &
317            run_off_lic, fqfonte, ffonte, evap, erod, fluxsens, fluxlat,dflux_s, dflux_l, &
318            tsurf_new, alb1, alb2, alb3, alb6, &
319            emis_new, z0m, z0h, qsurf)
320
321          debut_is=.false.
322
323
324        ! Treatment of snow melting and calving
325
326        ! for consistency with standard LMDZ, add calving to run_off_lic
327        run_off_lic(:)=run_off_lic(:) + to_ice(:)
328
329        DO i = 1, knon
330           ffonte_global(knindex(i),is_lic)    = ffonte(i)
331           fqfonte_global(knindex(i),is_lic)   = fqfonte(i)! net melting= melting - refreezing
332           fqcalving_global(knindex(i),is_lic) = to_ice(i) ! flux
333           runofflic_global(knindex(i)) = run_off_lic(i)
334        ENDDO
335        ! Here, we assume that the calving term is equal to the to_ice term
336        ! (no ice accumulation)
337
338
339#else
340       abort_message='Pb de coherence: landice_opt = 1  mais CPP_INLANDSIS = .false.'
341       CALL abort_physic(modname,abort_message,1)
342#endif
343
344
345    ELSE
346
347!****************************************************************************************
348! Soil calculations
349!
350!****************************************************************************************
351
352    ! EV: use calbeta
353    CALL calbeta(dtime, is_lic, knon, snow, qsol, beta, cal, dif_grnd)
354
355
356    ! use soil model and recalculate properly cal
357    IF (soil_model) THEN
358       CALL soil(dtime, is_lic, knon, snow, tsurf, qsol, &
359        & longitude(knindex(1:knon)), latitude(knindex(1:knon)), tsoil, soilcap, soilflux)
360       cal(1:knon) = RCPD / soilcap(1:knon)
361       radsol(1:knon)  = radsol(1:knon) + soilflux(1:knon)
362    ELSE
363       cal = RCPD * calice
364       WHERE (snow > 0.0) cal = RCPD * calsno
365    ENDIF
366
367
368!****************************************************************************************
369! Calulate fluxes
370!
371!****************************************************************************************
372!    beta(:) = 1.0
373!    dif_grnd(:) = 0.0
374
375! Suppose zero surface speed
376    u0(:)=0.0
377    v0(:)=0.0
378    u1_lay(:) = u1(:) - u0(:)
379    v1_lay(:) = v1(:) - v0(:)
380
381    CALL calcul_fluxs(knon, is_lic, dtime, &
382         tsurf, p1lay, cal, beta, cdragh, cdragh, ps, &
383         precip_rain, precip_snow, snow, qsurf,  &
384         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
385         1.,AcoefH, AcoefQ, BcoefH, BcoefQ, &
386         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
387
388#ifdef ISO
389#ifdef ISOVERIF
390     !write(*,*) 'surf_land_ice 1499'   
391     DO i=1,knon
392       IF (iso_eau > 0) THEN
393         IF (snow(i) > ridicule) THEN
394           CALL iso_verif_egalite_choix(xtsnow(iso_eau,i),snow(i), &
395    &                                   'surf_land_ice 1151',errmax,errmaxrel)
396         ENDIF !IF ((snow(i) > ridicule)) THEN
397       ENDIF !IF (iso_eau > 0) THEN
398     ENDDO !DO i=1,knon 
399#endif
400
401    DO i=1,knon
402      snow_prec(i)=snow(i)
403      DO ixt=1,niso
404        xtsnow_prec(ixt,i)=xtsnow(ixt,i)
405      ENDDO !DO ixt=1,niso
406      ! initialisation:
407      fq_fonte_diag(i)=0.0
408      fqfonte_diag(i)=0.0
409      snow_evap_diag(i)=0.0
410    ENDDO !DO i=1,knon
411#endif         
412
413    CALL calcul_flux_wind(knon, dtime, &
414         u0, v0, u1, v1, gustiness, cdragm, &
415         AcoefU, AcoefV, BcoefU, BcoefV, &
416         p1lay, temp_air, &
417         flux_u1, flux_v1)
418
419
420!****************************************************************************************
421! Calculate albedo
422!
423!****************************************************************************************
424
425!
426!IM: plusieurs choix/tests sur l'albedo des "glaciers continentaux"
427!       alb1(1 : knon)  = 0.6 !IM cf FH/GK
428!       alb1(1 : knon)  = 0.82
429!       alb1(1 : knon)  = 0.77 !211003 Ksta0.77
430!       alb1(1 : knon)  = 0.8 !KstaTER0.8 & LMD_ARMIP5
431!IM: KstaTER0.77 & LMD_ARMIP6   
432
433! Attantion: alb1 and alb2 are not the same!
434    alb1(1:knon)  = alb_vis_sno_lic
435    alb2(1:knon)  = alb_nir_sno_lic
436
437
438!****************************************************************************************
439! Rugosity
440!
441!****************************************************************************************
442
443if (z0m_landice .GT. 0.) then
444    z0m(1:knon) = z0m_landice
445    z0h(1:knon) = z0h_landice
446else
447    ! parameterization of z0=f(T) following measurements in Adelie Land by Amory et al 2018
448    coefa = 0.1658 !0.1862 !Ant
449    coefb = -50.3869 !-55.7718 !Ant
450    ta1 = 253.15 !255. Ant
451    ta2 = 273.15
452    ta3 = 273.15+3
453    z01 = exp(coefa*ta1 + coefb) !~0.2 ! ~0.25 mm
454    z02 = exp(coefa*ta2 + coefb) !~6  !~7 mm
455    z03 = z01
456    coefc = log(z03/z02)/(ta3-ta2)
457    coefd = log(z03)-coefc*ta3
458    do j=1,knon
459      if (temp_air(j) .lt. ta1) then
460        z0m(j) = z01
461      else if (temp_air(j).ge.ta1 .and. temp_air(j).lt.ta2) then
462        z0m(j) = exp(coefa*temp_air(j) + coefb)
463      else if (temp_air(j).ge.ta2 .and. temp_air(j).lt.ta3) then
464        ! if st > 0, melting induce smooth surface
465        z0m(j) = exp(coefc*temp_air(j) + coefd)
466      else
467        z0m(j) = z03
468      endif
469      z0h(j)=z0m(j)
470    enddo
471
472endif   
473 
474
475!****************************************************************************************
476! Simple blowing snow param
477!****************************************************************************************
478! we proceed in 2 steps:
479! first we erode - if possible -the accumulated snow during the time step
480! then we update the density of the underlying layer and see if we can also erode
481! this layer
482
483
484   if (ok_bs) then
485       fluxbs(:)=0.
486       do j=1,knon
487          ws1(j)=(u1(j)**2+v1(j)**2)**0.5
488          ustar(j)=(cdragm(j)*(u1(j)**2+v1(j)**2))**0.5
489          rhod(j)=p1lay(j)/RD/temp_air(j)
490          ustart0(j) =(log(2.868)-log(1.625))/0.085*sqrt(cdragm(j))
491       enddo
492
493       ! 1st step: erosion of fresh snow accumulated during the time step
494       do j=1, knon
495       if (precip_snow(j) .GT. 0.) then
496           rhos(j)=rhofresh_bs
497           ! blowing snow flux formula used in MAR
498           ustart(j)=ustart0(j)*exp(max(rhoice_bs/rhofresh_bs-rhoice_bs/rhos(j),0.))*exp(max(0.,rhos(j)-rhohard_bs))
499           ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs
500           ! computation of qbs at the top of the saltation layer
501           ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001)
502           esalt=1./(c_esalt_bs*max(1.e-6,ustar(j)))
503           hsalt(j)=0.08436*(max(1.e-6,ustar(j))**1.27)
504           qsalt(j)=(max(ustar(j)**2-ustart(j)**2,0.))/(RG*hsalt(j))*esalt
505           ! calculation of erosion (flux positive towards the surface here)
506           ! consistent with implicit resolution of turbulent mixing equation
507           ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s
508           ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level)
509           ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer
510           ! (rho*qsalt*hsalt)
511           ! during this first step we also lower bound the erosion to the amount of fresh snow accumulated during the time step
512           maxerosion=min(precip_snow(j),hsalt(j)*qsalt(j)*rhod(j)/tau_eqsalt_bs)
513
514           fluxbs_1(j)=rhod(j)*ws1(j)*cdragh(j)*zeta_bs*(AcoefQBS(j)-qsalt(j)) &
515                   / (1.-rhod(j)*ws1(j)*cdragh(j)*zeta_bs*BcoefQBS(j)*dtime)
516           fluxbs_1(j)=max(-maxerosion,fluxbs_1(j))
517
518           if (precip_snow(j) .gt. abs(fluxbs_1(j))) then
519               ok_remaining_freshsnow(j)=.true.
520               bsweight_fresh(j)=1.
521           else
522               ok_remaining_freshsnow(j)=.false.
523               bsweight_fresh(j)=exp(-(abs(fluxbs_1(j))-precip_snow(j))/precip_snow(j))
524           endif
525       else
526           ok_remaining_freshsnow(j)=.false.
527           fluxbs_1(j)=0.
528           bsweight_fresh(j)=0.
529       endif
530       enddo
531
532
533       ! we now compute the snow age of the overlying layer (snow surface after erosion of the fresh snow accumulated during the time step)
534       ! this is done through the routine albsno
535       CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)+fluxbs_1(:))
536
537       ! 2nd step:
538       ! computation of threshold friction velocity
539       ! which depends on surface snow density
540       do j = 1, knon
541        if (ok_remaining_freshsnow(j)) then
542           fluxbs_2(j)=0.
543        else
544           ! we start eroding the underlying layer
545           ! estimation of snow density
546           ! snow density increases with snow age and
547           ! increases even faster in case of sedimentation of blowing snow or rain
548           tau_dens=max(tau_densmin_bs, tau_dens0_bs*exp(-abs(precip_bs(j))/pbst_bs - &
549                    abs(precip_rain(j))/prt_bs)*exp(-max(tsurf(j)-RTT,0.)))
550           rhos(j)=rhofresh_bs+(rhohard_bs-rhofresh_bs)*(1.-exp(-agesno(j)*86400.0/tau_dens))
551           ! blowing snow flux formula used in MAR
552           ustart(j)=ustart0(j)*exp(max(rhoice_bs/rhofresh_bs-rhoice_bs/rhos(j),0.))*exp(max(0.,rhos(j)-rhohard_bs))
553           ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs
554           ! computation of qbs at the top of the saltation layer
555           ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001)
556           esalt=1./(c_esalt_bs*max(1.e-6,ustar(j)))
557           hsalt(j)=0.08436*(max(1.e-6,ustar(j))**1.27)
558           qsalt(j)=(max(ustar(j)**2-ustart(j)**2,0.))/(RG*hsalt(j))*esalt
559           ! calculation of erosion (flux positive towards the surface here)
560           ! consistent with implicit resolution of turbulent mixing equation
561           ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s
562           ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level)
563           ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer
564           ! (rho*qsalt*hsalt)
565           maxerosion=hsalt(j)*qsalt(j)*rhod(j)/tau_eqsalt_bs
566           fluxbs_2(j)=rhod(j)*ws1(j)*cdragh(j)*zeta_bs*(AcoefQBS(j)-qsalt(j)) &
567                   / (1.-rhod(j)*ws1(j)*cdragh(j)*zeta_bs*BcoefQBS(j)*dtime)
568           fluxbs_2(j)=max(-maxerosion,fluxbs_2(j))
569         endif
570       enddo
571
572
573
574
575       ! final flux and outputs       
576        do j=1, knon
577              ! total flux is the erosion of fresh snow +
578              ! a fraction of the underlying snow (if all the fresh snow has been eroded)
579              ! the calculation of the fraction is quite delicate since we do not know
580              ! how much time was needed to erode the fresh snow. We assume that this time
581              ! is dt*exp(-(abs(fluxbs1)-precipsnow)/precipsnow)=dt*bsweight_fresh
582
583              fluxbs(j)=fluxbs_1(j)+fluxbs_2(j)*(1.-bsweight_fresh(j))
584              i = knindex(j)
585              zxustartlic(i) = ustart(j)
586              zxrhoslic(i) = rhos(j)
587              zxqsaltlic(i)=qsalt(j)
588        enddo
589
590
591  else ! not ok_bs
592  ! those lines are useful to calculate the snow age
593       CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:))
594
595  endif ! if ok_bs
596
597
598
599!****************************************************************************************
600! Calculate snow amount
601!   
602!****************************************************************************************
603    IF (ok_bs) THEN
604      precip_totsnow(:)=precip_snow(:)+precip_bs(:)
605      evap_totsnow(:)=evap(:)-fluxbs(:) ! flux bs is positive towards the surface (snow erosion)
606    ELSE
607      precip_totsnow(:)=precip_snow(:)
608      evap_totsnow(:)=evap(:)
609    ENDIF
610   
611 
612    CALL fonte_neige(knon, is_lic, knindex, dtime, &
613         tsurf, precip_rain, precip_totsnow, &
614         snow, qsol, tsurf_new, evap_totsnow &
615#ifdef ISO   
616     &  ,fq_fonte_diag,fqfonte_diag,snow_evap_diag,fqcalving_diag     &
617     &  ,max_eau_sol_diag,runoff_diag,run_off_lic_diag,coeff_rel_diag &
618#endif
619     &   )
620
621
622#ifdef ISO
623#ifdef ISOVERIF
624    DO i=1,knon 
625      IF (iso_eau > 0) THEN 
626        CALL iso_verif_egalite_choix(Rland_ice(iso_eau,i),1.0, &
627     &                               'surf_landice_mod 217',errmax,errmaxrel)
628      ENDIF !IF (iso_eau > 0) THEN
629    ENDDO !DO i=1,knon
630#endif
631
632    CALL calcul_iso_surf_lic_vectall(klon,knon, &
633     &    evap,snow_evap_diag,Tsurf_new,snow, &
634     &    fq_fonte_diag,fqfonte_diag,dtime,t_coup, &
635     &    precip_snow,xtprecip_snow,precip_rain,xtprecip_rain, snow_prec,xtsnow_prec, &
636     &    xtspechum,spechum,ps,Rland_ice, &
637     &    xtevap,xtsnow,fqcalving_diag, &
638     &    knindex,is_lic,run_off_lic_diag,coeff_rel_diag &
639     &   )
640
641!        call fonte_neige_export_xtrun_off_lic_0(knon,xtrun_off_lic_0_diag)
642
643#endif
644   
645    WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.                                         
646    zfra(1:knon) = MAX(0.0,MIN(1.0,snow(1:knon)/(snow(1:knon)+10.0))) 
647
648
649    END IF ! landice_opt
650
651
652!****************************************************************************************
653! Send run-off on land-ice to coupler if coupled ocean.
654! run_off_lic has been calculated in fonte_neige or surf_inlandsis
655! If landice_opt>=2, corresponding call is done from surf_land_orchidee
656!****************************************************************************************
657    IF (type_ocean=='couple' .AND. landice_opt .LT. 2) THEN
658       ! Compress fraction where run_off_lic is active (here all pctsrf(is_lic))
659       run_off_lic_frac(:)=0.0
660       DO j = 1, knon
661          i = knindex(j)
662          run_off_lic_frac(j) = pctsrf(i,is_lic)
663       ENDDO
664
665       CALL cpl_send_landice_fields(itime, knon, knindex, run_off_lic, run_off_lic_frac)
666    ENDIF
667
668 ! transfer runoff rate [kg/m2/s](!) to physiq for output
669    runoff(1:knon)=run_off_lic(1:knon)/dtime
670
671       snow_o=0.
672       zfra_o = 0.
673       DO j = 1, knon
674           i = knindex(j)
675           snow_o(i) = snow(j)
676           zfra_o(i) = zfra(j)
677       ENDDO
678
679
680!albedo SB >>>
681     select case(NSW)
682     case(2)
683       alb_dir(1:knon,1)=alb1(1:knon)
684       alb_dir(1:knon,2)=alb2(1:knon)
685     case(4)
686       alb_dir(1:knon,1)=alb1(1:knon)
687       alb_dir(1:knon,2)=alb2(1:knon)
688       alb_dir(1:knon,3)=alb2(1:knon)
689       alb_dir(1:knon,4)=alb2(1:knon)
690     case(6)
691       alb_dir(1:knon,1)=alb1(1:knon)
692       alb_dir(1:knon,2)=alb1(1:knon)
693       alb_dir(1:knon,3)=alb1(1:knon)
694       alb_dir(1:knon,4)=alb2(1:knon)
695       alb_dir(1:knon,5)=alb2(1:knon)
696       alb_dir(1:knon,6)=alb2(1:knon)
697
698       IF ((landice_opt .EQ. 1) .AND. (iflag_albcalc .EQ. 2)) THEN
699       alb_dir(1:knon,1)=alb6(1:knon,1)
700       alb_dir(1:knon,2)=alb6(1:knon,2)
701       alb_dir(1:knon,3)=alb6(1:knon,3)
702       alb_dir(1:knon,4)=alb6(1:knon,4)
703       alb_dir(1:knon,5)=alb6(1:knon,5)
704       alb_dir(1:knon,6)=alb6(1:knon,6)
705       ENDIF
706
707     end select
708alb_dif=alb_dir
709!albedo SB <<<
710
711
712  END SUBROUTINE surf_landice
713!
714!****************************************************************************************
715!
716END MODULE surf_landice_mod
717
718
719
Note: See TracBrowser for help on using the repository browser.