LMDZ
surf_landice_mod.F90
Go to the documentation of this file.
1 !
3 
4  IMPLICIT NONE
5 
6 CONTAINS
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, temp_air, spechum, &
15  acoefh, acoefq, bcoefh, bcoefq, &
16  acoefu, acoefv, bcoefu, bcoefv, &
17  ps, u1, v1, gustiness, rugoro, pctsrf, &
18  snow, qsurf, qsol, agesno, &
19  tsoil, z0m, z0h, sfrwl, alb_dir, alb_dif, evap, fluxsens, fluxlat, &
20  tsurf_new, dflux_s, dflux_l, &
21  slope, cloudf, &
22  snowhgt, qsnow, to_ice, sissnow, &
23  alb3, runoff, &
24  flux_u1, flux_v1)
25 
26  USE dimphy
32 #ifdef CPP_SISVAT
33  USE surf_sisvat_mod, ONLY : surf_sisvat
34 #endif
35  USE indice_sol_mod
36 
37 ! INCLUDE "indicesol.h"
38  include "dimsoil.h"
39  include "YOMCST.h"
40  include "clesphys.h"
41 
42 ! Input variables
43 !****************************************************************************************
44  INTEGER, INTENT(IN) :: itime, knon
45  INTEGER, DIMENSION(klon), INTENT(in) :: knindex
46  REAL, INTENT(in) :: dtime
47  REAL, DIMENSION(klon), INTENT(IN) :: swnet ! net shortwave radiance
48  REAL, DIMENSION(klon), INTENT(IN) :: lwnet ! net longwave radiance
49  REAL, DIMENSION(klon), INTENT(IN) :: tsurf
50  REAL, DIMENSION(klon), INTENT(IN) :: p1lay
51  REAL, DIMENSION(klon), INTENT(IN) :: cdragh, cdragm
52  REAL, DIMENSION(klon), INTENT(IN) :: precip_rain, precip_snow
53  REAL, DIMENSION(klon), INTENT(IN) :: temp_air, spechum
54  REAL, DIMENSION(klon), INTENT(IN) :: AcoefH, AcoefQ
55  REAL, DIMENSION(klon), INTENT(IN) :: BcoefH, BcoefQ
56  REAL, DIMENSION(klon), INTENT(IN) :: AcoefU, AcoefV, BcoefU, BcoefV
57  REAL, DIMENSION(klon), INTENT(IN) :: ps
58  REAL, DIMENSION(klon), INTENT(IN) :: u1, v1, gustiness
59  REAL, DIMENSION(klon), INTENT(IN) :: rugoro
60  REAL, DIMENSION(klon,nbsrf), INTENT(IN) :: pctsrf
61 
62  LOGICAL, INTENT(IN) :: debut !true if first step
63  LOGICAL, INTENT(IN) :: lafin !true if last step
64  REAL, DIMENSION(klon), INTENT(IN) :: rlon, rlat
65  REAL, DIMENSION(klon), INTENT(IN) :: rmu0
66  REAL, DIMENSION(klon), INTENT(IN) :: lwdownm !ylwdown
67  REAL, DIMENSION(klon), INTENT(IN) :: albedo !mean albedo
68  REAL, DIMENSION(klon), INTENT(IN) :: pphi1
69  REAL, DIMENSION(klon), INTENT(IN) :: slope !mean slope in grid box
70  REAL, DIMENSION(klon), INTENT(IN) :: cloudf !total cloud fraction
71 
72 ! In/Output variables
73 !****************************************************************************************
74  REAL, DIMENSION(klon), INTENT(INOUT) :: snow, qsol
75  REAL, DIMENSION(klon), INTENT(INOUT) :: agesno
76  REAL, DIMENSION(klon, nsoilmx), INTENT(INOUT) :: tsoil
77 
78 ! Output variables
79 !****************************************************************************************
80  REAL, DIMENSION(klon), INTENT(OUT) :: qsurf
81  REAL, DIMENSION(klon), INTENT(OUT) :: z0m, z0h
82 !albedo SB >>>
83 ! REAL, DIMENSION(klon), INTENT(OUT) :: alb1 ! new albedo in visible SW interval
84 ! REAL, DIMENSION(klon), INTENT(OUT) :: alb2 ! new albedo in near IR interval
85  REAL, DIMENSION(6), INTENT(IN) ::SFRWL
86  REAL, DIMENSION(klon,nsw), INTENT(OUT) ::alb_dir,alb_dif
87 !albedo SB <<<
88  REAL, DIMENSION(klon), INTENT(OUT) :: evap, fluxsens, fluxlat
89  REAL, DIMENSION(klon), INTENT(OUT) :: tsurf_new
90  REAL, DIMENSION(klon), INTENT(OUT) :: dflux_s, dflux_l
91  REAL, DIMENSION(klon), INTENT(OUT) :: flux_u1, flux_v1
92 
93  REAL, DIMENSION(klon), INTENT(OUT) :: alb3
94  REAL, DIMENSION(klon), INTENT(OUT) :: qsnow !column water in snow [kg/m2]
95  REAL, DIMENSION(klon), INTENT(OUT) :: snowhgt !Snow height (m)
96  REAL, DIMENSION(klon), INTENT(OUT) :: to_ice
97  REAL, DIMENSION(klon), INTENT(OUT) :: sissnow
98  REAL, DIMENSION(klon), INTENT(OUT) :: runoff !Land ice runoff
99 
100 
101 ! Local variables
102 !****************************************************************************************
103  REAL, DIMENSION(klon) :: soilcap, soilflux
104  REAL, DIMENSION(klon) :: cal, beta, dif_grnd
105  REAL, DIMENSION(klon) :: zfra, alb_neig
106  REAL, DIMENSION(klon) :: radsol
107  REAL, DIMENSION(klon) :: u0, v0, u1_lay, v1_lay
108  INTEGER :: i,j
109 
110  REAL, DIMENSION(klon) :: emis_new !Emissivity
111  REAL, DIMENSION(klon) :: swdown,lwdown
112  REAL, DIMENSION(klon) :: precip_snow_adv, snow_adv !Snow Drift precip./advection
113  REAL, DIMENSION(klon) :: bl_height, wind_velo !height boundary layer, wind spd
114  REAL, DIMENSION(klon) :: dens_air, snow_cont_air !air density; snow content air
115  REAL, DIMENSION(klon) :: alb_soil !albedo of underlying ice
116  REAL, DIMENSION(klon) :: pexner !Exner potential
117  REAL :: pref
118  REAL, DIMENSION(klon,nsoilmx) :: tsoil0 !modfi
119 
120  CHARACTER (len = 20) :: modname = 'surf_landice'
121  CHARACTER (len = 80) :: abort_message
122 
123 !albedo SB >>>
124  real,dimension(klon) :: alb1,alb2
125 !albedo SB <<<
126 
127 ! End definition
128 !****************************************************************************************
129 !
130 ! Initialize output variables
131  alb3(:) = 999999.
132  alb2(:) = 999999.
133  alb1(:) = 999999.
134 
135  runoff(:) = 0.
136 !****************************************************************************************
137 ! Calculate total absorbed radiance at surface
138 !
139 !****************************************************************************************
140  radsol(:) = 0.0
141  radsol(1:knon) = swnet(1:knon) + lwnet(1:knon)
142 
143 !****************************************************************************************
144 ! ok_snow = TRUE : prepare and call SISVAT snow model
145 ! ok_snow = FALSE : soil_model, calcul_flux, fonte_neige, ...
146 !
147 !****************************************************************************************
148  IF (ok_snow) THEN
149 #ifdef CPP_SISVAT
150  ! Prepare for calling SISVAT
151 
152  ! Calculate incoming flux for SW and LW interval: swdown, lwdown
153  swdown(:) = 0.0
154  lwdown(:) = 0.0
155  DO i = 1, knon
156  swdown(i) = swnet(i)/(1-albedo(i))
157  lwdown(i) = lwdownm(i)
158  END DO
159 
160  ! Set constants and compute some input for SISVAT
161  snow_adv(:) = 0. ! no snow blown in for now
162  snow_cont_air(:) = 0.
163  alb_soil(:) = albedo(:)
164  pref = 100000. ! = 1000 hPa
165  DO i = 1, knon
166  wind_velo(i) = u1(i)**2 + v1(i)**2
167  wind_velo(i) = wind_velo(i)**0.5
168  pexner(i) = (p1lay(i)/pref)**(rd/rcpd)
169  dens_air(i) = p1lay(i)/rd/temp_air(i) ! dry air density
170  bl_height(i) = pphi1(i)/rg
171  END DO
172 
173 !****************************************************************************************
174 ! CALL to SISVAT interface
175 !
176 !****************************************************************************************
177  ! config: compute everything with SV but temperatures afterwards with soil/calculfluxs
178  DO i = 1, knon
179  tsoil0(i,:)=tsoil(i,:)
180  END DO
181  ! Martin
182  print*, 'on appelle surf_sisvat'
183  ! Martin
184  CALL surf_sisvat(knon, rlon, rlat, knindex, itime, dtime, debut, lafin, &
185  rmu0, swdown, lwdown, pexner, ps, p1lay, &
186  precip_rain, precip_snow, precip_snow_adv, snow_adv, &
187  bl_height, wind_velo, temp_air, dens_air, spechum, tsurf, &
188  rugoro, snow_cont_air, alb_soil, slope, cloudf, &
189  radsol, qsol, tsoil0, snow, snowhgt, qsnow, to_ice,sissnow, agesno, &
190  acoefh, acoefq, bcoefh, bcoefq, cdragh, &
191  run_off_lic, evap, fluxsens, fluxlat, dflux_s, dflux_l, &
192  tsurf_new, alb1, alb2, alb3, &
193  emis_new, z0m, qsurf)
194  z0h(1:knon)=z0m(1:knon) ! en attendant mieux
195 
196  ! Suppose zero surface speed
197  u0(:) = 0.0
198  v0(:) = 0.0
199  ! The calculation of heat/water fluxes, otherwise done by "CALL calcul_fluxs" is
200  ! integrated in SISVAT, using the same method. It can be found in "sisvat.f", in the
201  ! subroutine "SISVAT_TS2".
202  ! u0, v0=0., dif_grnd=0. and beta=1 are assumed there!
203 
204  CALL calcul_flux_wind(knon, dtime, &
205  u0, v0, u1, v1, gustiness, cdragm, &
206  acoefu, acoefv, bcoefu, bcoefv, &
207  p1lay, temp_air, &
208  flux_u1, flux_v1)
209 #else
210  abort_message='Pb de coherence: ok_snow = .true. mais CPP_SISVAT = .false.'
211  CALL abort_physic(modname,abort_message,1)
212 #endif
213  ELSE ! ok_snow=FALSE
214 
215 !****************************************************************************************
216 ! Soil calculations
217 !
218 !****************************************************************************************
219  IF (soil_model) THEN
220  CALL soil(dtime, is_lic, knon, snow, tsurf, tsoil, soilcap, soilflux)
221  cal(1:knon) = rcpd / soilcap(1:knon)
222  radsol(1:knon) = radsol(1:knon) + soilflux(1:knon)
223  ELSE
224  cal = rcpd * calice
225  WHERE (snow > 0.0) cal = rcpd * calsno
226  ENDIF
227 
228 
229 !****************************************************************************************
230 ! Calulate fluxes
231 !
232 !****************************************************************************************
233  beta(:) = 1.0
234  dif_grnd(:) = 0.0
235 
236 ! Suppose zero surface speed
237  u0(:)=0.0
238  v0(:)=0.0
239  u1_lay(:) = u1(:) - u0(:)
240  v1_lay(:) = v1(:) - v0(:)
241 
242  CALL calcul_fluxs(knon, is_lic, dtime, &
243  tsurf, p1lay, cal, beta, cdragh, cdragh, ps, &
244  precip_rain, precip_snow, snow, qsurf, &
245  radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, gustiness, &
246  1.,acoefh, acoefq, bcoefh, bcoefq, &
247  tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
248 
249  CALL calcul_flux_wind(knon, dtime, &
250  u0, v0, u1, v1, gustiness, cdragm, &
251  acoefu, acoefv, bcoefu, bcoefv, &
252  p1lay, temp_air, &
253  flux_u1, flux_v1)
254 
255 !****************************************************************************************
256 ! Calculate snow height, age, run-off,..
257 !
258 !****************************************************************************************
259  CALL fonte_neige( knon, is_lic, knindex, dtime, &
260  tsurf, precip_rain, precip_snow, &
261  snow, qsol, tsurf_new, evap)
262 
263 
264 !****************************************************************************************
265 ! Calculate albedo
266 !
267 !****************************************************************************************
268  CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:))
269  WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.
270  zfra(1:knon) = max(0.0,min(1.0,snow(1:knon)/(snow(1:knon)+10.0)))
271  alb1(1:knon) = alb_neig(1:knon)*zfra(1:knon) + &
272  0.6 * (1.0-zfra(1:knon))
273 !
274 !IM: plusieurs choix/tests sur l'albedo des "glaciers continentaux"
275 ! alb1(1 : knon) = 0.6 !IM cf FH/GK
276 ! alb1(1 : knon) = 0.82
277 ! alb1(1 : knon) = 0.77 !211003 Ksta0.77
278 ! alb1(1 : knon) = 0.8 !KstaTER0.8 & LMD_ARMIP5
279 !IM: KstaTER0.77 & LMD_ARMIP6
280 
281 ! Attantion: alb1 and alb2 are the same!
282  alb1(1:knon) = 0.77
283  alb2(1:knon) = alb1(1:knon)
284 
285 
286 !****************************************************************************************
287 ! Rugosity
288 !
289 !****************************************************************************************
290  z0m=1.e-3
291  z0h = z0m
292  z0m = sqrt(z0m**2+rugoro**2)
293 
294  END IF ! ok_snow
295 
296 
297 !****************************************************************************************
298 ! Send run-off on land-ice to coupler if coupled ocean.
299 ! run_off_lic has been calculated in fonte_neige or surf_sisvat
300 !
301 !****************************************************************************************
302  IF (type_ocean=='couple') THEN
303  CALL cpl_send_landice_fields(itime, knon, knindex, run_off_lic)
304  ENDIF
305 
306  ! transfer runoff rate [kg/m2/s](!) to physiq for output
307  runoff(1:knon)=run_off_lic(1:knon)/dtime
308 
309 
310 !****************************************************************************************
311  snow_o=0.
312  zfra_o = 0.
313  DO j = 1, knon
314  i = knindex(j)
315  snow_o(i) = snow(j)
316  zfra_o(i) = zfra(j)
317  ENDDO
318 
319 !****************************************************************************************
320  snow_o=0.
321  zfra_o = 0.
322  DO j = 1, knon
323  i = knindex(j)
324  snow_o(i) = snow(j)
325  zfra_o(i) = zfra(j)
326  ENDDO
327 
328 
329 !albedo SB >>>
330  select case(nsw)
331  case(2)
332  alb_dir(1:knon,1)=alb1(1:knon)
333  alb_dir(1:knon,2)=alb2(1:knon)
334  case(4)
335  alb_dir(1:knon,1)=alb1(1:knon)
336  alb_dir(1:knon,2)=alb2(1:knon)
337  alb_dir(1:knon,3)=alb2(1:knon)
338  alb_dir(1:knon,4)=alb2(1:knon)
339  case(6)
340  alb_dir(1:knon,1)=alb1(1:knon)
341  alb_dir(1:knon,2)=alb1(1:knon)
342  alb_dir(1:knon,3)=alb1(1:knon)
343  alb_dir(1:knon,4)=alb2(1:knon)
344  alb_dir(1:knon,5)=alb2(1:knon)
345  alb_dir(1:knon,6)=alb2(1:knon)
346  end select
347 alb_dif=alb_dir
348 !albedo SB <<<
349 
350 
351 
352 
353  END SUBROUTINE surf_landice
354 !
355 !****************************************************************************************
356 !
357 END MODULE surf_landice_mod
358 
359 
360 
subroutine albsno(klon, knon, dtime, agesno, alb_neig_grid, precip_snow)
Definition: albsno.F90:5
!$Header!c include clesph0 h c COMMON clesph0 soil_model
Definition: clesph0.h:6
subroutine fonte_neige(knon, nisurf, knindex, dtime, tsurf, precip_rain, precip_snow, snow, qsol, tsurf_new, evap)
real, parameter calsno
Definition: surface_data.F90:8
subroutine calcul_flux_wind(knon, dtime, u0, v0, u1, v1, gustiness, cdrag_m, AcoefU, AcoefV, BcoefU, BcoefV, p1lay, t1lay, flux_u1, flux_v1)
real, dimension(:), allocatable run_off_lic
integer, save klon
Definition: dimphy.F90:3
real, dimension(:), allocatable, save zfra_o
integer, parameter is_lic
subroutine, public cpl_send_landice_fields(itime, knon, knindex, rlic_in)
Definition: cpl_mod.F90:967
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, temp_air, spechum, AcoefH, AcoefQ, BcoefH, BcoefQ, AcoefU, AcoefV, BcoefU, BcoefV, ps, u1, v1, gustiness, rugoro, pctsrf, snow, qsurf, qsol, agesno, tsoil, z0m, z0h, SFRWL, alb_dir, alb_dif, evap, fluxsens, fluxlat, tsurf_new, dflux_s, dflux_l, slope, cloudf, snowhgt, qsnow, to_ice, sissnow, alb3, runoff, flux_u1, flux_v1)
c c $Id c nbregdyn DO klon c rlat(i) c ENDIF!lon c ENDIF!lat ENDIF!pctsrf ENDDO!klon ENDDO!nbregdyn cIM 190504 ENDIF!ok_regdyn cIM somme de toutes les nhistoW BEG IF(debut) THEN DO nreg
character(len=6), save type_ocean
Definition: albedo.F90:2
logical, save ok_snow
subroutine abort_physic(modname, message, ierr)
Definition: abort_physic.F90:3
real, parameter calice
Definition: surface_data.F90:6
subroutine calcul_fluxs(knon, nisurf, dtime, tsurf, p1lay, cal, beta, cdragh, cdragq, ps, precip_rain, precip_snow, snow, qsurf, radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, gustiness, fqsat, petAcoef, peqAcoef, petBcoef, peqBcoef, tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
subroutine surf_sisvat(knon, rlon, rlat, ikl2i, itime, dtime, debut, lafin, rmu0, swdown, lwdown, pexner, ps, p1lay, precip_rain, precip_snow, precip_snow_adv, snow_adv, bl_height, wind_velo, temp_air, dens_air, spechum, tsurf, rugos, snow_cont_air, alb_soil, slope, cloudf, radsol, qsol, tsoil, snow, snowhgt, qsnow, to_ice, sissnow, agesno, AcoefH, AcoefQ, BcoefH, BcoefQ, cdragh, runoff_lic, evap, fluxsens, fluxlat, dflux_s, dflux_l, tsurf_new, alb1, alb2, alb3, emis_new, z0_new, qsurf)
Definition: dimphy.F90:1
c c $Id c nbregdyn DO klon c rlon(i)
real, dimension(:), allocatable, save snow_o
real rg
Definition: comcstphy.h:1
subroutine soil(ptimestep, indice, knon, snow, ptsrf, ptsoil, pcapcal, pfluxgrd)
Definition: soil.F90:6