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

Last change on this file since 4835 was 4835, checked in by evignon, 3 months ago

commission pour la suite du travail sur la mise a jour
de la param de neige soufflee

  • 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: 22.3 KB
RevLine 
[781]1!
2MODULE surf_landice_mod
3 
4  IMPLICIT NONE
5
6CONTAINS
7!
8!****************************************************************************************
9!
10  SUBROUTINE surf_landice(itime, dtime, knon, knindex, &
[1865]11       rlon, rlat, debut, lafin, &
12       rmu0, lwdownm, albedo, pphi1, &
[888]13       swnet, lwnet, tsurf, p1lay, &
[4523]14       cdragh, cdragm, precip_rain, precip_snow, precip_bs, temp_air, spechum, &
[1067]15       AcoefH, AcoefQ, BcoefH, BcoefQ, &
16       AcoefU, AcoefV, BcoefU, BcoefV, &
[4529]17       AcoefQBS, BcoefQBS, &
[2240]18       ps, u1, v1, gustiness, rugoro, pctsrf, &
[4523]19       snow, qsurf, qsol, qbs1, agesno, &
20       tsoil, z0m, z0h, SFRWL, alb_dir, alb_dif, evap, fluxsens, fluxlat, fluxbs, &
[1067]21       tsurf_new, dflux_s, dflux_l, &
[3900]22       alt, slope, cloudf, &
[1865]23       snowhgt, qsnow, to_ice, sissnow, &
24       alb3, runoff, &
[1067]25       flux_u1, flux_v1)
[781]26
[1067]27    USE dimphy
[3974]28    USE geometry_mod,     ONLY : longitude,latitude
[3900]29    USE surface_data,     ONLY : type_ocean, calice, calsno, landice_opt, iflag_albcalc
30    USE fonte_neige_mod,  ONLY : fonte_neige,run_off_lic,fqcalving_global,ffonte_global,fqfonte_global,runofflic_global
[1067]31    USE cpl_mod,          ONLY : cpl_send_landice_fields
32    USE calcul_fluxs_mod
[4835]33    USE phys_local_var_mod, ONLY : zxrhoslic, zxustartlic, zxqsaltlic
[4414]34    USE phys_output_var_mod, ONLY : snow_o,zfra_o
[2728]35!FC
36    USE ioipsl_getin_p_mod, ONLY : getin_p
[4835]37    USE lmdz_blowing_snow_ini, ONLY : c_esalt_bs, zeta_bs, pbst_bs, prt_bs, rhoice_bs, rhohard_bs
38    USE lmdz_blowing_snow_ini, ONLY : rhofresh_bs, tau_eqsalt_bs, tau_dens0_bs, tau_densmin_bs
[3792]39#ifdef CPP_INLANDSIS
40    USE surf_inlandsis_mod,  ONLY : surf_inlandsis
41#endif
42
[1785]43    USE indice_sol_mod
[1067]44
[1785]45!    INCLUDE "indicesol.h"
[781]46    INCLUDE "dimsoil.h"
[793]47    INCLUDE "YOMCST.h"
48    INCLUDE "clesphys.h"
[781]49
50! Input variables
51!****************************************************************************************
52    INTEGER, INTENT(IN)                           :: itime, knon
53    INTEGER, DIMENSION(klon), INTENT(in)          :: knindex
54    REAL, INTENT(in)                              :: dtime
[888]55    REAL, DIMENSION(klon), INTENT(IN)             :: swnet ! net shortwave radiance
56    REAL, DIMENSION(klon), INTENT(IN)             :: lwnet ! net longwave radiance
[781]57    REAL, DIMENSION(klon), INTENT(IN)             :: tsurf
58    REAL, DIMENSION(klon), INTENT(IN)             :: p1lay
[1067]59    REAL, DIMENSION(klon), INTENT(IN)             :: cdragh, cdragm
[4523]60    REAL, DIMENSION(klon), INTENT(IN)             :: precip_rain, precip_snow, precip_bs
[781]61    REAL, DIMENSION(klon), INTENT(IN)             :: temp_air, spechum
[1067]62    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefH, AcoefQ
63    REAL, DIMENSION(klon), INTENT(IN)             :: BcoefH, BcoefQ
64    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefU, AcoefV, BcoefU, BcoefV
[4529]65    REAL, DIMENSION(klon), INTENT(IN)             :: AcoefQBS, BcoefQBS
[781]66    REAL, DIMENSION(klon), INTENT(IN)             :: ps
[4523]67    REAL, DIMENSION(klon), INTENT(IN)             :: u1, v1, gustiness, qbs1
[781]68    REAL, DIMENSION(klon), INTENT(IN)             :: rugoro
69    REAL, DIMENSION(klon,nbsrf), INTENT(IN)       :: pctsrf
70
[1865]71    LOGICAL,  INTENT(IN)                          :: debut   !true if first step
72    LOGICAL,  INTENT(IN)                          :: lafin   !true if last step
73    REAL, DIMENSION(klon), INTENT(IN)             :: rlon, rlat
74    REAL, DIMENSION(klon), INTENT(IN)             :: rmu0
75    REAL, DIMENSION(klon), INTENT(IN)             :: lwdownm !ylwdown
76    REAL, DIMENSION(klon), INTENT(IN)             :: albedo  !mean albedo
77    REAL, DIMENSION(klon), INTENT(IN)             :: pphi1   
[3900]78    REAL, DIMENSION(klon), INTENT(IN)             :: alt   !mean altitude of the grid box 
[1865]79    REAL, DIMENSION(klon), INTENT(IN)             :: slope   !mean slope in grid box 
80    REAL, DIMENSION(klon), INTENT(IN)             :: cloudf  !total cloud fraction
81
[781]82! In/Output variables
83!****************************************************************************************
84    REAL, DIMENSION(klon), INTENT(INOUT)          :: snow, qsol
85    REAL, DIMENSION(klon), INTENT(INOUT)          :: agesno
86    REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil
87
88! Output variables
89!****************************************************************************************
90    REAL, DIMENSION(klon), INTENT(OUT)            :: qsurf
[2243]91    REAL, DIMENSION(klon), INTENT(OUT)            :: z0m, z0h
[2227]92!albedo SB >>>
93!    REAL, DIMENSION(klon), INTENT(OUT)            :: alb1  ! new albedo in visible SW interval
94!    REAL, DIMENSION(klon), INTENT(OUT)            :: alb2  ! new albedo in near IR interval
[3792]95    REAL, DIMENSION(6), INTENT(IN)                :: SFRWL
96    REAL, DIMENSION(klon,nsw), INTENT(OUT)        :: alb_dir,alb_dif
[2227]97!albedo SB <<<
[781]98    REAL, DIMENSION(klon), INTENT(OUT)            :: evap, fluxsens, fluxlat
[4523]99    REAL, DIMENSION(klon), INTENT(OUT)            :: fluxbs
[888]100    REAL, DIMENSION(klon), INTENT(OUT)            :: tsurf_new
[781]101    REAL, DIMENSION(klon), INTENT(OUT)            :: dflux_s, dflux_l     
[1067]102    REAL, DIMENSION(klon), INTENT(OUT)            :: flux_u1, flux_v1
[781]103
[1865]104    REAL, DIMENSION(klon), INTENT(OUT)           :: alb3
105    REAL, DIMENSION(klon), INTENT(OUT)           :: qsnow   !column water in snow [kg/m2]
106    REAL, DIMENSION(klon), INTENT(OUT)           :: snowhgt !Snow height (m)
107    REAL, DIMENSION(klon), INTENT(OUT)           :: to_ice
108    REAL, DIMENSION(klon), INTENT(OUT)           :: sissnow
109    REAL, DIMENSION(klon), INTENT(OUT)           :: runoff  !Land ice runoff
110 
111
[781]112! Local variables
113!****************************************************************************************
114    REAL, DIMENSION(klon)    :: soilcap, soilflux
115    REAL, DIMENSION(klon)    :: cal, beta, dif_grnd
116    REAL, DIMENSION(klon)    :: zfra, alb_neig
[888]117    REAL, DIMENSION(klon)    :: radsol
[3792]118    REAL, DIMENSION(klon)    :: u0, v0, u1_lay, v1_lay, ustar
119    INTEGER                  :: i,j,nt
[3900]120    REAL, DIMENSION(klon)    :: fqfonte,ffonte
[4283]121    REAL, DIMENSION(klon)    :: run_off_lic_frac
[1865]122    REAL, DIMENSION(klon)    :: emis_new                  !Emissivity
123    REAL, DIMENSION(klon)    :: swdown,lwdown
[3900]124    REAL, DIMENSION(klon)    :: precip_snow_adv, snow_adv !Snow Drift precip./advection (not used in inlandsis)
125    REAL, DIMENSION(klon)    :: erod                      !erosion of surface snow (flux, kg/m2/s like evap)
126    REAL, DIMENSION(klon)    :: zsl_height, wind_velo     !surface layer height, wind spd
[1865]127    REAL, DIMENSION(klon)    :: dens_air,  snow_cont_air  !air density; snow content air
128    REAL, DIMENSION(klon)    :: alb_soil                  !albedo of underlying ice
129    REAL, DIMENSION(klon)    :: pexner                    !Exner potential
130    REAL                     :: pref
[3792]131    REAL, DIMENSION(klon,nsoilmx) :: tsoil0               !modif
132    REAL                          :: dtis                ! subtimestep
133    LOGICAL                       :: debut_is, lafin_is  ! debut and lafin for inlandsis
[1865]134
135    CHARACTER (len = 20)                      :: modname = 'surf_landice'
136    CHARACTER (len = 80)                      :: abort_message
137
[2728]138
[3900]139    REAL,DIMENSION(klon) :: alb1,alb2
[4523]140    REAL,DIMENSION(klon) :: precip_totsnow, evap_totsnow
[3900]141    REAL, DIMENSION (klon,6) :: alb6
[4835]142    REAL                   :: esalt
[4529]143    REAL                   :: lambdasalt,fluxsalt, csalt, nunu, aa, bb, cc
[4835]144    REAL                   :: tau_dens, maxerosion, fluxbs_2
145    REAL, DIMENSION(klon)  :: ws1, rhod, rhos, ustart0, ustart, qsalt, hsalt
[4672]146
[781]147! End definition
148!****************************************************************************************
[2728]149!FC
150!FC
151   REAL,SAVE :: alb_vis_sno_lic
152  !$OMP THREADPRIVATE(alb_vis_sno_lic)
153   REAL,SAVE :: alb_nir_sno_lic
154  !$OMP THREADPRIVATE(alb_nir_sno_lic)
155  LOGICAL, SAVE :: firstcall = .TRUE.
156  !$OMP THREADPRIVATE(firstcall)
157
158
[3792]159!FC firtscall initializations
160!******************************************************************************************
[2728]161  IF (firstcall) THEN
162  alb_vis_sno_lic=0.77
163  CALL getin_p('alb_vis_sno_lic',alb_vis_sno_lic)
164           PRINT*, 'alb_vis_sno_lic',alb_vis_sno_lic
165  alb_nir_sno_lic=0.77
166  CALL getin_p('alb_nir_sno_lic',alb_nir_sno_lic)
167           PRINT*, 'alb_nir_sno_lic',alb_nir_sno_lic
[3792]168 
[2728]169  firstcall=.false.
170  ENDIF
[3792]171!******************************************************************************************
172
[781]173! Initialize output variables
[1865]174    alb3(:) = 999999.
[888]175    alb2(:) = 999999.
176    alb1(:) = 999999.
[4523]177    fluxbs(:)=0. 
[1865]178    runoff(:) = 0.
[888]179!****************************************************************************************
180! Calculate total absorbed radiance at surface
181!
182!****************************************************************************************
183    radsol(:) = 0.0
184    radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
[781]185
186!****************************************************************************************
[4523]187
188!****************************************************************************************
[3792]189!  landice_opt = 0 : soil_model, calcul_flux, fonte_neige, ... 
[3901]190!  landice_opt = 1  : prepare and call INterace Lmdz SISvat (INLANDSIS)
[1865]191!****************************************************************************************
[3792]192
193
194    IF (landice_opt .EQ. 1) THEN
195
[3901]196!****************************************************************************************   
[3792]197! CALL to INLANDSIS interface
198!****************************************************************************************
199#ifdef CPP_INLANDSIS
200
201        debut_is=debut
202        lafin_is=.false.
203        ! Suppose zero surface speed
204        u0(:)            = 0.0
205        v0(:)            = 0.0
206
207
208        CALL calcul_flux_wind(knon, dtime, &
209         u0, v0, u1, v1, gustiness, cdragm, &
210         AcoefU, AcoefV, BcoefU, BcoefV, &
211         p1lay, temp_air, &
212         flux_u1, flux_v1)
213
214       
215       ! Set constants and compute some input for SISVAT
216       ! = 1000 hPa
217       ! and calculate incoming flux for SW and LW interval: swdown, lwdown
218       swdown(:)        = 0.0
219       lwdown(:)        = 0.0
[3900]220       snow_cont_air(:) = 0.  ! the snow content in air is not a prognostic variable of the model     
[3792]221       alb_soil(:)      = 0.4 ! before albedo(:) but here it is the ice albedo that we have to set
222       ustar(:)         = 0.
223       pref             = 100000.       
224       DO i = 1, knon
225          swdown(i)        = swnet(i)/(1-albedo(i))
226          lwdown(i)        = lwdownm(i)
227          wind_velo(i)     = u1(i)**2 + v1(i)**2
228          wind_velo(i)     = wind_velo(i)**0.5
229          pexner(i)        = (p1lay(i)/pref)**(RD/RCPD)
230          dens_air(i)      = p1lay(i)/RD/temp_air(i)  ! dry air density
231          zsl_height(i)    = pphi1(i)/RG     
232          tsoil0(i,:)      = tsoil(i,:) 
233          ustar(i)= (cdragm(i)*(wind_velo(i)**2))**0.5   
234       END DO
235       
236
237
[3900]238        dtis=dtime
[3792]239
[3900]240          IF (lafin) THEN
[3792]241            lafin_is=.true.
242          END IF
243
[3900]244          CALL surf_inlandsis(knon, rlon, rlat, knindex, itime, dtis, debut_is, lafin_is,&
245            rmu0, swdown, lwdown, albedo, pexner, ps, p1lay, precip_rain, precip_snow,   &
246            zsl_height, wind_velo, ustar, temp_air, dens_air, spechum, tsurf,&
247            rugoro, snow_cont_air, alb_soil, alt, slope, cloudf, &
248            radsol, qsol, tsoil0, snow, zfra, snowhgt, qsnow, to_ice, sissnow,agesno,   &
[3792]249            AcoefH, AcoefQ, BcoefH, BcoefQ, cdragm, cdragh, &
[3900]250            run_off_lic, fqfonte, ffonte, evap, erod, fluxsens, fluxlat,dflux_s, dflux_l, &
251            tsurf_new, alb1, alb2, alb3, alb6, &
252            emis_new, z0m, z0h, qsurf)
[3792]253
[3900]254          debut_is=.false.
[3792]255
256
[3900]257        ! Treatment of snow melting and calving
[3792]258
[3900]259        ! for consistency with standard LMDZ, add calving to run_off_lic
260        run_off_lic(:)=run_off_lic(:) + to_ice(:)
261
262        DO i = 1, knon
263           ffonte_global(knindex(i),is_lic)    = ffonte(i)
264           fqfonte_global(knindex(i),is_lic)   = fqfonte(i)! net melting= melting - refreezing
265           fqcalving_global(knindex(i),is_lic) = to_ice(i) ! flux
266           runofflic_global(knindex(i)) = run_off_lic(i)
267        ENDDO
268        ! Here, we assume that the calving term is equal to the to_ice term
269        ! (no ice accumulation)
270
271
[3792]272#else
[3901]273       abort_message='Pb de coherence: landice_opt = 1  mais CPP_INLANDSIS = .false.'
[3792]274       CALL abort_physic(modname,abort_message,1)
275#endif
276
277
278    ELSE
279
280!****************************************************************************************
[781]281! Soil calculations
282!
283!****************************************************************************************
[3780]284
285    ! EV: use calbeta
286    CALL calbeta(dtime, is_lic, knon, snow, qsol, beta, cal, dif_grnd)
287
288
289    ! use soil model and recalculate properly cal
[781]290    IF (soil_model) THEN
[3974]291       CALL soil(dtime, is_lic, knon, snow, tsurf, qsol, &
292        & longitude(knindex(1:knon)), latitude(knindex(1:knon)), tsoil, soilcap, soilflux)
[781]293       cal(1:knon) = RCPD / soilcap(1:knon)
294       radsol(1:knon)  = radsol(1:knon) + soilflux(1:knon)
295    ELSE
296       cal = RCPD * calice
297       WHERE (snow > 0.0) cal = RCPD * calsno
298    ENDIF
299
300
301!****************************************************************************************
302! Calulate fluxes
303!
304!****************************************************************************************
[3792]305!    beta(:) = 1.0
306!    dif_grnd(:) = 0.0
[781]307
[1067]308! Suppose zero surface speed
309    u0(:)=0.0
310    v0(:)=0.0
311    u1_lay(:) = u1(:) - u0(:)
312    v1_lay(:) = v1(:) - v0(:)
313
[781]314    CALL calcul_fluxs(knon, is_lic, dtime, &
[2254]315         tsurf, p1lay, cal, beta, cdragh, cdragh, ps, &
[781]316         precip_rain, precip_snow, snow, qsurf,  &
[2240]317         radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
[2254]318         1.,AcoefH, AcoefQ, BcoefH, BcoefQ, &
[781]319         tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
320
[1067]321    CALL calcul_flux_wind(knon, dtime, &
[2240]322         u0, v0, u1, v1, gustiness, cdragm, &
[1067]323         AcoefU, AcoefV, BcoefU, BcoefV, &
324         p1lay, temp_air, &
325         flux_u1, flux_v1)
[781]326
327
328!****************************************************************************************
329! Calculate albedo
330!
331!****************************************************************************************
[3780]332
[781]333!
334!IM: plusieurs choix/tests sur l'albedo des "glaciers continentaux"
[888]335!       alb1(1 : knon)  = 0.6 !IM cf FH/GK
336!       alb1(1 : knon)  = 0.82
337!       alb1(1 : knon)  = 0.77 !211003 Ksta0.77
338!       alb1(1 : knon)  = 0.8 !KstaTER0.8 & LMD_ARMIP5
[781]339!IM: KstaTER0.77 & LMD_ARMIP6   
340
[3780]341! Attantion: alb1 and alb2 are not the same!
[2728]342    alb1(1:knon)  = alb_vis_sno_lic
343    alb2(1:knon)  = alb_nir_sno_lic
[781]344
345
346!****************************************************************************************
347! Rugosity
348!
349!****************************************************************************************
[4245]350    z0m = z0m_landice
351    z0h = z0h_landice
352    !z0m = SQRT(z0m**2+rugoro**2)
[2243]353
[781]354
[4835]355!****************************************************************************************
356! Simple blowing snow param
357!****************************************************************************************
358! we proceed in 2 steps:
359! first we erode - if possible -the accumulated snow during the time step
360! then we update the density of the underlying layer and see if we can also erode
361! this layer
[1865]362
[4529]363
[4835]364   if (ok_bs) then
365       fluxbs(:)=0.
366       do j=1,klon
367          ws1(j)=(u1(j)**2+v1(j)**2)**0.5
368          ustar(j)=(cdragm(j)*(u1(j)**2+v1(j)**2))**0.5
369          rhod(j)=p1lay(j)/RD/temp_air(j)
370          ustart0(j) =(log(2.868)-log(1.625))/0.085*sqrt(cdragm(j))
371       enddo
372
373       ! 1st step: erosion of fresh snow accumulated during the time step
374       do j=1, knon
375
376           rhos(j)=rhofresh_bs
377           ! blowing snow flux formula used in MAR
378           ustart(j)=ustart0(j)*exp(max(rhoice_bs/rhofresh_bs-rhoice_bs/rhos(j),0.))*exp(max(0.,rhos(j)-rhohard_bs))
379           ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs
380           ! computation of qbs at the top of the saltation layer
381           ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001)
382           esalt=1./(c_esalt_bs*max(1.e-6,ustar(j)))
383           hsalt(j)=0.08436*(max(1.e-6,ustar(j))**1.27)
384           qsalt(j)=(max(ustar(j)**2-ustart(j)**2,0.))/(RG*hsalt(j))*esalt
385           ! calculation of erosion (flux positive towards the surface here)
386           ! consistent with implicit resolution of turbulent mixing equation
387           ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s
388           ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level)
389           ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer
390           ! (rho*qsalt*hsalt)
391           ! during this first step we also lower bound the erosion to the amount of fresh snow accumulated during the time step
392           maxerosion=min(precip_snow(j),hsalt(j)*qsalt(j)*rhod(j)/tau_eqsalt_bs)
393
394           fluxbs(j)=rhod(j)*ws1(j)*cdragm(j)*(AcoefQBS(j)-qsalt(j)) &
395                   / (1.-rhod(j)*ws1(j)*cdragm(j)*BcoefQBS(j)*dtime)
396           !fluxbs(j)= rhod*ws1(j)*cdragm(j)*(qbs1(j)-qsalt(j))
397           fluxbs(j)=max(-maxerosion,fluxbs(j))
398       enddo
399
400
401       ! we now compute the snow age of the overlying layer (snow surface after erosion of the fresh snow accumulated during the time step)
402       ! this is done through the routine albsno
403       CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)+fluxbs(:))
404
405       ! 2nd step:
[4529]406       ! computation of threshold friction velocity
407       ! which depends on surface snow density
[4672]408       do j = 1, knon
[4523]409           ! estimation of snow density
410           ! snow density increases with snow age and
[4672]411           ! increases even faster in case of sedimentation of blowing snow or rain
[4835]412           tau_dens=max(tau_densmin_bs, tau_dens0_bs*exp(-abs(precip_bs(j))/pbst_bs - &
413                    abs(precip_rain(j))/prt_bs)*exp(-max(tsurf(j)-RTT,0.)))
414           rhos(j)=rhofresh_bs+(rhohard_bs-rhofresh_bs)*(1.-exp(-agesno(j)*86400.0/tau_dens))
[4523]415           ! blowing snow flux formula used in MAR
[4835]416           ustart(j)=ustart0(j)*exp(max(rhoice_bs/rhofresh_bs-rhoice_bs/rhos(j),0.))*exp(max(0.,rhos(j)-rhohard_bs))
417           ! we have multiplied by exp to prevent erosion when rhos>rhohard_bs
418           ! computation of qbs at the top of the saltation layer
419           ! default formulation from MAR model (Amory et al. 2021, Gallee et al. 2001)
420           esalt=c_esalt_bs*ustar(j)
421           hsalt(j)=0.08436*(max(1.e-6,ustar(j))**1.27)
422           qsalt(j)=(max(ustar(j)**2-ustart(j)**2,0.))/(RG*hsalt(j))*esalt
423           ! calculation of erosion (flux positive towards the surface here)
424           ! consistent with implicit resolution of turbulent mixing equation
425           ! Nemoto and Nishimura 2004 show that steady-state saltation is achieved within a time tau_eqsalt_bs of about 10s
426           ! we thus prevent snowerosion (snow particle transfer from the saltation layer to the first model level)
427           ! integrated over tau_eqsalt_bs to exceed the total mass of snow particle in the saltation layer
428           ! (rho*qsalt*hsalt)
429           maxerosion=hsalt(j)*qsalt(j)*rhod(j)/tau_eqsalt_bs
430           fluxbs_2=rhod(j)*ws1(j)*cdragm(j)*(AcoefQBS(j)-qsalt(j)) &
431                   / (1.-rhod(j)*ws1(j)*cdragm(j)*BcoefQBS(j)*dtime)
432           !fluxbs_2= rhod*ws1(j)*cdragm(j)*(qbs1(j)-qsalt(j))
433           fluxbs_2=max(-maxerosion,fluxbs_2)
434           fluxbs(j)=fluxbs(j)+fluxbs_2
[4523]435       enddo
436
[4835]437       ! for outputs       
438        do j=1, knon
[4672]439              i = knindex(j)
440              zxustartlic(i) = ustart(j)
441              zxrhoslic(i) = rhos(j)
[4835]442              zxqsaltlic(i)=qsalt(j)
443        enddo
[4523]444
445
[4835]446  else
447  ! those lines are useful to calculate the snow age
448       CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:))
[4523]449
[4835]450  endif ! if ok_bs
[4523]451
[4835]452
453
[4523]454!****************************************************************************************
[4672]455! Calculate snow amount
[4523]456!   
457!****************************************************************************************
458    IF (ok_bs) THEN
[4526]459      precip_totsnow(:)=precip_snow(:)+precip_bs(:)
460      evap_totsnow(:)=evap(:)-fluxbs(:) ! flux bs is positive towards the surface (snow erosion)
[4523]461    ELSE
[4526]462      precip_totsnow(:)=precip_snow(:)
463      evap_totsnow(:)=evap(:)
[4523]464    ENDIF
[4672]465   
466 
[4523]467    CALL fonte_neige(knon, is_lic, knindex, dtime, &
468         tsurf, precip_rain, precip_totsnow,  &
469         snow, qsol, tsurf_new, evap_totsnow)
[4672]470   
471   
[4523]472    WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.                                         
473    zfra(1:knon) = MAX(0.0,MIN(1.0,snow(1:knon)/(snow(1:knon)+10.0))) 
474
475
[3792]476    END IF ! landice_opt
477
478
[781]479!****************************************************************************************
480! Send run-off on land-ice to coupler if coupled ocean.
[3903]481! run_off_lic has been calculated in fonte_neige or surf_inlandsis
[4283]482! If landice_opt>=2, corresponding call is done from surf_land_orchidee
[781]483!****************************************************************************************
[4283]484    IF (type_ocean=='couple' .AND. landice_opt .LT. 2) THEN
485       ! Compress fraction where run_off_lic is active (here all pctsrf(is_lic))
486       run_off_lic_frac(:)=0.0
487       DO j = 1, knon
488          i = knindex(j)
489          run_off_lic_frac(j) = pctsrf(i,is_lic)
490       ENDDO
491
492       CALL cpl_send_landice_fields(itime, knon, knindex, run_off_lic, run_off_lic_frac)
[781]493    ENDIF
[1865]494
495 ! transfer runoff rate [kg/m2/s](!) to physiq for output
496    runoff(1:knon)=run_off_lic(1:knon)/dtime
497
[1403]498       snow_o=0.
499       zfra_o = 0.
500       DO j = 1, knon
501           i = knindex(j)
502           snow_o(i) = snow(j)
503           zfra_o(i) = zfra(j)
504       ENDDO
505
506
[2227]507!albedo SB >>>
508     select case(NSW)
509     case(2)
510       alb_dir(1:knon,1)=alb1(1:knon)
511       alb_dir(1:knon,2)=alb2(1:knon)
512     case(4)
513       alb_dir(1:knon,1)=alb1(1:knon)
514       alb_dir(1:knon,2)=alb2(1:knon)
515       alb_dir(1:knon,3)=alb2(1:knon)
516       alb_dir(1:knon,4)=alb2(1:knon)
517     case(6)
518       alb_dir(1:knon,1)=alb1(1:knon)
519       alb_dir(1:knon,2)=alb1(1:knon)
520       alb_dir(1:knon,3)=alb1(1:knon)
521       alb_dir(1:knon,4)=alb2(1:knon)
522       alb_dir(1:knon,5)=alb2(1:knon)
523       alb_dir(1:knon,6)=alb2(1:knon)
[3900]524
[3901]525       IF ((landice_opt .EQ. 1) .AND. (iflag_albcalc .EQ. 2)) THEN
[3900]526       alb_dir(1:knon,1)=alb6(1:knon,1)
527       alb_dir(1:knon,2)=alb6(1:knon,2)
528       alb_dir(1:knon,3)=alb6(1:knon,3)
529       alb_dir(1:knon,4)=alb6(1:knon,4)
530       alb_dir(1:knon,5)=alb6(1:knon,5)
531       alb_dir(1:knon,6)=alb6(1:knon,6)
532       ENDIF
533
[2227]534     end select
535alb_dif=alb_dir
536!albedo SB <<<
537
538
[781]539  END SUBROUTINE surf_landice
540!
541!****************************************************************************************
542!
543END MODULE surf_landice_mod
544
545
546
Note: See TracBrowser for help on using the repository browser.