source: LMDZ4/branches/IPSL-CM4_IPCC_branch/libf/phylmd/interface_surf.F90 @ 634

Last change on this file since 634 was 628, checked in by Laurent Fairhead, 20 years ago

Modif pour compatibilité avec OASIS3 AC
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 95.0 KB
Line 
1!
2! $Header$
3!
4
5  MODULE interface_surf
6
7! Ce module regroupe toutes les routines gerant l'interface entre le modele
8! atmospherique et les modeles de surface (sols continentaux, oceans, glaces)
9! Les routines sont les suivantes:
10!
11!   interfsurf_*: routines d'aiguillage vers les interfaces avec les
12!                 differents modeles de surface
13!   interfsol\
14!             > routines d'interface proprement dite
15!   interfoce/
16!
17!   interfstart: routine d'initialisation et de lecture de l'etat initial
18!                "interface"
19!   interffin  : routine d'ecriture de l'etat de redemmarage de l'interface
20!
21!
22! L. Fairhead, LMD, 02/2000
23
24  USE ioipsl
25
26  IMPLICIT none
27
28  PRIVATE
29  PUBLIC :: interfsurf,interfsurf_hq, gath2cpl
30
31  INTERFACE interfsurf
32    module procedure interfsurf_hq, interfsurf_vent
33  END INTERFACE
34
35  INTERFACE interfoce
36    module procedure interfoce_cpl, interfoce_slab, interfoce_lim
37  END INTERFACE
38
39#include "YOMCST.inc"
40#include "indicesol.inc"
41
42
43! run_off      ruissellement total
44  REAL, ALLOCATABLE, DIMENSION(:),SAVE    :: run_off, run_off_lic
45  real, allocatable, dimension(:),save    :: coastalflow, riverflow
46!!$PB
47  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE :: tmp_rriv, tmp_rcoa,tmp_rlic
48!! pour simuler la fonte des glaciers antarctiques
49  REAL, ALLOCATABLE, DIMENSION(:,:), SAVE :: coeff_iceberg
50  real, save                              :: surf_maille
51  real, save                              :: cte_flux_iceberg = 6.3e7
52  integer, save                           :: num_antarctic = 1
53  REAL, save                              :: tau_calv
54!!$
55  CONTAINS
56!
57!############################################################################
58!
59  SUBROUTINE interfsurf_hq(itime, dtime, date0, jour, rmu0, &
60      & klon, iim, jjm, nisurf, knon, knindex, pctsrf, &
61      & rlon, rlat, cufi, cvfi,&
62      & debut, lafin, ok_veget, soil_model, nsoilmx, tsoil, qsol,&
63      & zlev,  u1_lay, v1_lay, temp_air, spechum, epot_air, ccanopy, &
64      & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
65      & precip_rain, precip_snow, sollw, sollwdown, swnet, swdown, &
66      & fder, taux, tauy, rugos, rugoro, &
67      & albedo, snow, qsurf, &
68      & tsurf, p1lay, ps, radsol, &
69      & ocean, npas, nexca, zmasq, &
70      & evap, fluxsens, fluxlat, dflux_l, dflux_s, &             
71      & tsol_rad, tsurf_new, alb_new, alblw, emis_new, &
72      & z0_new, pctsrf_new, agesno,fqcalving,ffonte, run_off_lic_0)
73
74
75! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
76! (sols continentaux, oceans, glaces) pour les fluxs de chaleur et d'humidite.
77! En pratique l'interface se fait entre la couche limite du modele
78! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
79!
80!
81! L.Fairhead 02/2000
82!
83! input:
84!   itime        numero du pas de temps
85!   klon         nombre total de points de grille
86!   iim, jjm     nbres de pts de grille
87!   dtime        pas de temps de la physique (en s)
88!   date0        jour initial
89!   jour         jour dans l'annee en cours,
90!   rmu0         cosinus de l'angle solaire zenithal
91!   nexca        pas de temps couplage
92!   nisurf       index de la surface a traiter (1 = sol continental)
93!   knon         nombre de points de la surface a traiter
94!   knindex      index des points de la surface a traiter
95!   pctsrf       tableau des pourcentages de surface de chaque maille
96!   rlon         longitudes
97!   rlat         latitudes
98!   cufi,cvfi    resolution des mailles en x et y (m)
99!   debut        logical: 1er appel a la physique
100!   lafin        logical: dernier appel a la physique
101!   ok_veget     logical: appel ou non au schema de surface continental
102!                     (si false calcul simplifie des fluxs sur les continents)
103!   zlev         hauteur de la premiere couche
104!   u1_lay       vitesse u 1ere couche
105!   v1_lay       vitesse v 1ere couche
106!   temp_air     temperature de l'air 1ere couche
107!   spechum      humidite specifique 1ere couche
108!   epot_air     temp potentielle de l'air
109!   ccanopy      concentration CO2 canopee
110!   tq_cdrag     cdrag
111!   petAcoef     coeff. A de la resolution de la CL pour t
112!   peqAcoef     coeff. A de la resolution de la CL pour q
113!   petBcoef     coeff. B de la resolution de la CL pour t
114!   peqBcoef     coeff. B de la resolution de la CL pour q
115!   precip_rain  precipitation liquide
116!   precip_snow  precipitation solide
117!   sollw        flux IR net a la surface
118!   sollwdown    flux IR descendant a la surface
119!   swnet        flux solaire net
120!   swdown       flux solaire entrant a la surface
121!   albedo       albedo de la surface
122!   tsurf        temperature de surface
123!   p1lay        pression 1er niveau (milieu de couche)
124!   ps           pression au sol
125!   radsol       rayonnement net aus sol (LW + SW)
126!   ocean        type d'ocean utilise (force, slab, couple)
127!   fder         derivee des flux (pour le couplage)
128!   taux, tauy   tension de vents
129!   rugos        rugosite
130!   zmasq        masque terre/ocean
131!   rugoro       rugosite orographique
132!   run_off_lic_0 runoff glacier du pas de temps precedent
133!
134! output:
135!   evap         evaporation totale
136!   fluxsens     flux de chaleur sensible
137!   fluxlat      flux de chaleur latente
138!   tsol_rad     
139!   tsurf_new    temperature au sol
140!   alb_new      albedo
141!   emis_new     emissivite
142!   z0_new       surface roughness
143!   pctsrf_new   nouvelle repartition des surfaces
144
145#include "iniprint.h"
146
147
148! Parametres d'entree
149  integer, intent(IN) :: itime
150  integer, intent(IN) :: iim, jjm
151  integer, intent(IN) :: klon
152  real, intent(IN) :: dtime
153  real, intent(IN) :: date0
154  integer, intent(IN) :: jour
155  real, intent(IN)    :: rmu0(klon)
156  integer, intent(IN) :: nisurf
157  integer, intent(IN) :: knon
158  integer, dimension(klon), intent(in) :: knindex
159  real, dimension(klon,nbsrf), intent(IN) :: pctsrf
160  logical, intent(IN) :: debut, lafin, ok_veget
161  real, dimension(klon), intent(IN) :: rlon, rlat
162  real, dimension(klon), intent(IN) :: cufi, cvfi
163  real, dimension(klon), intent(INOUT) :: tq_cdrag
164  real, dimension(klon), intent(IN) :: zlev
165  real, dimension(klon), intent(IN) :: u1_lay, v1_lay
166  real, dimension(klon), intent(IN) :: temp_air, spechum
167  real, dimension(klon), intent(IN) :: epot_air, ccanopy
168  real, dimension(klon), intent(IN) :: petAcoef, peqAcoef
169  real, dimension(klon), intent(IN) :: petBcoef, peqBcoef
170  real, dimension(klon), intent(IN) :: precip_rain, precip_snow
171  real, dimension(klon), intent(IN) :: sollw, sollwdown, swnet, swdown
172  real, dimension(klon), intent(IN) :: ps, albedo
173  real, dimension(klon), intent(IN) :: tsurf, p1lay
174  REAL, DIMENSION(klon), INTENT(INOUT) :: radsol,fder
175  real, dimension(klon), intent(IN) :: zmasq
176  real, dimension(klon), intent(IN) :: taux, tauy, rugos, rugoro
177  character (len = 6)  :: ocean
178  integer              :: npas, nexca ! nombre et pas de temps couplage
179  real, dimension(klon), intent(INOUT) :: evap, snow, qsurf
180!! PB ajout pour soil
181  logical          :: soil_model
182  integer          :: nsoilmx
183  REAL, DIMENSION(klon, nsoilmx) :: tsoil
184  REAL, dimension(klon), intent(INOUT) :: qsol
185  REAL, dimension(klon)          :: soilcap
186  REAL, dimension(klon)          :: soilflux
187! Parametres de sortie
188  real, dimension(klon), intent(OUT):: fluxsens, fluxlat
189  real, dimension(klon), intent(OUT):: tsol_rad, tsurf_new, alb_new
190  real, dimension(klon), intent(OUT):: alblw
191  real, dimension(klon), intent(OUT):: emis_new, z0_new
192  real, dimension(klon), intent(OUT):: dflux_l, dflux_s
193  real, dimension(klon,nbsrf), intent(OUT) :: pctsrf_new
194  real, dimension(klon), intent(INOUT):: agesno
195  real, dimension(klon), intent(INOUT):: run_off_lic_0
196
197! Flux thermique utiliser pour fondre la neige
198!jld a rajouter   real, dimension(klon), intent(INOUT):: ffonte
199  real, dimension(klon), intent(INOUT):: ffonte
200! Flux d'eau "perdue" par la surface et nécessaire pour que limiter la
201! hauteur de neige, en kg/m2/s
202!jld a rajouter   real, dimension(klon), intent(INOUT):: fqcalving
203  real, dimension(klon), intent(INOUT):: fqcalving
204
205! Local
206  character (len = 20),save :: modname = 'interfsurf_hq'
207  character (len = 80) :: abort_message
208  logical, save        :: first_call = .true.
209  integer, save        :: error
210  integer              :: ii, index
211  logical,save              :: check = .false.
212  real, dimension(klon):: cal, beta, dif_grnd, capsol
213!!$PB  real, parameter      :: calice=1.0/(5.1444e+06*0.15), tau_gl=86400.*5.
214  real, parameter      :: calice=1.0/(5.1444e+06*0.15), tau_gl=86400.*5.
215  real, parameter      :: calsno=1./(2.3867e+06*.15)
216  real, dimension(klon):: alb_ice
217  real, dimension(klon):: tsurf_temp
218  real, dimension(klon):: qsurf_new
219!!  real, allocatable, dimension(:), save :: alb_neig_grid
220  real, dimension(klon):: alb_neig, alb_eau
221  real, DIMENSION(klon):: zfra
222  logical              :: cumul = .false.
223  INTEGER,dimension(1) :: iloc
224  INTEGER                 :: isize
225  real, dimension(klon):: fder_prev
226  REAL, dimension(klon) :: bidule
227
228  if (check) write(*,*) 'Entree ', modname
229!
230! On doit commencer par appeler les schemas de surfaces continentales
231! car l'ocean a besoin du ruissellement qui est y calcule
232!
233  if (first_call) then
234    call conf_interface(tau_calv)
235    if (nisurf /= is_ter .and. klon > 1) then
236      write(*,*)' *** Warning ***'
237      write(*,*)' nisurf = ',nisurf,' /= is_ter = ',is_ter
238      write(*,*)'or on doit commencer par les surfaces continentales'
239      abort_message='voir ci-dessus'
240      call abort_gcm(modname,abort_message,1)
241    endif
242    if (ocean /= 'slab  ' .and. ocean /= 'force ' .and. ocean /= 'couple') then
243      write(*,*)' *** Warning ***'
244      write(*,*)'Option couplage pour l''ocean = ', ocean
245      abort_message='option pour l''ocean non valable'
246      call abort_gcm(modname,abort_message,1)
247    endif
248    if ( is_oce > is_sic ) then
249      write(*,*)' *** Warning ***'
250      write(*,*)' Pour des raisons de sequencement dans le code'
251      write(*,*)' l''ocean doit etre traite avant la banquise'
252      write(*,*)' or is_oce = ',is_oce, '> is_sic = ',is_sic
253      abort_message='voir ci-dessus'
254      call abort_gcm(modname,abort_message,1)
255    endif
256!    allocate(alb_neig_grid(klon), stat = error)
257!    if (error /= 0) then
258!      abort_message='Pb allocation alb_neig_grid'
259!      call abort_gcm(modname,abort_message,1)
260!    endif
261  endif
262  first_call = .false.
263 
264! Initialisations diverses
265!
266!!$  cal=0.; beta=1.; dif_grnd=0.; capsol=0.
267!!$  alb_new = 0.; z0_new = 0.; alb_neig = 0.0
268!!$! PB
269!!$  tsurf_new = 0.
270
271!IM cf JLD
272  ffonte(1:knon)=0.
273  fqcalving(1:knon)=0.
274
275  cal = 999999. ; beta = 999999. ; dif_grnd = 999999. ; capsol = 999999.
276  alb_new = 999999. ; z0_new = 999999. ; alb_neig = 999999.
277  tsurf_new = 999999.
278  alblw = 999999.
279! Aiguillage vers les differents schemas de surface
280
281  if (nisurf == is_ter) then
282!
283! Surface "terre" appel a l'interface avec les sols continentaux
284!
285! allocation du run-off
286    if (.not. allocated(coastalflow)) then
287      allocate(coastalflow(knon), stat = error)
288      if (error /= 0) then
289        abort_message='Pb allocation coastalflow'
290        call abort_gcm(modname,abort_message,1)
291      endif
292      allocate(riverflow(knon), stat = error)
293      if (error /= 0) then
294        abort_message='Pb allocation riverflow'
295        call abort_gcm(modname,abort_message,1)
296      endif
297      allocate(run_off(knon), stat = error)
298      if (error /= 0) then
299        abort_message='Pb allocation run_off'
300        call abort_gcm(modname,abort_message,1)
301      endif
302!!$PB
303      ALLOCATE (tmp_rriv(iim,jjm+1), stat=error)
304      if (error /= 0) then
305        abort_message='Pb allocation tmp_rriv'
306        call abort_gcm(modname,abort_message,1)
307      endif
308      ALLOCATE (tmp_rcoa(iim,jjm+1), stat=error)
309      if (error /= 0) then
310        abort_message='Pb allocation tmp_rcoa'
311        call abort_gcm(modname,abort_message,1)
312      endif
313      ALLOCATE (tmp_rlic(iim,jjm+1), stat=error)
314      if (error /= 0) then
315        abort_message='Pb allocation tmp_rlic'
316        call abort_gcm(modname,abort_message,1)
317      endif
318
319!!$
320    else if (size(coastalflow) /= knon) then
321      write(*,*)'Bizarre, le nombre de points continentaux'
322      write(*,*)'a change entre deux appels. J''arrete ...'
323      abort_message='voir ci-dessus'
324      call abort_gcm(modname,abort_message,1)
325    endif
326    coastalflow = 0.
327    riverflow = 0.
328!
329! Calcul age de la neige
330!
331!!$ PB ATTENTION changement ordre des appels
332!!$    CALL albsno(klon,agesno,alb_neig_grid) 
333
334    if (.not. ok_veget) then
335!
336! calcul albedo: lecture albedo fichier CL puis ajout albedo neige
337!
338       call interfsur_lim(itime, dtime, jour, &
339     & klon, nisurf, knon, knindex, debut,  &
340     & alb_new, z0_new)
341
342! calcul snow et qsurf, hydrol adapté
343!
344       CALL calbeta(dtime, nisurf, knon, snow, qsol, beta, capsol, dif_grnd)
345
346       IF (soil_model) THEN
347           CALL soil(dtime, nisurf, knon,snow, tsurf, tsoil,soilcap, soilflux)
348           cal(1:knon) = RCPD / soilcap(1:knon)
349           radsol(1:knon) = radsol(1:knon)  + soilflux(1:knon)
350       ELSE
351           cal = RCPD * capsol
352!!$      cal = capsol
353       ENDIF
354       CALL calcul_fluxs( klon, knon, nisurf, dtime, &
355     &   tsurf, p1lay, cal, beta, tq_cdrag, ps, &
356     &   precip_rain, precip_snow, snow, qsurf,  &
357     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
358     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
359     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
360
361       CALL fonte_neige( klon, knon, nisurf, dtime, &
362     &   tsurf, p1lay, cal, beta, tq_cdrag, ps, &
363     &   precip_rain, precip_snow, snow, qsol,  &
364     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
365     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
366     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
367     &   fqcalving,ffonte, run_off_lic_0)
368
369
370     call albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)) 
371     where (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.
372     zfra(1:knon) = max(0.0,min(1.0,snow(1:knon)/(snow(1:knon)+10.0)))
373     alb_new(1 : knon)  = alb_neig(1 : knon) *zfra(1:knon) + &
374    &                     alb_new(1 : knon)*(1.0-zfra(1:knon))
375     z0_new = sqrt(z0_new**2+rugoro**2)
376     alblw(1 : knon) = alb_new(1 : knon)
377
378    else
379!!      CALL albsno(klon,agesno,alb_neig_grid) 
380!
381!  appel a sechiba
382!
383#ifdef CPP_VEGET
384      call interfsol(itime, klon, dtime, date0, nisurf, knon, &
385     &  knindex, rlon, rlat, cufi, cvfi, iim, jjm, pctsrf, &
386     &  debut, lafin, ok_veget, &
387     &  zlev,  u1_lay, v1_lay, temp_air, spechum, epot_air, ccanopy, &
388     &  tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
389     &  precip_rain, precip_snow, sollwdown, swnet, swdown, &
390     &  tsurf, p1lay/100., ps/100., radsol, &
391     &  evap, fluxsens, fluxlat, &             
392     &  tsol_rad, tsurf_new, alb_new, alblw, &
393     &  emis_new, z0_new, dflux_l, dflux_s, qsurf_new)
394
395
396! ajout de la contribution du relief
397
398      z0_new = SQRT(z0_new**2+rugoro**2)
399!
400! mise a jour de l'humidite saturante calculee par ORCHIDEE
401      qsurf(1:knon) = qsurf_new(1:knon)
402#endif
403
404    endif   
405!
406! Remplissage des pourcentages de surface
407!
408    pctsrf_new(:,nisurf) = pctsrf(:,nisurf)
409
410  else if (nisurf == is_oce) then
411
412    if (check) write(*,*)'ocean, nisurf = ',nisurf
413
414
415!
416! Surface "ocean" appel a l'interface avec l'ocean
417!
418    if (ocean == 'couple') then
419      if (nexca == 0) then
420        abort_message='nexca = 0 dans interfoce_cpl'
421        call abort_gcm(modname,abort_message,1)
422      endif
423
424      cumul = .false.
425
426      iloc = maxloc(fder(1:klon))
427      if (check) then
428        if (fder(iloc(1))> 0.) then
429          WRITE(*,*)'**** Debug fder ****'
430          WRITE(*,*)'max fder(',iloc(1),') = ',fder(iloc(1))
431        endif
432      endif
433!!$
434!!$      where(fder.gt.0.)
435!!$        fder = 0.
436!!$      endwhere
437
438      call interfoce(itime, dtime, cumul, &
439      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
440      & ocean, npas, nexca, debut, lafin, &
441      & swdown, sollw, precip_rain, precip_snow, evap, tsurf, &
442      & fluxlat, fluxsens, fder, albedo, taux, tauy, zmasq, &
443      & tsurf_new, alb_new, pctsrf_new)
444
445!    else if (ocean == 'slab  ') then
446!      call interfoce(nisurf)
447    else                              ! lecture conditions limites
448      call interfoce(itime, dtime, jour, &
449     &  klon, nisurf, knon, knindex, &
450     &  debut, &
451     &  tsurf_new, pctsrf_new)
452
453    endif
454
455    tsurf_temp = tsurf_new
456    cal = 0.
457    beta = 1.
458    dif_grnd = 0.
459    alb_neig(:) = 0.
460    agesno(:) = 0.
461
462    call calcul_fluxs( klon, knon, nisurf, dtime, &
463     &   tsurf_temp, p1lay, cal, beta, tq_cdrag, ps, &
464     &   precip_rain, precip_snow, snow, qsurf,  &
465     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
466     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
467     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
468
469    fder_prev = fder   
470    fder = fder_prev + dflux_s + dflux_l
471
472      iloc = maxloc(fder(1:klon))
473        if (check.and.fder(iloc(1))> 0.) then
474          WRITE(*,*)'**** Debug fder****'
475          WRITE(*,*)'max fder(',iloc(1),') = ',fder(iloc(1))
476          WRITE(*,*)'fder_prev, dflux_s, dflux_l',fder_prev(iloc(1)), &
477     &                        dflux_s(iloc(1)), dflux_l(iloc(1))
478        endif
479!!$
480!!$      where(fder.gt.0.)
481!!$        fder = 0.
482!!$      endwhere
483
484!
485! 2eme appel a interfoce pour le cumul des champs (en particulier
486! fluxsens et fluxlat calcules dans calcul_fluxs)
487!
488    if (ocean == 'couple') then
489
490      cumul = .true.
491
492      call interfoce(itime, dtime, cumul, &
493      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
494      & ocean, npas, nexca, debut, lafin, &
495      & swdown, sollw, precip_rain, precip_snow, evap, tsurf, &
496      & fluxlat, fluxsens, fder, albedo, taux, tauy, zmasq, &
497      & tsurf_new, alb_new, pctsrf_new)
498
499!    else if (ocean == 'slab  ') then
500!      call interfoce(nisurf)
501
502    endif
503
504!
505! calcul albedo
506!
507
508    if ( minval(rmu0) == maxval(rmu0) .and. minval(rmu0) == -999.999 ) then
509      CALL alboc(FLOAT(jour),rlat,alb_eau)
510    else  ! cycle diurne
511      CALL alboc_cd(rmu0,alb_eau)
512    endif
513    DO ii =1, knon
514      alb_new(ii) = alb_eau(knindex(ii))
515    enddo
516
517    z0_new = sqrt(rugos**2 + rugoro**2)
518    alblw(1:knon) = alb_new(1:knon)
519
520!
521  else if (nisurf == is_sic) then
522
523    if (check) write(*,*)'sea ice, nisurf = ',nisurf
524
525!
526! Surface "glace de mer" appel a l'interface avec l'ocean
527!
528!
529    if (ocean == 'couple') then
530
531      cumul =.false.
532
533      iloc = maxloc(fder(1:klon))
534      if (check.and.fder(iloc(1))> 0.) then
535        WRITE(*,*)'**** Debug fder ****'
536        WRITE(*,*)'max fder(',iloc(1),') = ',fder(iloc(1))
537      endif
538!!$
539!!$      where(fder.gt.0.)
540!!$        fder = 0.
541!!$      endwhere
542
543      call interfoce(itime, dtime, cumul, &
544      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
545      & ocean, npas, nexca, debut, lafin, &
546      & swdown, sollw, precip_rain, precip_snow, evap, tsurf, &
547      & fluxlat, fluxsens, fder, albedo, taux, tauy, zmasq, &
548      & tsurf_new, alb_new, pctsrf_new)
549
550      tsurf_temp = tsurf_new
551      cal = 0.
552      dif_grnd = 0.
553      beta = 1.0
554
555!    else if (ocean == 'slab  ') then
556!      call interfoce(nisurf)
557    ELSE
558!                              ! lecture conditions limites
559      CALL interfoce(itime, dtime, jour, &
560             &  klon, nisurf, knon, knindex, &
561             &  debut, &
562             &  tsurf_new, pctsrf_new)
563
564!IM cf LF
565      DO ii = 1, knon
566       IF (pctsrf_new(ii,nisurf) < EPSFRA) then
567          snow(ii) = 0.0
568!IM cf LF/JLD         tsurf(ii) = RTT - 1.8
569          tsurf_new(ii) = RTT - 1.8
570          IF (soil_model) tsoil(ii,:) = RTT -1.8
571        endif
572      enddo
573
574      CALL calbeta(dtime, nisurf, knon, snow, qsol, beta, capsol, dif_grnd)
575     
576      IF (soil_model) THEN
577!IM cf LF/JLD        CALL soil(dtime, nisurf, knon,snow, tsurf, tsoil,soilcap, soilflux)
578         CALL soil(dtime, nisurf, knon,snow, tsurf_new, tsoil,soilcap, soilflux)
579         cal(1:knon) = RCPD / soilcap(1:knon)
580         radsol(1:knon) = radsol(1:knon)  + soilflux(1:knon)
581         dif_grnd = 0.
582      ELSE
583         dif_grnd = 1.0 / tau_gl
584         cal = RCPD * calice
585         WHERE (snow > 0.0) cal = RCPD * calsno
586      ENDIF
587      tsurf_temp = tsurf
588      beta = 1.0
589    ENDIF
590
591    CALL calcul_fluxs( klon, knon, nisurf, dtime, &
592         &   tsurf_temp, p1lay, cal, beta, tq_cdrag, ps, &
593         &   precip_rain, precip_snow, snow, qsurf,  &
594         &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
595         &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
596         &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
597
598    IF (ocean /= 'couple') THEN
599      CALL fonte_neige( klon, knon, nisurf, dtime, &
600             &   tsurf_temp, p1lay, cal, beta, tq_cdrag, ps, &
601             &   precip_rain, precip_snow, snow, qsol,  &
602             &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
603             &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
604             &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
605             &   fqcalving,ffonte, run_off_lic_0)
606
607!     calcul albedo
608
609      CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)) 
610      WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.
611      zfra(1:knon) = MAX(0.0,MIN(1.0,snow(1:knon)/(snow(1:knon)+10.0)))
612      alb_new(1 : knon) = alb_neig(1 : knon) *zfra(1:knon) + &
613     &                    0.6 * (1.0-zfra(1:knon))
614!!      alb_new(1 : knon) = 0.6
615    ENDIF
616
617    fder_prev = fder   
618    fder = fder_prev + dflux_s + dflux_l
619
620      iloc = maxloc(fder(1:klon))
621      if (check.and.fder(iloc(1))> 0.) then
622        WRITE(*,*)'**** Debug fder ****'
623        WRITE(*,*)'max fder(',iloc(1),') = ',fder(iloc(1))
624        WRITE(*,*)'fder_prev, dflux_s, dflux_l',fder_prev(iloc(1)), &
625     &                        dflux_s(iloc(1)), dflux_l(iloc(1))
626      endif
627!!$      where(fder.gt.0.)
628!!$        fder = 0.
629!!$      endwhere
630
631!
632! 2eme appel a interfoce pour le cumul et le passage des flux a l'ocean
633!
634    if (ocean == 'couple') then
635
636      cumul =.true.
637
638      call interfoce(itime, dtime, cumul, &
639      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
640      & ocean, npas, nexca, debut, lafin, &
641      & swdown, sollw, precip_rain, precip_snow, evap, tsurf, &
642      & fluxlat, fluxsens, fder, albedo, taux, tauy, zmasq, &
643      & tsurf_new, alb_new, pctsrf_new)
644
645!    else if (ocean == 'slab  ') then
646!      call interfoce(nisurf)
647
648    endif
649
650     
651    z0_new = 0.002
652    z0_new = SQRT(z0_new**2+rugoro**2)
653    alblw(1:knon) = alb_new(1:knon)
654
655  else if (nisurf == is_lic) then
656
657    if (check) write(*,*)'glacier, nisurf = ',nisurf
658
659    if (.not. allocated(run_off_lic)) then
660      allocate(run_off_lic(knon), stat = error)
661      if (error /= 0) then
662        abort_message='Pb allocation run_off_lic'
663        call abort_gcm(modname,abort_message,1)
664      endif
665      run_off_lic = 0.
666    endif
667!
668! Surface "glacier continentaux" appel a l'interface avec le sol
669!
670!    call interfsol(nisurf)
671    IF (soil_model) THEN
672        CALL soil(dtime, nisurf, knon, snow, tsurf, tsoil,soilcap, soilflux)
673        cal(1:knon) = RCPD / soilcap(1:knon)
674        radsol(1:knon)  = radsol(1:knon) + soilflux(1:knon)
675    ELSE
676        cal = RCPD * calice
677        WHERE (snow > 0.0) cal = RCPD * calsno
678    ENDIF
679    beta = 1.0
680    dif_grnd = 0.0
681
682    call calcul_fluxs( klon, knon, nisurf, dtime, &
683     &   tsurf, p1lay, cal, beta, tq_cdrag, ps, &
684     &   precip_rain, precip_snow, snow, qsurf,  &
685     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
686     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
687     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
688
689    call fonte_neige( klon, knon, nisurf, dtime, &
690     &   tsurf, p1lay, cal, beta, tq_cdrag, ps, &
691     &   precip_rain, precip_snow, snow, qsol,  &
692     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
693     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
694     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
695     &   fqcalving,ffonte, run_off_lic_0)
696
697! passage du run-off des glaciers calcule dans fonte_neige au coupleur
698    bidule=0.
699    bidule(1:knon)= run_off_lic(1:knon)   
700    call gath2cpl(bidule, tmp_rlic, klon, knon,iim,jjm,knindex)
701!
702! calcul albedo
703!
704     CALL albsno(klon,knon,dtime,agesno(:),alb_neig(:), precip_snow(:)) 
705     WHERE (snow(1 : knon) .LT. 0.0001) agesno(1 : knon) = 0.
706     zfra(1:knon) = MAX(0.0,MIN(1.0,snow(1:knon)/(snow(1:knon)+10.0)))
707     alb_new(1 : knon)  = alb_neig(1 : knon)*zfra(1:knon) + &
708    &                     0.6 * (1.0-zfra(1:knon))
709!IM cf FH/GK     alb_new(1 : knon)  = 0.6
710!       alb_new(1 : knon)  = 0.82
711!IM cf JLD/ GK
712!IM: 211003 Ksta0.77      alb_new(1 : knon)  = 0.77
713!IM: KstaTER0.8 & LMD_ARMIP5    alb_new(1 : knon)  = 0.8
714!IM: KstaTER0.77 & LMD_ARMIP6   
715        alb_new(1 : knon)  = 0.77
716
717!
718! Rugosite
719!
720    z0_new = rugoro
721!
722! Remplissage des pourcentages de surface
723!
724    pctsrf_new(:,nisurf) = pctsrf(:,nisurf)
725
726    alblw(1:knon) = alb_new(1:knon)
727  else
728    write(*,*)'Index surface = ',nisurf
729    abort_message = 'Index surface non valable'
730    call abort_gcm(modname,abort_message,1)
731  endif
732
733  END SUBROUTINE interfsurf_hq
734
735!
736!#########################################################################
737!
738  SUBROUTINE interfsurf_vent(nisurf, knon   &         
739  &                     )
740!
741! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
742! (sols continentaux, oceans, glaces) pour les tensions de vents.
743! En pratique l'interface se fait entre la couche limite du modele
744! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
745!
746!
747! L.Fairhead 02/2000
748!
749! input:
750!   nisurf       index de la surface a traiter (1 = sol continental)
751!   knon         nombre de points de la surface a traiter
752
753! Parametres d'entree
754  integer, intent(IN) :: nisurf
755  integer, intent(IN) :: knon
756
757
758  return
759  END SUBROUTINE interfsurf_vent
760!
761!#########################################################################
762!
763#ifdef CPP_VEGET
764  SUBROUTINE interfsol(itime, klon, dtime, date0, nisurf, knon, &
765     & knindex, rlon, rlat, cufi, cvfi, iim, jjm, pctsrf, &
766     & debut, lafin, ok_veget, &
767     & plev,  u1_lay, v1_lay, temp_air, spechum, epot_air, ccanopy, &
768     & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
769     & precip_rain, precip_snow, lwdown, swnet, swdown, &
770     & tsurf, p1lay, ps, radsol, &
771     & evap, fluxsens, fluxlat, &             
772     & tsol_rad, tsurf_new, alb_new, alblw, &
773     & emis_new, z0_new, dflux_l, dflux_s, qsurf)
774
775  USE intersurf
776
777! Cette routine sert d'interface entre le modele atmospherique et le
778! modele de sol continental. Appel a sechiba
779!
780! L. Fairhead 02/2000
781!
782! input:
783!   itime        numero du pas de temps
784!   klon         nombre total de points de grille
785!   dtime        pas de temps de la physique (en s)
786!   nisurf       index de la surface a traiter (1 = sol continental)
787!   knon         nombre de points de la surface a traiter
788!   knindex      index des points de la surface a traiter
789!   rlon         longitudes de la grille entiere
790!   rlat         latitudes de la grille entiere
791!   pctsrf       tableau des fractions de surface de chaque maille
792!   debut        logical: 1er appel a la physique (lire les restart)
793!   lafin        logical: dernier appel a la physique (ecrire les restart)
794!   ok_veget     logical: appel ou non au schema de surface continental
795!                     (si false calcul simplifie des fluxs sur les continents)
796!   plev         hauteur de la premiere couche (Pa)     
797!   u1_lay       vitesse u 1ere couche
798!   v1_lay       vitesse v 1ere couche
799!   temp_air     temperature de l'air 1ere couche
800!   spechum      humidite specifique 1ere couche
801!   epot_air     temp pot de l'air
802!   ccanopy      concentration CO2 canopee
803!   tq_cdrag     cdrag
804!   petAcoef     coeff. A de la resolution de la CL pour t
805!   peqAcoef     coeff. A de la resolution de la CL pour q
806!   petBcoef     coeff. B de la resolution de la CL pour t
807!   peqBcoef     coeff. B de la resolution de la CL pour q
808!   precip_rain  precipitation liquide
809!   precip_snow  precipitation solide
810!   lwdown       flux IR descendant a la surface
811!   swnet        flux solaire net
812!   swdown       flux solaire entrant a la surface
813!   tsurf        temperature de surface
814!   p1lay        pression 1er niveau (milieu de couche)
815!   ps           pression au sol
816!   radsol       rayonnement net aus sol (LW + SW)
817!   
818!
819! input/output
820!   run_off      ruissellement total
821!
822! output:
823!   evap         evaporation totale
824!   fluxsens     flux de chaleur sensible
825!   fluxlat      flux de chaleur latente
826!   tsol_rad     
827!   tsurf_new    temperature au sol
828!   alb_new      albedo
829!   emis_new     emissivite
830!   z0_new       surface roughness
831!   qsurf        air moisture at surface
832
833! Parametres d'entree
834  integer, intent(IN) :: itime
835  integer, intent(IN) :: klon
836  real, intent(IN)    :: dtime
837  real, intent(IN)    :: date0
838  integer, intent(IN) :: nisurf
839  integer, intent(IN) :: knon
840  integer, intent(IN) :: iim, jjm
841  integer, dimension(klon), intent(IN) :: knindex
842  logical, intent(IN) :: debut, lafin, ok_veget
843  real, dimension(klon,nbsrf), intent(IN) :: pctsrf
844  real, dimension(klon), intent(IN) :: rlon, rlat
845  real, dimension(klon), intent(IN) :: cufi, cvfi
846  real, dimension(klon), intent(IN) :: plev
847  real, dimension(klon), intent(IN) :: u1_lay, v1_lay
848  real, dimension(klon), intent(IN) :: temp_air, spechum
849  real, dimension(klon), intent(IN) :: epot_air, ccanopy
850  real, dimension(klon), intent(INOUT) :: tq_cdrag
851  real, dimension(klon), intent(IN) :: petAcoef, peqAcoef
852  real, dimension(klon), intent(IN) :: petBcoef, peqBcoef
853  real, dimension(klon), intent(IN) :: precip_rain, precip_snow
854  real, dimension(klon), intent(IN) :: lwdown, swnet, swdown, ps
855!IM cf. JP +++
856  real, dimension(klon) :: swdown_vrai
857!IM cf. JP ---
858  real, dimension(klon), intent(IN) :: tsurf, p1lay
859  real, dimension(klon), intent(IN) :: radsol
860! Parametres de sortie
861  real, dimension(klon), intent(OUT):: evap, fluxsens, fluxlat, qsurf
862  real, dimension(klon), intent(OUT):: tsol_rad, tsurf_new, alb_new, alblw
863  real, dimension(klon), intent(OUT):: emis_new, z0_new
864  real, dimension(klon), intent(OUT):: dflux_s, dflux_l
865
866! Local
867!
868  integer              :: ii, ij, jj, igrid, ireal, i, index, iglob
869  integer              :: error
870  character (len = 20) :: modname = 'interfsol'
871  character (len = 80) :: abort_message
872  logical,save              :: check = .FALSE.
873  real, dimension(klon) :: cal, beta, dif_grnd, capsol
874! type de couplage dans sechiba
875!  character (len=10)   :: coupling = 'implicit'
876! drapeaux controlant les appels dans SECHIBA
877!  type(control_type), save   :: control_in
878! Preserved albedo
879!IM cf. JP +++
880  real, allocatable, dimension(:), save :: albedo_keep, zlev
881!IM cf. JP ---
882! coordonnees geographiques
883  real, allocatable, dimension(:,:), save :: lalo
884! pts voisins
885  integer,allocatable, dimension(:,:), save :: neighbours
886! fractions continents
887  real,allocatable, dimension(:), save :: contfrac
888! resolution de la grille
889  real, allocatable, dimension (:,:), save :: resolution
890! correspondance point n -> indices (i,j)
891  integer, allocatable, dimension(:,:), save :: correspond
892! offset pour calculer les point voisins
893  integer, dimension(8,3), save :: off_ini
894  integer, dimension(8), save :: offset
895! Identifieurs des fichiers restart et histoire
896  integer, save          :: rest_id, hist_id
897  integer, save          :: rest_id_stom, hist_id_stom
898!
899  real, allocatable, dimension (:,:), save :: lon_scat, lat_scat 
900
901  logical, save          :: lrestart_read = .true. , lrestart_write = .false.
902
903  real, dimension(klon):: snow
904  real, dimension(knon,2) :: albedo_out
905! Pb de nomenclature
906  real, dimension(klon) :: petA_orc, peqA_orc
907  real, dimension(klon) :: petB_orc, peqB_orc
908! Pb de correspondances de grilles
909  integer, dimension(:), save, allocatable :: ig, jg
910  integer :: indi, indj
911  integer, dimension(klon) :: ktindex
912  REAL, dimension(klon) :: bidule
913! Essai cdrag
914  real, dimension(klon) :: cdrag
915
916#include "temps.inc"
917#include "YOMCST.inc"
918#include "iniprint.h"
919
920  if (check) write(lunout,*)'Entree ', modname
921  if (check) write(lunout,*)'ok_veget = ',ok_veget
922
923  ktindex(:) = knindex(:) + iim - 1
924
925! initialisation
926  if (debut) then
927
928  IF ( .NOT. allocated(albedo_keep)) THEN
929     ALLOCATE(albedo_keep(klon))
930     ALLOCATE(zlev(klon))
931  ENDIF
932! Pb de correspondances de grilles
933   allocate(ig(klon))
934   allocate(jg(klon))
935   ig(1) = 1
936   jg(1) = 1
937   indi = 0
938   indj = 2
939   do igrid = 2, klon - 1
940     indi = indi + 1
941     if ( indi > iim) then
942       indi = 1
943       indj = indj + 1
944     endif
945     ig(igrid) = indi
946     jg(igrid) = indj
947   enddo
948   ig(klon) = 1
949   jg(klon) = jjm + 1
950!
951!  Initialisation des offset   
952!
953! offset bord ouest
954   off_ini(1,1) = - iim  ; off_ini(2,1) = - iim + 1; off_ini(3,1) = 1
955   off_ini(4,1) = iim + 1; off_ini(5,1) = iim      ; off_ini(6,1) = 2 * iim - 1
956   off_ini(7,1) = iim -1 ; off_ini(8,1) = - 1
957! offset point normal
958   off_ini(1,2) = - iim  ; off_ini(2,2) = - iim + 1; off_ini(3,2) = 1
959   off_ini(4,2) = iim + 1; off_ini(5,2) = iim      ; off_ini(6,2) = iim - 1
960   off_ini(7,2) = -1     ; off_ini(8,2) = - iim - 1
961! offset bord   est
962   off_ini(1,3) = - iim; off_ini(2,3) = - 2 * iim + 1; off_ini(3,3) = - iim + 1
963   off_ini(4,3) =  1   ; off_ini(5,3) = iim          ; off_ini(6,3) = iim - 1
964   off_ini(7,3) = -1   ; off_ini(8,3) = - iim - 1
965!
966! Initialisation des correspondances point -> indices i,j
967!
968    if (( .not. allocated(correspond))) then
969      allocate(correspond(iim,jjm+1), stat = error)
970      if (error /= 0) then
971        abort_message='Pb allocation correspond'
972        call abort_gcm(modname,abort_message,1)
973      endif     
974    endif
975!
976! Attention aux poles
977!
978    do igrid = 1, knon
979      index = ktindex(igrid)
980          jj = int((index - 1)/iim) + 1
981          ij = index - (jj - 1) * iim
982      correspond(ij,jj) = igrid
983    enddo
984
985! Allouer et initialiser le tableau de coordonnees du sol
986!
987    if ((.not. allocated(lalo))) then
988      allocate(lalo(knon,2), stat = error)
989      if (error /= 0) then
990        abort_message='Pb allocation lalo'
991        call abort_gcm(modname,abort_message,1)
992      endif     
993    endif
994    if ((.not. allocated(lon_scat))) then
995      allocate(lon_scat(iim,jjm+1), stat = error)
996      if (error /= 0) then
997        abort_message='Pb allocation lon_scat'
998        call abort_gcm(modname,abort_message,1)
999      endif     
1000    endif
1001    if ((.not. allocated(lat_scat))) then
1002      allocate(lat_scat(iim,jjm+1), stat = error)
1003      if (error /= 0) then
1004        abort_message='Pb allocation lat_scat'
1005        call abort_gcm(modname,abort_message,1)
1006      endif     
1007    endif
1008    lon_scat = 0.
1009    lat_scat = 0.
1010    do igrid = 1, knon
1011      index = knindex(igrid)
1012      lalo(igrid,2) = rlon(index)
1013      lalo(igrid,1) = rlat(index)
1014      ij = index - int((index-1)/iim)*iim - 1
1015      jj = 2 + int((index-1)/iim)
1016      if (mod(index,iim) == 1 ) then
1017        jj = 1 + int((index-1)/iim)
1018        ij = iim
1019      endif
1020!      lon_scat(ij,jj) = rlon(index)
1021!      lat_scat(ij,jj) = rlat(index)
1022    enddo
1023    index = 1
1024    do jj = 2, jjm
1025      do ij = 1, iim
1026        index = index + 1
1027        lon_scat(ij,jj) = rlon(index)
1028        lat_scat(ij,jj) = rlat(index)
1029      enddo
1030    enddo
1031    lon_scat(:,1) = lon_scat(:,2)
1032    lat_scat(:,1) = rlat(1)
1033    lon_scat(:,jjm+1) = lon_scat(:,2)
1034    lat_scat(:,jjm+1) = rlat(klon)
1035! Pb de correspondances de grilles!
1036!    do igrid = 1, knon
1037!      index = ktindex(igrid)
1038!      ij = ig(index)
1039!      jj = jg(index)
1040!      lon_scat(ij,jj) = rlon(index)
1041!      lat_scat(ij,jj) = rlat(index)
1042!    enddo
1043
1044!
1045! Allouer et initialiser le tableau des voisins et des fraction de continents
1046!
1047    if ( (.not.allocated(neighbours))) THEN
1048      allocate(neighbours(knon,8), stat = error)
1049      if (error /= 0) then
1050        abort_message='Pb allocation neighbours'
1051        call abort_gcm(modname,abort_message,1)
1052      endif
1053    endif
1054    neighbours = -1.
1055    if (( .not. allocated(contfrac))) then
1056      allocate(contfrac(knon), stat = error)
1057      if (error /= 0) then
1058        abort_message='Pb allocation contfrac'
1059        call abort_gcm(modname,abort_message,1)
1060      endif     
1061    endif
1062
1063    do igrid = 1, knon
1064      ireal = knindex(igrid)
1065      contfrac(igrid) = pctsrf(ireal,is_ter)
1066    enddo
1067
1068    do igrid = 1, knon
1069      iglob = ktindex(igrid)
1070      if (mod(iglob, iim) == 1) then
1071        offset = off_ini(:,1)
1072      else if(mod(iglob, iim) == 0) then
1073        offset = off_ini(:,3)
1074      else
1075        offset = off_ini(:,2)
1076      endif
1077      do i = 1, 8
1078        index = iglob + offset(i)
1079        ireal = (min(max(1, index - iim + 1), klon))
1080        if (pctsrf(ireal, is_ter) > EPSFRA) then
1081          jj = int((index - 1)/iim) + 1
1082          ij = index - (jj - 1) * iim
1083            neighbours(igrid, i) = correspond(ij, jj)
1084        endif
1085      enddo
1086    enddo
1087
1088!
1089!  Allocation et calcul resolutions
1090    IF ( (.NOT.ALLOCATED(resolution))) THEN
1091      ALLOCATE(resolution(knon,2), stat = error)
1092      if (error /= 0) then
1093        abort_message='Pb allocation resolution'
1094        call abort_gcm(modname,abort_message,1)
1095      endif
1096    ENDIF
1097    do igrid = 1, knon
1098      ij = knindex(igrid)
1099      resolution(igrid,1) = cufi(ij)
1100      resolution(igrid,2) = cvfi(ij)
1101    enddo 
1102
1103  endif                          ! (fin debut)
1104
1105!
1106! Appel a la routine sols continentaux
1107!
1108  if (lafin) lrestart_write = .true.
1109  if (check) write(lunout,*)'lafin ',lafin,lrestart_write
1110
1111  petA_orc = petBcoef * dtime
1112  petB_orc = petAcoef
1113  peqA_orc = peqBcoef * dtime
1114  peqB_orc = peqAcoef
1115
1116  cdrag = 0.
1117  cdrag(1:knon) = tq_cdrag(1:knon)
1118
1119!IM cf. JP +++
1120! zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/287.05*temp_air(1:knon))*9.80665)
1121  zlev(1:knon) = (100.*plev(1:knon))/((ps(1:knon)/RD*temp_air(1:knon))*RG)
1122!IM cf. JP ---
1123
1124
1125! PF et PASB
1126!   where(cdrag > 0.01)
1127!     cdrag = 0.01
1128!   endwhere
1129!  write(*,*)'Cdrag = ',minval(cdrag),maxval(cdrag)
1130
1131!
1132! Init Orchidee
1133!
1134  if (debut) then
1135    call intersurf_main (itime+itau_phy-1, iim, jjm+1, knon, ktindex, dtime, &
1136     & lrestart_read, lrestart_write, lalo, &
1137     & contfrac, neighbours, resolution, date0, &
1138     & zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
1139     & cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
1140     & precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
1141     & evap, fluxsens, fluxlat, coastalflow, riverflow, &
1142     & tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
1143     & lon_scat, lat_scat)
1144
1145!IM cf. JP +++
1146    albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
1147!IM cf. JP ---
1148
1149  endif
1150
1151!IM cf. JP +++
1152  swdown_vrai(1:knon) = swnet(1:knon)/(1. - albedo_keep(1:knon))
1153!IM cf. JP ---
1154
1155  call intersurf_main (itime+itau_phy, iim, jjm+1, knon, ktindex, dtime, &
1156     & lrestart_read, lrestart_write, lalo, &
1157     & contfrac, neighbours, resolution, date0, &
1158     & zlev,  u1_lay, v1_lay, spechum, temp_air, epot_air, ccanopy, &
1159     & cdrag, petA_orc, peqA_orc, petB_orc, peqB_orc, &
1160     & precip_rain, precip_snow, lwdown, swnet, swdown_vrai, ps, &
1161     & evap, fluxsens, fluxlat, coastalflow, riverflow, &
1162     & tsol_rad, tsurf_new, qsurf, albedo_out, emis_new, z0_new, &
1163     & lon_scat, lat_scat)
1164
1165!IM cf. JP +++
1166    albedo_keep(1:knon) = (albedo_out(1:knon,1)+albedo_out(1:knon,2))/2.
1167!IM cf. JP ---
1168
1169    bidule=0.
1170    bidule(1:knon)=riverflow(1:knon)
1171    call gath2cpl(bidule, tmp_rriv, klon, knon,iim,jjm,knindex)
1172    bidule=0.
1173    bidule(1:knon)=coastalflow(1:knon)
1174    call gath2cpl(bidule, tmp_rcoa, klon, knon,iim,jjm,knindex)
1175    alb_new(1:knon) = albedo_out(1:knon,1)
1176    alblw(1:knon) = albedo_out(1:knon,2)
1177
1178
1179! Convention orchidee: positif vers le haut
1180  fluxsens(1:knon) = -1. * fluxsens(1:knon)
1181  fluxlat(1:knon)  = -1. * fluxlat(1:knon)
1182
1183!  evap     = -1. * evap
1184
1185  if (debut) lrestart_read = .false.
1186
1187  END SUBROUTINE interfsol
1188#endif
1189!
1190!#########################################################################
1191!
1192  SUBROUTINE interfoce_cpl(itime, dtime, cumul, &
1193      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
1194      & ocean, npas, nexca, debut, lafin, &
1195      & swdown, lwdown, precip_rain, precip_snow, evap, tsurf, &
1196      & fluxlat, fluxsens, fder, albsol, taux, tauy, zmasq, &
1197      & tsurf_new, alb_new, pctsrf_new)
1198
1199! Cette routine sert d'interface entre le modele atmospherique et un
1200! coupleur avec un modele d'ocean 'complet' derriere
1201!
1202! Le modele de glace qu'il est prevu d'utiliser etant couple directement a
1203! l'ocean presentement, on va passer deux fois dans cette routine par pas de
1204! temps physique, une fois avec les points oceans et l'autre avec les points
1205! glace. A chaque pas de temps de couplage, la lecture des champs provenant
1206! du coupleur se fera "dans" l'ocean et l'ecriture des champs a envoyer
1207! au coupleur "dans" la glace. Il faut donc des tableaux de travail "tampons"
1208! dimensionnes sur toute la grille qui remplissent les champs sur les
1209! domaines ocean/glace quand il le faut. Il est aussi necessaire que l'index
1210! ocean soit traiter avant l'index glace (sinon tout intervertir)
1211!
1212!
1213! L. Fairhead 02/2000
1214!
1215! input:
1216!   itime        numero du pas de temps
1217!   iim, jjm     nbres de pts de grille
1218!   dtime        pas de temps de la physique
1219!   klon         nombre total de points de grille
1220!   nisurf       index de la surface a traiter (1 = sol continental)
1221!   pctsrf       tableau des fractions de surface de chaque maille
1222!   knon         nombre de points de la surface a traiter
1223!   knindex      index des points de la surface a traiter
1224!   rlon         longitudes
1225!   rlat         latitudes
1226!   debut        logical: 1er appel a la physique
1227!   lafin        logical: dernier appel a la physique
1228!   ocean        type d'ocean
1229!   nexca        frequence de couplage
1230!   swdown       flux solaire entrant a la surface
1231!   lwdown       flux IR net a la surface
1232!   precip_rain  precipitation liquide
1233!   precip_snow  precipitation solide
1234!   evap         evaporation
1235!   tsurf        temperature de surface
1236!   fder         derivee dF/dT
1237!   albsol       albedo du sol (coherent avec swdown)
1238!   taux         tension de vent en x
1239!   tauy         tension de vent en y
1240!   nexca        frequence de couplage
1241!   zmasq        masque terre/ocean
1242!
1243!
1244! output:
1245!   tsurf_new    temperature au sol
1246!   alb_new      albedo
1247!   pctsrf_new   nouvelle repartition des surfaces
1248!   alb_ice      albedo de la glace
1249!
1250#ifdef CPP_PSMILE 
1251  USE oasis
1252  integer :: il_time_secs !time in seconds
1253#endif
1254
1255! Parametres d'entree
1256  integer, intent(IN) :: itime
1257  integer, intent(IN) :: iim, jjm
1258  real, intent(IN) :: dtime
1259  integer, intent(IN) :: klon
1260  integer, intent(IN) :: nisurf
1261  integer, intent(IN) :: knon
1262  real, dimension(klon,nbsrf), intent(IN) :: pctsrf
1263  integer, dimension(klon), intent(in) :: knindex
1264  logical, intent(IN) :: debut, lafin
1265  real, dimension(klon), intent(IN) :: rlon, rlat
1266  character (len = 6)  :: ocean
1267  real, dimension(klon), intent(IN) :: lwdown, swdown
1268  real, dimension(klon), intent(IN) :: precip_rain, precip_snow
1269  real, dimension(klon), intent(IN) :: tsurf, fder, albsol, taux, tauy
1270  INTEGER              :: nexca, npas, kstep
1271  real, dimension(klon), intent(IN) :: zmasq
1272  real, dimension(klon), intent(IN) :: fluxlat, fluxsens
1273  logical, intent(IN)               :: cumul
1274  real, dimension(klon), intent(INOUT) :: evap
1275
1276! Parametres de sortie
1277  real, dimension(klon), intent(OUT):: tsurf_new, alb_new
1278  real, dimension(klon,nbsrf), intent(OUT) :: pctsrf_new
1279
1280! Variables locales
1281  integer                    :: j, error, sum_error, ig, cpl_index,i
1282  character (len = 20) :: modname = 'interfoce_cpl'
1283  character (len = 80) :: abort_message
1284  logical,save              :: check = .FALSE.
1285! variables pour moyenner les variables de couplage
1286  real, allocatable, dimension(:,:),save :: cpl_sols, cpl_nsol, cpl_rain
1287  real, allocatable, dimension(:,:),save :: cpl_snow, cpl_evap, cpl_tsol
1288  real, allocatable, dimension(:,:),save :: cpl_fder, cpl_albe, cpl_taux
1289  real, allocatable, dimension(:,:),save :: cpl_tauy
1290  REAL, ALLOCATABLE, DIMENSION(:,:),SAVE :: cpl_rriv, cpl_rcoa, cpl_rlic
1291!!$
1292! variables tampons avant le passage au coupleur
1293  real, allocatable, dimension(:,:,:),save :: tmp_sols, tmp_nsol, tmp_rain
1294  real, allocatable, dimension(:,:,:),save :: tmp_snow, tmp_evap, tmp_tsol
1295  real, allocatable, dimension(:,:,:),save :: tmp_fder, tmp_albe, tmp_taux
1296!!$  real, allocatable, dimension(:,:,:),save :: tmp_tauy, tmp_rriv, tmp_rcoa
1297  REAL, ALLOCATABLE, DIMENSION(:,:,:),SAVE :: tmp_tauy
1298! variables a passer au coupleur
1299  real, dimension(iim, jjm+1) :: wri_sol_ice, wri_sol_sea, wri_nsol_ice
1300  real, dimension(iim, jjm+1) :: wri_nsol_sea, wri_fder_ice, wri_evap_ice
1301  REAL, DIMENSION(iim, jjm+1) :: wri_evap_sea, wri_rcoa, wri_rriv
1302  REAL, DIMENSION(iim, jjm+1) :: wri_rain, wri_snow, wri_taux, wri_tauy
1303  REAL, DIMENSION(iim, jjm+1) :: wri_calv
1304  REAL, DIMENSION(iim, jjm+1) :: wri_tauxx, wri_tauyy, wri_tauzz
1305  REAL, DIMENSION(iim, jjm+1) :: tmp_lon, tmp_lat
1306! variables relues par le coupleur
1307! read_sic = fraction de glace
1308! read_sit = temperature de glace
1309  real, allocatable, dimension(:,:),save :: read_sst, read_sic, read_sit
1310  real, allocatable, dimension(:,:),save :: read_alb_sic
1311! variable tampon
1312  real, dimension(klon)       :: tamp_sic
1313! sauvegarde des fractions de surface d'un pas de temps a l'autre apres
1314! l'avoir lu
1315  real, allocatable,dimension(:,:),save :: pctsrf_sav
1316  real, dimension(iim, jjm+1, 2) :: tamp_srf
1317  integer, allocatable, dimension(:), save :: tamp_ind
1318  real, allocatable, dimension(:,:),save :: tamp_zmasq
1319  real, dimension(iim, jjm+1) :: deno
1320  integer                     :: idtime
1321  integer, allocatable,dimension(:),save :: unity
1322!
1323  logical, save    :: first_appel = .true.
1324  logical,save          :: print
1325!maf
1326! variables pour avoir une sortie IOIPSL des champs echanges
1327  CHARACTER*80,SAVE :: clintocplnam, clfromcplnam
1328  INTEGER, SAVE :: jf,nhoridct,nidct
1329  INTEGER, SAVE :: nhoridcs,nidcs
1330  INTEGER :: ndexct(iim*(jjm+1)),ndexcs(iim*(jjm+1))
1331  REAL :: zx_lon(iim,jjm+1), zx_lat(iim,jjm+1), zjulian
1332  integer :: idayref, itau_w
1333#include "param_cou.h"
1334#include "inc_cpl.h"
1335#include "temps.inc"
1336!
1337! Initialisation
1338!
1339  if (check) write(*,*)'Entree ',modname,'nisurf = ',nisurf
1340 
1341  if (first_appel) then
1342    error = 0
1343    allocate(unity(klon), stat = error)
1344    if ( error  /=0) then
1345      abort_message='Pb allocation variable unity'
1346      call abort_gcm(modname,abort_message,1)
1347    endif
1348    allocate(pctsrf_sav(klon,nbsrf), stat = error)
1349    if ( error  /=0) then
1350      abort_message='Pb allocation variable pctsrf_sav'
1351      call abort_gcm(modname,abort_message,1)
1352    endif
1353    pctsrf_sav = 0.
1354
1355    do ig = 1, klon
1356      unity(ig) = ig
1357    enddo
1358    sum_error = 0
1359    allocate(cpl_sols(klon,2), stat = error); sum_error = sum_error + error
1360    allocate(cpl_nsol(klon,2), stat = error); sum_error = sum_error + error
1361    allocate(cpl_rain(klon,2), stat = error); sum_error = sum_error + error
1362    allocate(cpl_snow(klon,2), stat = error); sum_error = sum_error + error
1363    allocate(cpl_evap(klon,2), stat = error); sum_error = sum_error + error
1364    allocate(cpl_tsol(klon,2), stat = error); sum_error = sum_error + error
1365    allocate(cpl_fder(klon,2), stat = error); sum_error = sum_error + error
1366    allocate(cpl_albe(klon,2), stat = error); sum_error = sum_error + error
1367    allocate(cpl_taux(klon,2), stat = error); sum_error = sum_error + error
1368    allocate(cpl_tauy(klon,2), stat = error); sum_error = sum_error + error
1369    ALLOCATE(cpl_rriv(iim,jjm+1), stat=error); sum_error = sum_error + error
1370    ALLOCATE(cpl_rcoa(iim,jjm+1), stat=error); sum_error = sum_error + error
1371    ALLOCATE(cpl_rlic(iim,jjm+1), stat=error); sum_error = sum_error + error
1372!!
1373    allocate(read_sst(iim, jjm+1), stat = error); sum_error = sum_error + error
1374    allocate(read_sic(iim, jjm+1), stat = error); sum_error = sum_error + error
1375    allocate(read_sit(iim, jjm+1), stat = error); sum_error = sum_error + error
1376    allocate(read_alb_sic(iim, jjm+1), stat = error); sum_error = sum_error + error
1377
1378    if (sum_error /= 0) then
1379      abort_message='Pb allocation variables couplees'
1380      call abort_gcm(modname,abort_message,1)
1381    endif
1382    cpl_sols = 0.; cpl_nsol = 0.; cpl_rain = 0.; cpl_snow = 0.
1383    cpl_evap = 0.; cpl_tsol = 0.; cpl_fder = 0.; cpl_albe = 0.
1384    cpl_taux = 0.; cpl_tauy = 0.; cpl_rriv = 0.; cpl_rcoa = 0.; cpl_rlic = 0.
1385
1386    sum_error = 0
1387    allocate(tamp_ind(klon), stat = error); sum_error = sum_error + error
1388    allocate(tamp_zmasq(iim, jjm+1), stat = error); sum_error = sum_error + error   
1389    do ig = 1, klon
1390      tamp_ind(ig) = ig
1391    enddo
1392    call gath2cpl(zmasq, tamp_zmasq, klon, klon, iim, jjm, tamp_ind)
1393!
1394! initialisation couplage
1395!
1396    idtime = int(dtime)
1397#ifdef CPP_PSMILE
1398    CALL inicma(iim, (jjm+1))
1399#else
1400   call inicma(npas , nexca, idtime,(jjm+1)*iim)
1401#endif
1402
1403!
1404! initialisation sorties netcdf
1405!
1406    idayref = day_ini
1407    CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
1408    CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlon,zx_lon)
1409    DO i = 1, iim
1410      zx_lon(i,1) = rlon(i+1)
1411      zx_lon(i,jjm+1) = rlon(i+1)
1412    ENDDO
1413    CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlat,zx_lat)
1414    clintocplnam="cpl_atm_tauflx"
1415    CALL histbeg(clintocplnam, iim,zx_lon(:,1),jjm+1,zx_lat(1,:),1,iim,1,jjm+1, &
1416       & itau_phy,zjulian,dtime,nhoridct,nidct)
1417! no vertical axis
1418    CALL histdef(nidct, 'tauxe','tauxe', &
1419         & "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1420    CALL histdef(nidct, 'tauyn','tauyn', &
1421         & "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1422    CALL histdef(nidct, 'tmp_lon','tmp_lon', &
1423         & "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1424    CALL histdef(nidct, 'tmp_lat','tmp_lat', &
1425         & "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1426    DO jf=1,jpflda2o1 + jpflda2o2
1427      CALL histdef(nidct, cl_writ(jf),cl_writ(jf), &
1428         & "-",iim, jjm+1, nhoridct, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1429    END DO
1430    CALL histend(nidct)
1431    CALL histsync(nidct)
1432
1433    clfromcplnam="cpl_atm_sst"
1434    CALL histbeg(clfromcplnam, iim,zx_lon(:,1),jjm+1,zx_lat(1,:),1,iim,1,jjm+1, &
1435       & 0,zjulian,dtime,nhoridcs,nidcs)
1436! no vertical axis
1437    DO jf=1,jpfldo2a
1438      CALL histdef(nidcs, cl_read(jf),cl_read(jf), &
1439         & "-",iim, jjm+1, nhoridcs, 1, 1, 1, -99, 32, "inst", dtime,dtime)
1440    END DO
1441    CALL histend(nidcs)
1442    CALL histsync(nidcs)
1443
1444! pour simuler la fonte des glaciers antarctiques
1445!
1446    surf_maille = (4. * rpi * ra**2) / (iim * (jjm +1))
1447    ALLOCATE(coeff_iceberg(iim,jjm+1), stat=error)
1448    if (error /= 0) then
1449      abort_message='Pb allocation variable coeff_iceberg'
1450      call abort_gcm(modname,abort_message,1)
1451    endif
1452    open (12,file='flux_iceberg',form='formatted',status='old')
1453    read (12,*) coeff_iceberg
1454    close (12)
1455    num_antarctic = max(1, count(coeff_iceberg > 0))
1456   
1457    first_appel = .false.
1458  endif ! fin if (first_appel)
1459
1460! Initialisations
1461
1462! calcul des fluxs a passer
1463
1464  cpl_index = 1
1465  if (nisurf == is_sic) cpl_index = 2
1466  if (cumul) then
1467    if (check) write(*,*) modname, 'cumul des champs'
1468    do ig = 1, knon
1469      cpl_sols(ig,cpl_index) = cpl_sols(ig,cpl_index) &
1470       &                          + swdown(ig)      / FLOAT(nexca)
1471      cpl_nsol(ig,cpl_index) = cpl_nsol(ig,cpl_index) &
1472       &                          + (lwdown(ig) + fluxlat(ig) +fluxsens(ig))&
1473       &                                / FLOAT(nexca)
1474      cpl_rain(ig,cpl_index) = cpl_rain(ig,cpl_index) &
1475       &                          + precip_rain(ig) / FLOAT(nexca)
1476      cpl_snow(ig,cpl_index) = cpl_snow(ig,cpl_index) &
1477       &                          + precip_snow(ig) / FLOAT(nexca)
1478      cpl_evap(ig,cpl_index) = cpl_evap(ig,cpl_index) &
1479       &                          + evap(ig)        / FLOAT(nexca)
1480      cpl_tsol(ig,cpl_index) = cpl_tsol(ig,cpl_index) &
1481       &                          + tsurf(ig)       / FLOAT(nexca)
1482      cpl_fder(ig,cpl_index) = cpl_fder(ig,cpl_index) &
1483       &                          + fder(ig)        / FLOAT(nexca)
1484      cpl_albe(ig,cpl_index) = cpl_albe(ig,cpl_index) &
1485       &                          + albsol(ig)      / FLOAT(nexca)
1486      cpl_taux(ig,cpl_index) = cpl_taux(ig,cpl_index) &
1487       &                          + taux(ig)        / FLOAT(nexca)
1488      cpl_tauy(ig,cpl_index) = cpl_tauy(ig,cpl_index) &
1489       &                          + tauy(ig)        / FLOAT(nexca)
1490    enddo
1491    IF (cpl_index .EQ. 1) THEN
1492        cpl_rriv(:,:) = cpl_rriv(:,:) + tmp_rriv(:,:) / FLOAT(nexca)
1493        cpl_rcoa(:,:) = cpl_rcoa(:,:) + tmp_rcoa(:,:) / FLOAT(nexca)
1494        cpl_rlic(:,:) = cpl_rlic(:,:) + tmp_rlic(:,:) / FLOAT(nexca)
1495    ENDIF
1496  endif
1497
1498  if (mod(itime, nexca) == 1) then
1499!
1500! Demande des champs au coupleur
1501!
1502! Si le domaine considere est l'ocean, on lit les champs venant du coupleur
1503!
1504    if (nisurf == is_oce .and. .not. cumul) then
1505      if (check) write(*,*)'rentree fromcpl, itime-1 = ',itime-1
1506#ifdef CPP_PSMILE
1507      il_time_secs=(itime-1)*dtime
1508      CALL fromcpl(il_time_secs, iim, (jjm+1),                           &
1509     &        read_sst, read_sic, read_sit, read_alb_sic)
1510#else
1511      call fromcpl(itime-1,(jjm+1)*iim,                                  &
1512     &        read_sst, read_sic, read_sit, read_alb_sic)
1513#endif
1514!
1515! sorties NETCDF des champs recus
1516!
1517      ndexcs(:)=0
1518      itau_w = itau_phy + itime
1519      CALL histwrite(nidcs,cl_read(1),itau_w,read_sst,iim*(jjm+1),ndexcs)
1520      CALL histwrite(nidcs,cl_read(2),itau_w,read_sic,iim*(jjm+1),ndexcs)
1521      CALL histwrite(nidcs,cl_read(3),itau_w,read_alb_sic,iim*(jjm+1),ndexcs)
1522      CALL histwrite(nidcs,cl_read(4),itau_w,read_sit,iim*(jjm+1),ndexcs)
1523      CALL histsync(nidcs)
1524! pas utile      IF (npas-itime.LT.nexca )CALL histclo(nidcs)
1525
1526      do j = 1, jjm + 1
1527        do ig = 1, iim
1528          if (abs(1. - read_sic(ig,j)) < 0.00001) then
1529            read_sst(ig,j) = RTT - 1.8
1530            read_sit(ig,j) = read_sit(ig,j) / read_sic(ig,j)
1531            read_alb_sic(ig,j) = read_alb_sic(ig,j) / read_sic(ig,j)
1532          else if (abs(read_sic(ig,j)) < 0.00001) then
1533            read_sst(ig,j) = read_sst(ig,j) / (1. - read_sic(ig,j))
1534            read_sit(ig,j) = read_sst(ig,j)
1535            read_alb_sic(ig,j) =  0.6
1536          else
1537            read_sst(ig,j) = read_sst(ig,j) / (1. - read_sic(ig,j))
1538            read_sit(ig,j) = read_sit(ig,j) / read_sic(ig,j)
1539            read_alb_sic(ig,j) = read_alb_sic(ig,j) / read_sic(ig,j)
1540          endif
1541        enddo
1542      enddo
1543!
1544! transformer read_sic en pctsrf_sav
1545!
1546      call cpl2gath(read_sic, tamp_sic , klon, klon,iim,jjm, unity)
1547      do ig = 1, klon
1548        IF (pctsrf(ig,is_oce) > epsfra .OR.            &
1549     &             pctsrf(ig,is_sic) > epsfra) THEN
1550          pctsrf_sav(ig,is_sic) = (pctsrf(ig,is_oce) + pctsrf(ig,is_sic)) &
1551     &                               * tamp_sic(ig)
1552          pctsrf_sav(ig,is_oce) = (pctsrf(ig,is_oce) + pctsrf(ig,is_sic)) &
1553     &                        - pctsrf_sav(ig,is_sic)
1554        endif
1555      enddo
1556!
1557! Pour rattraper des erreurs d'arrondis
1558!
1559      where (abs(pctsrf_sav(:,is_sic)) .le. 2.*epsilon(pctsrf_sav(1,is_sic)))
1560        pctsrf_sav(:,is_sic) = 0.
1561        pctsrf_sav(:,is_oce) = pctsrf(:,is_oce) + pctsrf(:,is_sic)
1562      endwhere
1563      where (abs(pctsrf_sav(:,is_oce)) .le. 2.*epsilon(pctsrf_sav(1,is_oce)))
1564        pctsrf_sav(:,is_sic) = pctsrf(:,is_oce) + pctsrf(:,is_sic)
1565        pctsrf_sav(:,is_oce) = 0.
1566      endwhere
1567      if (minval(pctsrf_sav(:,is_oce)) < 0.) then
1568        write(*,*)'Pb fraction ocean inferieure a 0'
1569        write(*,*)'au point ',minloc(pctsrf_sav(:,is_oce))
1570        write(*,*)'valeur = ',minval(pctsrf_sav(:,is_oce))
1571        abort_message = 'voir ci-dessus'
1572        call abort_gcm(modname,abort_message,1)
1573      endif
1574      if (minval(pctsrf_sav(:,is_sic)) < 0.) then
1575        write(*,*)'Pb fraction glace inferieure a 0'
1576        write(*,*)'au point ',minloc(pctsrf_sav(:,is_sic))
1577        write(*,*)'valeur = ',minval(pctsrf_sav(:,is_sic))
1578        abort_message = 'voir ci-dessus'
1579        call abort_gcm(modname,abort_message,1)
1580      endif
1581    endif
1582  endif                         ! fin mod(itime, nexca) == 1
1583
1584  if (mod(itime, nexca) == 0) then
1585!
1586! allocation memoire
1587    if (nisurf == is_oce .and. (.not. cumul) ) then
1588      sum_error = 0
1589      allocate(tmp_sols(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1590      allocate(tmp_nsol(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1591      allocate(tmp_rain(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1592      allocate(tmp_snow(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1593      allocate(tmp_evap(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1594      allocate(tmp_tsol(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1595      allocate(tmp_fder(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1596      allocate(tmp_albe(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1597      allocate(tmp_taux(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1598      allocate(tmp_tauy(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1599!!$      allocate(tmp_rriv(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1600!!$      allocate(tmp_rcoa(iim,jjm+1,2), stat=error); sum_error = sum_error + error
1601      if (sum_error /= 0) then
1602        abort_message='Pb allocation variables couplees pour l''ecriture'
1603        call abort_gcm(modname,abort_message,1)
1604      endif
1605    endif
1606
1607!
1608! Mise sur la bonne grille des champs a passer au coupleur
1609!
1610    cpl_index = 1
1611    if (nisurf == is_sic) cpl_index = 2
1612    call gath2cpl(cpl_sols(1,cpl_index), tmp_sols(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1613    call gath2cpl(cpl_nsol(1,cpl_index), tmp_nsol(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1614    call gath2cpl(cpl_rain(1,cpl_index), tmp_rain(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1615    call gath2cpl(cpl_snow(1,cpl_index), tmp_snow(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1616    call gath2cpl(cpl_evap(1,cpl_index), tmp_evap(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1617    call gath2cpl(cpl_tsol(1,cpl_index), tmp_tsol(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1618    call gath2cpl(cpl_fder(1,cpl_index), tmp_fder(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1619    call gath2cpl(cpl_albe(1,cpl_index), tmp_albe(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1620    call gath2cpl(cpl_taux(1,cpl_index), tmp_taux(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1621    call gath2cpl(cpl_tauy(1,cpl_index), tmp_tauy(1,1,cpl_index), klon, knon,iim,jjm,                  knindex)
1622
1623!
1624! Si le domaine considere est la banquise, on envoie les champs au coupleur
1625!
1626    if (nisurf == is_sic .and. cumul) then
1627      wri_rain = 0.; wri_snow = 0.; wri_rcoa = 0.; wri_rriv = 0.
1628      wri_taux = 0.; wri_tauy = 0.
1629      call gath2cpl(pctsrf(1,is_oce), tamp_srf(1,1,1), klon, klon, iim, jjm, tamp_ind)
1630      call gath2cpl(pctsrf(1,is_sic), tamp_srf(1,1,2), klon, klon, iim, jjm, tamp_ind)
1631
1632      wri_sol_ice = tmp_sols(:,:,2)
1633      wri_sol_sea = tmp_sols(:,:,1)
1634      wri_nsol_ice = tmp_nsol(:,:,2)
1635      wri_nsol_sea = tmp_nsol(:,:,1)
1636      wri_fder_ice = tmp_fder(:,:,2)
1637      wri_evap_ice = tmp_evap(:,:,2)
1638      wri_evap_sea = tmp_evap(:,:,1)
1639!!$PB
1640      wri_rriv = cpl_rriv(:,:)
1641      wri_rcoa = cpl_rcoa(:,:)
1642      DO j = 1, jjm + 1
1643        wri_calv(:,j) = sum(cpl_rlic(:,j)) / iim
1644      enddo
1645
1646      where (tamp_zmasq /= 1.)
1647        deno =  tamp_srf(:,:,1) + tamp_srf(:,:,2)
1648        wri_rain = tmp_rain(:,:,1) * tamp_srf(:,:,1) / deno +    &
1649      &            tmp_rain(:,:,2) * tamp_srf(:,:,2) / deno
1650        wri_snow = tmp_snow(:,:,1) * tamp_srf(:,:,1) / deno +    &
1651      &            tmp_snow(:,:,2) * tamp_srf(:,:,2) / deno
1652        wri_taux = tmp_taux(:,:,1) * tamp_srf(:,:,1) / deno +    &
1653      &            tmp_taux(:,:,2) * tamp_srf(:,:,2) / deno
1654        wri_tauy = tmp_tauy(:,:,1) * tamp_srf(:,:,1) / deno +    &
1655      &            tmp_tauy(:,:,2) * tamp_srf(:,:,2) / deno
1656      endwhere
1657!
1658! pour simuler la fonte des glaciers antarctiques
1659!
1660!$$$        wri_rain = wri_rain      &
1661!$$$      &     + coeff_iceberg * cte_flux_iceberg / (num_antarctic * surf_maille)
1662!      wri_calv = coeff_iceberg * cte_flux_iceberg / (num_antarctic * surf_maille)
1663!
1664! on passe les coordonnées de la grille
1665!
1666
1667      CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlon,tmp_lon)
1668      CALL gr_fi_ecrit(1,klon,iim,jjm+1,rlat,tmp_lat)
1669
1670      DO i = 1, iim
1671        tmp_lon(i,1) = rlon(i+1)
1672        tmp_lon(i,jjm + 1) = rlon(i+1)
1673      ENDDO
1674!
1675! sortie netcdf des champs pour le changement de repere
1676!
1677      ndexct(:)=0
1678      CALL histwrite(nidct,'tauxe',itau_w,wri_taux,iim*(jjm+1),ndexct)
1679      CALL histwrite(nidct,'tauyn',itau_w,wri_tauy,iim*(jjm+1),ndexct)
1680      CALL histwrite(nidct,'tmp_lon',itau_w,tmp_lon,iim*(jjm+1),ndexct)
1681      CALL histwrite(nidct,'tmp_lat',itau_w,tmp_lat,iim*(jjm+1),ndexct)
1682
1683!
1684! calcul 3 coordonnées du vent
1685!
1686      CALL atm2geo (iim , jjm + 1, wri_taux, wri_tauy, tmp_lon, tmp_lat, &
1687         & wri_tauxx, wri_tauyy, wri_tauzz )
1688!
1689! sortie netcdf des champs apres changement de repere et juste avant
1690! envoi au coupleur
1691!
1692      CALL histwrite(nidct,cl_writ(1),itau_w,wri_sol_ice,iim*(jjm+1),ndexct)
1693      CALL histwrite(nidct,cl_writ(2),itau_w,wri_sol_sea,iim*(jjm+1),ndexct)
1694      CALL histwrite(nidct,cl_writ(3),itau_w,wri_nsol_ice,iim*(jjm+1),ndexct)
1695      CALL histwrite(nidct,cl_writ(4),itau_w,wri_nsol_sea,iim*(jjm+1),ndexct)
1696      CALL histwrite(nidct,cl_writ(5),itau_w,wri_fder_ice,iim*(jjm+1),ndexct)
1697      CALL histwrite(nidct,cl_writ(6),itau_w,wri_evap_ice,iim*(jjm+1),ndexct)
1698      CALL histwrite(nidct,cl_writ(7),itau_w,wri_evap_sea,iim*(jjm+1),ndexct)
1699      CALL histwrite(nidct,cl_writ(8),itau_w,wri_rain,iim*(jjm+1),ndexct)
1700      CALL histwrite(nidct,cl_writ(9),itau_w,wri_snow,iim*(jjm+1),ndexct)
1701      CALL histwrite(nidct,cl_writ(10),itau_w,wri_rcoa,iim*(jjm+1),ndexct)
1702      CALL histwrite(nidct,cl_writ(11),itau_w,wri_rriv,iim*(jjm+1),ndexct)
1703      CALL histwrite(nidct,cl_writ(12),itau_w,wri_calv,iim*(jjm+1),ndexct)
1704      CALL histwrite(nidct,cl_writ(13),itau_w,wri_tauxx,iim*(jjm+1),ndexct)
1705      CALL histwrite(nidct,cl_writ(14),itau_w,wri_tauyy,iim*(jjm+1),ndexct)
1706      CALL histwrite(nidct,cl_writ(15),itau_w,wri_tauzz,iim*(jjm+1),ndexct)
1707      CALL histwrite(nidct,cl_writ(16),itau_w,wri_tauxx,iim*(jjm+1),ndexct)
1708      CALL histwrite(nidct,cl_writ(17),itau_w,wri_tauyy,iim*(jjm+1),ndexct)
1709      CALL histwrite(nidct,cl_writ(18),itau_w,wri_tauzz,iim*(jjm+1),ndexct)
1710      CALL histsync(nidct)
1711! pas utile      IF (lafin) CALL histclo(nidct)
1712#ifdef CPP_PSMILE
1713      il_time_secs=(itime-1)*dtime
1714
1715      CALL intocpl(il_time_secs, iim, jjm+1, wri_sol_ice, wri_sol_sea, wri_nsol_ice,&
1716      & wri_nsol_sea, wri_fder_ice, wri_evap_ice, wri_evap_sea, wri_rain, &
1717      & wri_snow, wri_rcoa, wri_rriv, wri_calv, wri_tauxx, wri_tauyy,     &
1718      & wri_tauzz, wri_tauxx, wri_tauyy, wri_tauzz,lafin )
1719#else
1720      call intocpl(itime, (jjm+1)*iim, wri_sol_ice, wri_sol_sea, wri_nsol_ice,&
1721      & wri_nsol_sea, wri_fder_ice, wri_evap_ice, wri_evap_sea, wri_rain, &
1722      & wri_snow, wri_rcoa, wri_rriv, wri_calv, wri_tauxx, wri_tauyy,     &
1723      & wri_tauzz, wri_tauxx, wri_tauyy, wri_tauzz,lafin )
1724#endif
1725!
1726      cpl_sols = 0.; cpl_nsol = 0.; cpl_rain = 0.; cpl_snow = 0.
1727      cpl_evap = 0.; cpl_tsol = 0.; cpl_fder = 0.; cpl_albe = 0.
1728      cpl_taux = 0.; cpl_tauy = 0.; cpl_rriv = 0.; cpl_rcoa = 0.; cpl_rlic = 0.
1729!
1730! deallocation memoire variables temporaires
1731!
1732      sum_error = 0
1733      deallocate(tmp_sols, stat=error); sum_error = sum_error + error
1734      deallocate(tmp_nsol, stat=error); sum_error = sum_error + error
1735      deallocate(tmp_rain, stat=error); sum_error = sum_error + error
1736      deallocate(tmp_snow, stat=error); sum_error = sum_error + error
1737      deallocate(tmp_evap, stat=error); sum_error = sum_error + error
1738      deallocate(tmp_fder, stat=error); sum_error = sum_error + error
1739      deallocate(tmp_tsol, stat=error); sum_error = sum_error + error
1740      deallocate(tmp_albe, stat=error); sum_error = sum_error + error
1741      deallocate(tmp_taux, stat=error); sum_error = sum_error + error
1742      deallocate(tmp_tauy, stat=error); sum_error = sum_error + error
1743!!$PB
1744!!$      deallocate(tmp_rriv, stat=error); sum_error = sum_error + error
1745!!$      deallocate(tmp_rcoa, stat=error); sum_error = sum_error + error
1746      if (sum_error /= 0) then
1747        abort_message='Pb deallocation variables couplees'
1748        call abort_gcm(modname,abort_message,1)
1749      endif
1750
1751    endif
1752
1753  endif            ! fin (mod(itime, nexca) == 0)
1754!
1755! on range les variables lues/sauvegardees dans les bonnes variables de sortie
1756!
1757  if (nisurf == is_oce) then
1758    call cpl2gath(read_sst, tsurf_new, klon, knon,iim,jjm, knindex)
1759  else if (nisurf == is_sic) then
1760    call cpl2gath(read_sit, tsurf_new, klon, knon,iim,jjm, knindex)
1761    call cpl2gath(read_alb_sic, alb_new, klon, knon,iim,jjm, knindex)
1762  endif
1763  pctsrf_new(:,nisurf) = pctsrf_sav(:,nisurf)
1764 
1765!  if (lafin) call quitcpl
1766
1767  END SUBROUTINE interfoce_cpl
1768!
1769!#########################################################################
1770!
1771
1772  SUBROUTINE interfoce_slab(nisurf)
1773
1774! Cette routine sert d'interface entre le modele atmospherique et un
1775! modele de 'slab' ocean
1776!
1777! L. Fairhead 02/2000
1778!
1779! input:
1780!   nisurf       index de la surface a traiter (1 = sol continental)
1781!
1782! output:
1783!
1784
1785! Parametres d'entree
1786  integer, intent(IN) :: nisurf
1787
1788  END SUBROUTINE interfoce_slab
1789!
1790!#########################################################################
1791!
1792  SUBROUTINE interfoce_lim(itime, dtime, jour, &
1793     & klon, nisurf, knon, knindex, &
1794     & debut,  &
1795     & lmt_sst, pctsrf_new)
1796
1797! Cette routine sert d'interface entre le modele atmospherique et un fichier
1798! de conditions aux limites
1799!
1800! L. Fairhead 02/2000
1801!
1802! input:
1803!   itime        numero du pas de temps courant
1804!   dtime        pas de temps de la physique (en s)
1805!   jour         jour a lire dans l'annee
1806!   nisurf       index de la surface a traiter (1 = sol continental)
1807!   knon         nombre de points dans le domaine a traiter
1808!   knindex      index des points de la surface a traiter
1809!   klon         taille de la grille
1810!   debut        logical: 1er appel a la physique (initialisation)
1811!
1812! output:
1813!   lmt_sst      SST lues dans le fichier de CL
1814!   pctsrf_new   sous-maille fractionnelle
1815!
1816
1817
1818! Parametres d'entree
1819  integer, intent(IN) :: itime
1820  real   , intent(IN) :: dtime
1821  integer, intent(IN) :: jour
1822  integer, intent(IN) :: nisurf
1823  integer, intent(IN) :: knon
1824  integer, intent(IN) :: klon
1825  integer, dimension(klon), intent(in) :: knindex
1826  logical, intent(IN) :: debut
1827
1828! Parametres de sortie
1829  real, intent(out), dimension(klon) :: lmt_sst
1830  real, intent(out), dimension(klon,nbsrf) :: pctsrf_new
1831
1832! Variables locales
1833  integer     :: ii
1834  INTEGER,save :: lmt_pas     ! frequence de lecture des conditions limites
1835                             ! (en pas de physique)
1836  logical,save :: deja_lu    ! pour indiquer que le jour a lire a deja
1837                             ! lu pour une surface precedente
1838  integer,save :: jour_lu
1839  integer      :: ierr
1840  character (len = 20) :: modname = 'interfoce_lim'
1841  character (len = 80) :: abort_message
1842  character (len = 20),save :: fich ='limit.nc'
1843  logical, save     :: newlmt = .TRUE.
1844  logical, save     :: check = .FALSE.
1845! Champs lus dans le fichier de CL
1846  real, allocatable , save, dimension(:) :: sst_lu, rug_lu, nat_lu
1847  real, allocatable , save, dimension(:,:) :: pct_tmp
1848!
1849! quelques variables pour netcdf
1850!
1851#include "netcdf.inc"
1852  integer              :: nid, nvarid
1853  integer, dimension(2) :: start, epais
1854!
1855! Fin déclaration
1856!
1857   
1858  if (debut .and. .not. allocated(sst_lu)) then
1859    lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
1860    jour_lu = jour - 1
1861    allocate(sst_lu(klon))
1862    allocate(nat_lu(klon))
1863    allocate(pct_tmp(klon,nbsrf))
1864  endif
1865
1866  if ((jour - jour_lu) /= 0) deja_lu = .false.
1867 
1868  if (check) write(*,*)modname,' :: jour, jour_lu, deja_lu', jour, jour_lu, deja_lu
1869  if (check) write(*,*)modname,' :: itime, lmt_pas ', itime, lmt_pas,dtime
1870
1871! Tester d'abord si c'est le moment de lire le fichier
1872  if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu) then
1873!
1874! Ouverture du fichier
1875!
1876    fich = trim(fich)
1877    ierr = NF_OPEN (fich, NF_NOWRITE,nid)
1878    if (ierr.NE.NF_NOERR) then
1879      abort_message = 'Pb d''ouverture du fichier de conditions aux limites'
1880      call abort_gcm(modname,abort_message,1)
1881    endif
1882!
1883! La tranche de donnees a lire:
1884!
1885    start(1) = 1
1886    start(2) = jour
1887    epais(1) = klon
1888    epais(2) = 1
1889!
1890    if (newlmt) then
1891!
1892! Fraction "ocean"
1893!
1894      ierr = NF_INQ_VARID(nid, 'FOCE', nvarid)
1895      if (ierr /= NF_NOERR) then
1896        abort_message = 'Le champ <FOCE> est absent'
1897        call abort_gcm(modname,abort_message,1)
1898      endif
1899#ifdef NC_DOUBLE
1900      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_oce))
1901#else
1902      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_oce))
1903#endif
1904      if (ierr /= NF_NOERR) then
1905        abort_message = 'Lecture echouee pour <FOCE>'
1906        call abort_gcm(modname,abort_message,1)
1907      endif
1908!
1909! Fraction "glace de mer"
1910!
1911      ierr = NF_INQ_VARID(nid, 'FSIC', nvarid)
1912      if (ierr /= NF_NOERR) then
1913        abort_message = 'Le champ <FSIC> est absent'
1914        call abort_gcm(modname,abort_message,1)
1915      endif
1916#ifdef NC_DOUBLE
1917      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_sic))
1918#else
1919      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_sic))
1920#endif
1921      if (ierr /= NF_NOERR) then
1922        abort_message = 'Lecture echouee pour <FSIC>'
1923        call abort_gcm(modname,abort_message,1)
1924      endif
1925!
1926! Fraction "terre"
1927!
1928      ierr = NF_INQ_VARID(nid, 'FTER', nvarid)
1929      if (ierr /= NF_NOERR) then
1930        abort_message = 'Le champ <FTER> est absent'
1931        call abort_gcm(modname,abort_message,1)
1932      endif
1933#ifdef NC_DOUBLE
1934      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_ter))
1935#else
1936      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_ter))
1937#endif
1938      if (ierr /= NF_NOERR) then
1939        abort_message = 'Lecture echouee pour <FTER>'
1940        call abort_gcm(modname,abort_message,1)
1941      endif
1942!
1943! Fraction "glacier terre"
1944!
1945      ierr = NF_INQ_VARID(nid, 'FLIC', nvarid)
1946      if (ierr /= NF_NOERR) then
1947        abort_message = 'Le champ <FLIC> est absent'
1948        call abort_gcm(modname,abort_message,1)
1949      endif
1950#ifdef NC_DOUBLE
1951      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_lic))
1952#else
1953      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_lic))
1954#endif
1955      if (ierr /= NF_NOERR) then
1956        abort_message = 'Lecture echouee pour <FLIC>'
1957        call abort_gcm(modname,abort_message,1)
1958      endif
1959!
1960    else  ! on en est toujours a rnatur
1961!
1962      ierr = NF_INQ_VARID(nid, 'NAT', nvarid)
1963      if (ierr /= NF_NOERR) then
1964        abort_message = 'Le champ <NAT> est absent'
1965        call abort_gcm(modname,abort_message,1)
1966      endif
1967#ifdef NC_DOUBLE
1968      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, nat_lu)
1969#else
1970      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, nat_lu)
1971#endif
1972      if (ierr /= NF_NOERR) then
1973        abort_message = 'Lecture echouee pour <NAT>'
1974        call abort_gcm(modname,abort_message,1)
1975      endif
1976!
1977! Remplissage des fractions de surface
1978! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice
1979!
1980      pct_tmp = 0.0
1981      do ii = 1, klon
1982        pct_tmp(ii,nint(nat_lu(ii)) + 1) = 1.
1983      enddo
1984
1985!
1986!  On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire
1987!
1988      pctsrf_new = pct_tmp
1989      pctsrf_new (:,2)= pct_tmp (:,1)
1990      pctsrf_new (:,1)= pct_tmp (:,2)
1991      pct_tmp = pctsrf_new
1992    endif ! fin test sur newlmt
1993!
1994! Lecture SST
1995!
1996    ierr = NF_INQ_VARID(nid, 'SST', nvarid)
1997    if (ierr /= NF_NOERR) then
1998      abort_message = 'Le champ <SST> est absent'
1999      call abort_gcm(modname,abort_message,1)
2000    endif
2001#ifdef NC_DOUBLE
2002    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, sst_lu)
2003#else
2004    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, sst_lu)
2005#endif
2006    if (ierr /= NF_NOERR) then
2007      abort_message = 'Lecture echouee pour <SST>'
2008      call abort_gcm(modname,abort_message,1)
2009    endif   
2010
2011!
2012! Fin de lecture
2013!
2014    ierr = NF_CLOSE(nid)
2015    deja_lu = .true.
2016    jour_lu = jour
2017  endif
2018!
2019! Recopie des variables dans les champs de sortie
2020!
2021  lmt_sst = 999999999.
2022  do ii = 1, knon
2023    lmt_sst(ii) = sst_lu(knindex(ii))
2024  enddo
2025
2026  pctsrf_new(:,is_oce) = pct_tmp(:,is_oce)
2027  pctsrf_new(:,is_sic) = pct_tmp(:,is_sic)
2028
2029  END SUBROUTINE interfoce_lim
2030
2031!
2032!#########################################################################
2033!
2034  SUBROUTINE interfsur_lim(itime, dtime, jour, &
2035     & klon, nisurf, knon, knindex, &
2036     & debut,  &
2037     & lmt_alb, lmt_rug)
2038
2039! Cette routine sert d'interface entre le modele atmospherique et un fichier
2040! de conditions aux limites
2041!
2042! L. Fairhead 02/2000
2043!
2044! input:
2045!   itime        numero du pas de temps courant
2046!   dtime        pas de temps de la physique (en s)
2047!   jour         jour a lire dans l'annee
2048!   nisurf       index de la surface a traiter (1 = sol continental)
2049!   knon         nombre de points dans le domaine a traiter
2050!   knindex      index des points de la surface a traiter
2051!   klon         taille de la grille
2052!   debut        logical: 1er appel a la physique (initialisation)
2053!
2054! output:
2055!   lmt_sst      SST lues dans le fichier de CL
2056!   lmt_alb      Albedo lu
2057!   lmt_rug      longueur de rugosité lue
2058!   pctsrf_new   sous-maille fractionnelle
2059!
2060
2061
2062! Parametres d'entree
2063  integer, intent(IN) :: itime
2064  real   , intent(IN) :: dtime
2065  integer, intent(IN) :: jour
2066  integer, intent(IN) :: nisurf
2067  integer, intent(IN) :: knon
2068  integer, intent(IN) :: klon
2069  integer, dimension(klon), intent(in) :: knindex
2070  logical, intent(IN) :: debut
2071
2072! Parametres de sortie
2073  real, intent(out), dimension(klon) :: lmt_alb
2074  real, intent(out), dimension(klon) :: lmt_rug
2075
2076! Variables locales
2077  integer     :: ii
2078  integer,save :: lmt_pas     ! frequence de lecture des conditions limites
2079                             ! (en pas de physique)
2080  logical,save :: deja_lu_sur! pour indiquer que le jour a lire a deja
2081                             ! lu pour une surface precedente
2082  integer,save :: jour_lu_sur
2083  integer      :: ierr
2084  character (len = 20) :: modname = 'interfsur_lim'
2085  character (len = 80) :: abort_message
2086  character (len = 20),save :: fich ='limit.nc'
2087  logical,save     :: newlmt = .false.
2088  logical,save     :: check = .false.
2089! Champs lus dans le fichier de CL
2090  real, allocatable , save, dimension(:) :: alb_lu, rug_lu
2091!
2092! quelques variables pour netcdf
2093!
2094#include "netcdf.inc"
2095  integer ,save             :: nid, nvarid
2096  integer, dimension(2),save :: start, epais
2097!
2098! Fin déclaration
2099!
2100   
2101  if (debut) then
2102    lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
2103    jour_lu_sur = jour - 1
2104    allocate(alb_lu(klon))
2105    allocate(rug_lu(klon))
2106  endif
2107
2108  if ((jour - jour_lu_sur) /= 0) deja_lu_sur = .false.
2109 
2110  if (check) write(*,*)modname,':: jour_lu_sur, deja_lu_sur', jour_lu_sur, deja_lu_sur
2111  if (check) write(*,*)modname,':: itime, lmt_pas', itime, lmt_pas
2112  if (check) call flush(6)
2113
2114! Tester d'abord si c'est le moment de lire le fichier
2115  if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu_sur) then
2116!
2117! Ouverture du fichier
2118!
2119    fich = trim(fich)
2120    IF (check) WRITE(*,*)modname,' ouverture fichier ',fich
2121    if (check) CALL flush(6)
2122    ierr = NF_OPEN (fich, NF_NOWRITE,nid)
2123    if (ierr.NE.NF_NOERR) then
2124      abort_message = 'Pb d''ouverture du fichier de conditions aux limites'
2125      call abort_gcm(modname,abort_message,1)
2126    endif
2127!
2128! La tranche de donnees a lire:
2129 
2130    start(1) = 1
2131    start(2) = jour
2132    epais(1) = klon
2133    epais(2) = 1
2134!
2135! Lecture Albedo
2136!
2137    ierr = NF_INQ_VARID(nid, 'ALB', nvarid)
2138    if (ierr /= NF_NOERR) then
2139      abort_message = 'Le champ <ALB> est absent'
2140      call abort_gcm(modname,abort_message,1)
2141    endif
2142#ifdef NC_DOUBLE
2143    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, alb_lu)
2144#else
2145    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, alb_lu)
2146#endif
2147    if (ierr /= NF_NOERR) then
2148      abort_message = 'Lecture echouee pour <ALB>'
2149      call abort_gcm(modname,abort_message,1)
2150    endif
2151!
2152! Lecture rugosité
2153!
2154    ierr = NF_INQ_VARID(nid, 'RUG', nvarid)
2155    if (ierr /= NF_NOERR) then
2156      abort_message = 'Le champ <RUG> est absent'
2157      call abort_gcm(modname,abort_message,1)
2158    endif
2159#ifdef NC_DOUBLE
2160    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, rug_lu)
2161#else
2162    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, rug_lu)
2163#endif
2164    if (ierr /= NF_NOERR) then
2165      abort_message = 'Lecture echouee pour <RUG>'
2166      call abort_gcm(modname,abort_message,1)
2167    endif
2168
2169!
2170! Fin de lecture
2171!
2172    ierr = NF_CLOSE(nid)
2173    deja_lu_sur = .true.
2174    jour_lu_sur = jour
2175  endif
2176!
2177! Recopie des variables dans les champs de sortie
2178!
2179!!$  lmt_alb(:) = 0.0
2180!!$  lmt_rug(:) = 0.0
2181  lmt_alb(:) = 999999.
2182  lmt_rug(:) = 999999.
2183  DO ii = 1, knon
2184    lmt_alb(ii) = alb_lu(knindex(ii))
2185    lmt_rug(ii) = rug_lu(knindex(ii))
2186  enddo
2187
2188  END SUBROUTINE interfsur_lim
2189
2190!
2191!#########################################################################
2192!
2193
2194  SUBROUTINE calcul_fluxs( klon, knon, nisurf, dtime, &
2195     & tsurf, p1lay, cal, beta, coef1lay, ps, &
2196     & precip_rain, precip_snow, snow, qsurf, &
2197     & radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, &
2198     & petAcoef, peqAcoef, petBcoef, peqBcoef, &
2199     & tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
2200
2201! Cette routine calcule les fluxs en h et q a l'interface et eventuellement
2202! une temperature de surface (au cas ou ok_veget = false)
2203!
2204! L. Fairhead 4/2000
2205!
2206! input:
2207!   knon         nombre de points a traiter
2208!   nisurf       surface a traiter
2209!   tsurf        temperature de surface
2210!   p1lay        pression 1er niveau (milieu de couche)
2211!   cal          capacite calorifique du sol
2212!   beta         evap reelle
2213!   coef1lay     coefficient d'echange
2214!   ps           pression au sol
2215!   precip_rain  precipitations liquides
2216!   precip_snow  precipitations solides
2217!   snow         champs hauteur de neige
2218!   runoff       runoff en cas de trop plein
2219!   petAcoef     coeff. A de la resolution de la CL pour t
2220!   peqAcoef     coeff. A de la resolution de la CL pour q
2221!   petBcoef     coeff. B de la resolution de la CL pour t
2222!   peqBcoef     coeff. B de la resolution de la CL pour q
2223!   radsol       rayonnement net aus sol (LW + SW)
2224!   dif_grnd     coeff. diffusion vers le sol profond
2225!
2226! output:
2227!   tsurf_new    temperature au sol
2228!   qsurf        humidite de l'air au dessus du sol
2229!   fluxsens     flux de chaleur sensible
2230!   fluxlat      flux de chaleur latente
2231!   dflux_s      derivee du flux de chaleur sensible / Ts
2232!   dflux_l      derivee du flux de chaleur latente  / Ts
2233!
2234
2235#include "YOETHF.inc"
2236#include "FCTTRE.inc"
2237#include "indicesol.inc"
2238
2239! Parametres d'entree
2240  integer, intent(IN) :: knon, nisurf, klon
2241  real   , intent(IN) :: dtime
2242  real, dimension(klon), intent(IN) :: petAcoef, peqAcoef
2243  real, dimension(klon), intent(IN) :: petBcoef, peqBcoef
2244  real, dimension(klon), intent(IN) :: ps, q1lay
2245  real, dimension(klon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
2246  real, dimension(klon), intent(IN) :: precip_rain, precip_snow
2247  real, dimension(klon), intent(IN) :: radsol, dif_grnd
2248  real, dimension(klon), intent(IN) :: t1lay, u1lay, v1lay
2249  real, dimension(klon), intent(INOUT) :: snow, qsurf
2250
2251! Parametres sorties
2252  real, dimension(klon), intent(OUT):: tsurf_new, evap, fluxsens, fluxlat
2253  real, dimension(klon), intent(OUT):: dflux_s, dflux_l
2254
2255! Variables locales
2256  integer :: i
2257  real, dimension(klon) :: zx_mh, zx_nh, zx_oh
2258  real, dimension(klon) :: zx_mq, zx_nq, zx_oq
2259  real, dimension(klon) :: zx_pkh, zx_dq_s_dt, zx_qsat, zx_coef
2260  real, dimension(klon) :: zx_sl, zx_k1
2261  real, dimension(klon) :: zx_q_0 , d_ts
2262  real                  :: zdelta, zcvm5, zx_qs, zcor, zx_dq_s_dh
2263  real                  :: bilan_f, fq_fonte
2264  REAL                  :: subli, fsno
2265  REAL                  :: qsat_new, q1_new
2266  real, parameter :: t_grnd = 271.35, t_coup = 273.15
2267!! PB temporaire en attendant mieux pour le modele de neige
2268  REAL, parameter :: chasno = 3.334E+05/(2.3867E+06*0.15)
2269!
2270  logical, save         :: check = .false.
2271  character (len = 20)  :: modname = 'calcul_fluxs'
2272  logical, save         :: fonte_neige = .false.
2273  real, save            :: max_eau_sol = 150.0
2274  character (len = 80) :: abort_message
2275  logical,save         :: first = .true.,second=.false.
2276
2277  if (check) write(*,*)'Entree ', modname,' surface = ',nisurf
2278
2279  IF (check) THEN
2280      WRITE(*,*)' radsol (min, max)' &
2281         &     , MINVAL(radsol(1:knon)), MAXVAL(radsol(1:knon))
2282      CALL flush(6)
2283  ENDIF
2284
2285  if (size(coastalflow) /= knon .AND. nisurf == is_ter) then
2286    write(*,*)'Bizarre, le nombre de points continentaux'
2287    write(*,*)'a change entre deux appels. J''arrete ...'
2288    abort_message='Pb run_off'
2289    call abort_gcm(modname,abort_message,1)
2290  endif
2291!
2292! Traitement neige et humidite du sol
2293!
2294!!$  WRITE(*,*)'test calcul_flux, surface ', nisurf
2295!!PB test
2296!!$    if (nisurf == is_oce) then
2297!!$      snow = 0.
2298!!$      qsol = max_eau_sol
2299!!$    else
2300!!$      where (precip_snow > 0.) snow = snow + (precip_snow * dtime)
2301!!$      where (snow > epsilon(snow)) snow = max(0.0, snow - (evap * dtime))
2302!!$!      snow = max(0.0, snow + (precip_snow - evap) * dtime)
2303!!$      where (precip_rain > 0.) qsol = qsol + (precip_rain - evap) * dtime
2304!!$    endif
2305!!$    IF (nisurf /= is_ter) qsol = max_eau_sol
2306
2307
2308!
2309! Initialisation
2310!
2311  evap = 0.
2312  fluxsens=0.
2313  fluxlat=0.
2314  dflux_s = 0.
2315  dflux_l = 0. 
2316!
2317! zx_qs = qsat en kg/kg
2318!
2319  DO i = 1, knon
2320    zx_pkh(i) = (ps(i)/ps(i))**RKAPPA
2321  IF (thermcep) THEN
2322      zdelta=MAX(0.,SIGN(1.,rtt-tsurf(i)))
2323      zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta
2324      zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q1lay(i))
2325      zx_qs= r2es * FOEEW(tsurf(i),zdelta)/ps(i)
2326      zx_qs=MIN(0.5,zx_qs)
2327      zcor=1./(1.-retv*zx_qs)
2328      zx_qs=zx_qs*zcor
2329      zx_dq_s_dh = FOEDE(tsurf(i),zdelta,zcvm5,zx_qs,zcor) &
2330     &                 /RLVTT / zx_pkh(i)
2331    ELSE
2332      IF (tsurf(i).LT.t_coup) THEN
2333        zx_qs = qsats(tsurf(i)) / ps(i)
2334        zx_dq_s_dh = dqsats(tsurf(i),zx_qs)/RLVTT &
2335     &                    / zx_pkh(i)
2336      ELSE
2337        zx_qs = qsatl(tsurf(i)) / ps(i)
2338        zx_dq_s_dh = dqsatl(tsurf(i),zx_qs)/RLVTT &
2339     &               / zx_pkh(i)
2340      ENDIF
2341    ENDIF
2342    zx_dq_s_dt(i) = RCPD * zx_pkh(i) * zx_dq_s_dh
2343    zx_qsat(i) = zx_qs
2344    zx_coef(i) = coef1lay(i) &
2345     & * (1.0+SQRT(u1lay(i)**2+v1lay(i)**2)) &
2346     & * p1lay(i)/(RD*t1lay(i))
2347
2348  ENDDO
2349
2350
2351! === Calcul de la temperature de surface ===
2352!
2353! zx_sl = chaleur latente d'evaporation ou de sublimation
2354!
2355  do i = 1, knon
2356    zx_sl(i) = RLVTT
2357    if (tsurf(i) .LT. RTT) zx_sl(i) = RLSTT
2358    zx_k1(i) = zx_coef(i)
2359  enddo
2360
2361
2362  do i = 1, knon
2363! Q
2364    zx_oq(i) = 1. - (beta(i) * zx_k1(i) * peqBcoef(i) * dtime)
2365    zx_mq(i) = beta(i) * zx_k1(i) * &
2366     &             (peqAcoef(i) - zx_qsat(i) &
2367     &                          + zx_dq_s_dt(i) * tsurf(i)) &
2368     &             / zx_oq(i)
2369    zx_nq(i) = beta(i) * zx_k1(i) * (-1. * zx_dq_s_dt(i)) &
2370     &                              / zx_oq(i)
2371
2372! H
2373    zx_oh(i) = 1. - (zx_k1(i) * petBcoef(i) * dtime)
2374    zx_mh(i) = zx_k1(i) * petAcoef(i) / zx_oh(i)
2375    zx_nh(i) = - (zx_k1(i) * RCPD * zx_pkh(i))/ zx_oh(i)
2376
2377! Tsurface
2378    tsurf_new(i) = (tsurf(i) + cal(i)/(RCPD * zx_pkh(i)) * dtime * &
2379     &             (radsol(i) + zx_mh(i) + zx_sl(i) * zx_mq(i)) &
2380     &                 + dif_grnd(i) * t_grnd * dtime)/ &
2381     &          ( 1. - dtime * cal(i)/(RCPD * zx_pkh(i)) * ( &
2382     &                       zx_nh(i) + zx_sl(i) * zx_nq(i)) & 
2383     &                     + dtime * dif_grnd(i))
2384
2385!
2386! Y'a-t-il fonte de neige?
2387!
2388!    fonte_neige = (nisurf /= is_oce) .AND. &
2389!     & (snow(i) > epsfra .OR. nisurf == is_sic .OR. nisurf == is_lic) &
2390!     & .AND. (tsurf_new(i) >= RTT)
2391!    if (fonte_neige) tsurf_new(i) = RTT 
2392    d_ts(i) = tsurf_new(i) - tsurf(i)
2393!    zx_h_ts(i) = tsurf_new(i) * RCPD * zx_pkh(i)
2394!    zx_q_0(i) = zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)
2395!== flux_q est le flux de vapeur d'eau: kg/(m**2 s)  positive vers bas
2396!== flux_t est le flux de cpt (energie sensible): j/(m**2 s)
2397    evap(i) = - zx_mq(i) - zx_nq(i) * tsurf_new(i)
2398    fluxlat(i) = - evap(i) * zx_sl(i)
2399    fluxsens(i) = zx_mh(i) + zx_nh(i) * tsurf_new(i)
2400! Derives des flux dF/dTs (W m-2 K-1):
2401    dflux_s(i) = zx_nh(i)
2402    dflux_l(i) = (zx_sl(i) * zx_nq(i))
2403! Nouvelle valeure de l'humidite au dessus du sol
2404    qsat_new=zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)
2405    q1_new = peqAcoef(i) - peqBcoef(i)*evap(i)*dtime
2406    qsurf(i)=q1_new*(1.-beta(i)) + beta(i)*qsat_new
2407!
2408! en cas de fonte de neige
2409!
2410!    if (fonte_neige) then
2411!      bilan_f = radsol(i) + fluxsens(i) - (zx_sl(i) * evap (i)) - &
2412!     &          dif_grnd(i) * (tsurf_new(i) - t_grnd) - &
2413!     &          RCPD * (zx_pkh(i))/cal(i)/dtime * (tsurf_new(i) - tsurf(i))
2414!      bilan_f = max(0., bilan_f)
2415!      fq_fonte = bilan_f / zx_sl(i)
2416!      snow(i) = max(0., snow(i) - fq_fonte * dtime)
2417!      qsol(i) = qsol(i) + (fq_fonte * dtime)
2418!    endif
2419!!$    if (nisurf == is_ter)  &
2420!!$     &  run_off(i) = run_off(i) + max(qsol(i) - max_eau_sol, 0.0)
2421!!$    qsol(i) = min(qsol(i), max_eau_sol)
2422  ENDDO
2423
2424  END SUBROUTINE calcul_fluxs
2425!
2426!#########################################################################
2427!
2428  SUBROUTINE gath2cpl(champ_in, champ_out, klon, knon, iim, jjm, knindex)
2429
2430! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
2431! au coupleur.
2432!
2433!
2434! input:         
2435!   champ_in     champ sur la grille gathere       
2436!   knon         nombre de points dans le domaine a traiter
2437!   knindex      index des points de la surface a traiter
2438!   klon         taille de la grille
2439!   iim,jjm      dimension de la grille 2D
2440!
2441! output:
2442!   champ_out    champ sur la grille 2D
2443!
2444! input
2445  integer                   :: klon, knon, iim, jjm
2446  real, dimension(klon)     :: champ_in
2447  integer, dimension(klon)  :: knindex
2448! output
2449  real, dimension(iim,jjm+1)  :: champ_out
2450! local
2451  integer                   :: i, ig, j
2452  real, dimension(klon)     :: tamp
2453
2454  tamp = 0.
2455  do i = 1, knon
2456    ig = knindex(i)
2457    tamp(ig) = champ_in(i)
2458  enddo   
2459  ig = 1
2460  champ_out(:,1) = tamp(ig)
2461  do j = 2, jjm
2462    do i = 1, iim
2463      ig = ig + 1
2464      champ_out(i,j) = tamp(ig)
2465    enddo
2466  enddo
2467  ig = ig + 1
2468  champ_out(:,jjm+1) = tamp(ig)
2469
2470  END SUBROUTINE gath2cpl
2471!
2472!#########################################################################
2473!
2474  SUBROUTINE cpl2gath(champ_in, champ_out, klon, knon, iim, jjm, knindex)
2475
2476! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
2477! au coupleur.
2478!
2479!
2480! input:         
2481!   champ_in     champ sur la grille gathere       
2482!   knon         nombre de points dans le domaine a traiter
2483!   knindex      index des points de la surface a traiter
2484!   klon         taille de la grille
2485!   iim,jjm      dimension de la grille 2D
2486!
2487! output:
2488!   champ_out    champ sur la grille 2D
2489!
2490! input
2491  integer                   :: klon, knon, iim, jjm
2492  real, dimension(iim,jjm+1)     :: champ_in
2493  integer, dimension(klon)  :: knindex
2494! output
2495  real, dimension(klon)  :: champ_out
2496! local
2497  integer                   :: i, ig, j
2498  real, dimension(klon)     :: tamp
2499  logical ,save                  :: check = .false.
2500
2501  ig = 1
2502  tamp(ig) = champ_in(1,1)
2503  do j = 2, jjm
2504    do i = 1, iim
2505      ig = ig + 1
2506      tamp(ig) = champ_in(i,j)
2507    enddo
2508  enddo
2509  ig = ig + 1
2510  tamp(ig) = champ_in(1,jjm+1)
2511
2512  do i = 1, knon
2513    ig = knindex(i)
2514    champ_out(i) = tamp(ig)
2515  enddo   
2516
2517  END SUBROUTINE cpl2gath
2518!
2519!#########################################################################
2520!
2521  SUBROUTINE albsno(klon, knon,dtime,agesno,alb_neig_grid, precip_snow)
2522  IMPLICIT none
2523 
2524  INTEGER :: klon, knon
2525  INTEGER, PARAMETER :: nvm = 8
2526  REAL   :: dtime
2527  REAL, dimension(klon,nvm) :: veget
2528  REAL, DIMENSION(klon) :: alb_neig_grid, agesno, precip_snow
2529 
2530  INTEGER :: i, nv
2531 
2532  REAL, DIMENSION(nvm),SAVE :: init, decay
2533  REAL :: as
2534  DATA init /0.55, 0.14, 0.18, 0.29, 0.15, 0.15, 0.14, 0./
2535  DATA decay/0.30, 0.67, 0.63, 0.45, 0.40, 0.14, 0.06, 1./
2536 
2537  veget = 0.
2538  veget(:,1) = 1.     ! desert partout
2539  DO i = 1, knon
2540    alb_neig_grid(i) = 0.0
2541  ENDDO
2542  DO nv = 1, nvm
2543    DO i = 1, knon
2544      as = init(nv)+decay(nv)*EXP(-agesno(i)/5.)
2545      alb_neig_grid(i) = alb_neig_grid(i) + veget(i,nv)*as
2546    ENDDO
2547  ENDDO
2548!
2549!! modilation en fonction de l'age de la neige
2550!
2551  DO i = 1, knon
2552    agesno(i)  = (agesno(i) + (1.-agesno(i)/50.)*dtime/86400.)&
2553            &             * EXP(-1.*MAX(0.0,precip_snow(i))*dtime/0.3)
2554    agesno(i) =  MAX(agesno(i),0.0)
2555  ENDDO
2556 
2557  END SUBROUTINE albsno
2558!
2559!#########################################################################
2560!
2561
2562  SUBROUTINE fonte_neige( klon, knon, nisurf, dtime, &
2563     & tsurf, p1lay, cal, beta, coef1lay, ps, &
2564     & precip_rain, precip_snow, snow, qsol, &
2565     & radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, &
2566     & petAcoef, peqAcoef, petBcoef, peqBcoef, &
2567     & tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l, &
2568     & fqcalving,ffonte,run_off_lic_0)
2569
2570! Routine de traitement de la fonte de la neige dans le cas du traitement
2571! de sol simplifié
2572!
2573! LF 03/2001
2574! input:
2575!   knon         nombre de points a traiter
2576!   nisurf       surface a traiter
2577!   tsurf        temperature de surface
2578!   p1lay        pression 1er niveau (milieu de couche)
2579!   cal          capacite calorifique du sol
2580!   beta         evap reelle
2581!   coef1lay     coefficient d'echange
2582!   ps           pression au sol
2583!   precip_rain  precipitations liquides
2584!   precip_snow  precipitations solides
2585!   snow         champs hauteur de neige
2586!   qsol         hauteur d'eau contenu dans le sol
2587!   runoff       runoff en cas de trop plein
2588!   petAcoef     coeff. A de la resolution de la CL pour t
2589!   peqAcoef     coeff. A de la resolution de la CL pour q
2590!   petBcoef     coeff. B de la resolution de la CL pour t
2591!   peqBcoef     coeff. B de la resolution de la CL pour q
2592!   radsol       rayonnement net aus sol (LW + SW)
2593!   dif_grnd     coeff. diffusion vers le sol profond
2594!
2595! output:
2596!   tsurf_new    temperature au sol
2597!   fluxsens     flux de chaleur sensible
2598!   fluxlat      flux de chaleur latente
2599!   dflux_s      derivee du flux de chaleur sensible / Ts
2600!   dflux_l      derivee du flux de chaleur latente  / Ts
2601! in/out:
2602!   run_off_lic_0 run off glacier du pas de temps précedent
2603!
2604
2605#include "YOETHF.inc"
2606#include "FCTTRE.inc"
2607#include "indicesol.inc"
2608!IM cf JLD
2609#include "YOMCST.inc"
2610
2611! Parametres d'entree
2612  integer, intent(IN) :: knon, nisurf, klon
2613  real   , intent(IN) :: dtime
2614  real, dimension(klon), intent(IN) :: petAcoef, peqAcoef
2615  real, dimension(klon), intent(IN) :: petBcoef, peqBcoef
2616  real, dimension(klon), intent(IN) :: ps, q1lay
2617  real, dimension(klon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
2618  real, dimension(klon), intent(IN) :: precip_rain, precip_snow
2619  real, dimension(klon), intent(IN) :: radsol, dif_grnd
2620  real, dimension(klon), intent(IN) :: t1lay, u1lay, v1lay
2621  real, dimension(klon), intent(INOUT) :: snow, qsol
2622
2623! Parametres sorties
2624  real, dimension(klon), intent(INOUT):: tsurf_new, evap, fluxsens, fluxlat
2625  real, dimension(klon), intent(INOUT):: dflux_s, dflux_l
2626! Flux thermique utiliser pour fondre la neige
2627  real, dimension(klon), intent(INOUT):: ffonte
2628! Flux d'eau "perdue" par la surface et necessaire pour que limiter la
2629! hauteur de neige, en kg/m2/s
2630  real, dimension(klon), intent(INOUT):: fqcalving
2631  real, dimension(klon), intent(INOUT):: run_off_lic_0
2632! Variables locales
2633! Masse maximum de neige (kg/m2). Au dessus de ce seuil, la neige
2634! en exces "s'ecoule" (calving)
2635!  real, parameter :: snow_max=1.
2636!IM cf JLD/GK
2637  real, parameter :: snow_max=3000.
2638  integer :: i
2639  real, dimension(klon) :: zx_mh, zx_nh, zx_oh
2640  real, dimension(klon) :: zx_mq, zx_nq, zx_oq
2641  real, dimension(klon) :: zx_pkh, zx_dq_s_dt, zx_qsat, zx_coef
2642  real, dimension(klon) :: zx_sl, zx_k1
2643  real, dimension(klon) :: zx_q_0 , d_ts
2644  real                  :: zdelta, zcvm5, zx_qs, zcor, zx_dq_s_dh
2645  real                  :: bilan_f, fq_fonte
2646  REAL                  :: subli, fsno
2647  REAL, DIMENSION(klon) :: bil_eau_s, snow_evap
2648  real, parameter :: t_grnd = 271.35, t_coup = 273.15
2649!! PB temporaire en attendant mieux pour le modele de neige
2650! REAL, parameter :: chasno = RLMLT/(2.3867E+06*0.15)
2651  REAL, parameter :: chasno = 3.334E+05/(2.3867E+06*0.15)
2652!IM cf JLD/ GKtest
2653  REAL, parameter :: chaice = 3.334E+05/(2.3867E+06*0.15)
2654! fin GKtest
2655!
2656  logical, save         :: check = .FALSE.
2657  character (len = 20)  :: modname = 'fonte_neige'
2658  logical, save         :: neige_fond = .false.
2659  real, save            :: max_eau_sol = 150.0
2660  character (len = 80) :: abort_message
2661  logical,save         :: first = .true.,second=.false.
2662  real                 :: coeff_rel
2663
2664
2665  if (check) write(*,*)'Entree ', modname,' surface = ',nisurf
2666
2667! Initialisations
2668  coeff_rel = dtime/(tau_calv * rday)
2669  bil_eau_s(:) = 0.
2670  DO i = 1, knon
2671    zx_pkh(i) = (ps(i)/ps(i))**RKAPPA
2672  IF (thermcep) THEN
2673      zdelta=MAX(0.,SIGN(1.,rtt-tsurf(i)))
2674      zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta
2675      zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q1lay(i))
2676      zx_qs= r2es * FOEEW(tsurf(i),zdelta)/ps(i)
2677      zx_qs=MIN(0.5,zx_qs)
2678      zcor=1./(1.-retv*zx_qs)
2679      zx_qs=zx_qs*zcor
2680      zx_dq_s_dh = FOEDE(tsurf(i),zdelta,zcvm5,zx_qs,zcor) &
2681     &                 /RLVTT / zx_pkh(i)
2682    ELSE
2683      IF (tsurf(i).LT.t_coup) THEN
2684        zx_qs = qsats(tsurf(i)) / ps(i)
2685        zx_dq_s_dh = dqsats(tsurf(i),zx_qs)/RLVTT &
2686     &                    / zx_pkh(i)
2687      ELSE
2688        zx_qs = qsatl(tsurf(i)) / ps(i)
2689        zx_dq_s_dh = dqsatl(tsurf(i),zx_qs)/RLVTT &
2690     &               / zx_pkh(i)
2691      ENDIF
2692    ENDIF
2693    zx_dq_s_dt(i) = RCPD * zx_pkh(i) * zx_dq_s_dh
2694    zx_qsat(i) = zx_qs
2695    zx_coef(i) = coef1lay(i) &
2696     & * (1.0+SQRT(u1lay(i)**2+v1lay(i)**2)) &
2697     & * p1lay(i)/(RD*t1lay(i))
2698  ENDDO
2699
2700
2701! === Calcul de la temperature de surface ===
2702!
2703! zx_sl = chaleur latente d'evaporation ou de sublimation
2704!
2705  do i = 1, knon
2706    zx_sl(i) = RLVTT
2707    if (tsurf(i) .LT. RTT) zx_sl(i) = RLSTT
2708    zx_k1(i) = zx_coef(i)
2709  enddo
2710
2711
2712  do i = 1, knon
2713! Q
2714    zx_oq(i) = 1. - (beta(i) * zx_k1(i) * peqBcoef(i) * dtime)
2715    zx_mq(i) = beta(i) * zx_k1(i) * &
2716     &             (peqAcoef(i) - zx_qsat(i) &
2717     &                          + zx_dq_s_dt(i) * tsurf(i)) &
2718     &             / zx_oq(i)
2719    zx_nq(i) = beta(i) * zx_k1(i) * (-1. * zx_dq_s_dt(i)) &
2720     &                              / zx_oq(i)
2721
2722! H
2723    zx_oh(i) = 1. - (zx_k1(i) * petBcoef(i) * dtime)
2724    zx_mh(i) = zx_k1(i) * petAcoef(i) / zx_oh(i)
2725    zx_nh(i) = - (zx_k1(i) * RCPD * zx_pkh(i))/ zx_oh(i)
2726  enddo
2727
2728
2729  WHERE (precip_snow > 0.) snow = snow + (precip_snow * dtime)
2730  snow_evap = 0.
2731  WHERE (evap > 0. )
2732    snow_evap = MIN (snow / dtime, evap)
2733    snow = snow - snow_evap * dtime
2734    snow = MAX(0.0, snow)
2735  end where
2736 
2737!  bil_eau_s = bil_eau_s + (precip_rain * dtime) - (evap - snow_evap) * dtime
2738  bil_eau_s = (precip_rain * dtime) - (evap - snow_evap) * dtime
2739
2740!
2741! Y'a-t-il fonte de neige?
2742!
2743  ffonte=0.
2744  do i = 1, knon
2745    neige_fond = ((snow(i) > epsfra .OR. nisurf == is_sic .OR. nisurf == is_lic) &
2746     & .AND. tsurf_new(i) >= RTT)
2747    if (neige_fond) then
2748      fq_fonte = MIN( MAX((tsurf_new(i)-RTT )/chasno,0.0),snow(i))
2749      ffonte(i) = fq_fonte * RLMLT/dtime
2750      snow(i) = max(0., snow(i) - fq_fonte)
2751      bil_eau_s(i) = bil_eau_s(i) + fq_fonte
2752      tsurf_new(i) = tsurf_new(i) - fq_fonte * chasno 
2753!IM cf JLD OK     
2754!IM cf JLD/ GKtest fonte aussi pour la glace
2755      IF (nisurf == is_sic .OR. nisurf == is_lic ) THEN
2756        fq_fonte = MAX((tsurf_new(i)-RTT )/chaice,0.0)
2757        ffonte(i) = ffonte(i) + fq_fonte * RLMLT/dtime
2758        bil_eau_s(i) = bil_eau_s(i) + fq_fonte
2759        tsurf_new(i) = RTT
2760      ENDIF
2761      d_ts(i) = tsurf_new(i) - tsurf(i)
2762    endif
2763!
2764!   s'il y a une hauteur trop importante de neige, elle s'coule
2765    fqcalving(i) = max(0., snow(i) - snow_max)/dtime
2766    snow(i)=min(snow(i),snow_max)
2767!
2768    IF (nisurf == is_ter) then
2769      qsol(i) = qsol(i) + bil_eau_s(i)
2770      run_off(i) = run_off(i) + MAX(qsol(i) - max_eau_sol, 0.0)
2771      qsol(i) = MIN(qsol(i), max_eau_sol)
2772    else if (nisurf == is_lic) then
2773      run_off_lic(i) = (coeff_rel *  fqcalving(i)) + &
2774 &                        (1. - coeff_rel) * run_off_lic_0(i)
2775      run_off_lic_0(i) = run_off_lic(i)
2776      run_off_lic(i) = run_off_lic(i) + bil_eau_s(i)/dtime
2777    endif
2778  enddo
2779
2780  END SUBROUTINE fonte_neige
2781!
2782!#########################################################################
2783!
2784  END MODULE interface_surf
Note: See TracBrowser for help on using the repository browser.