source: LMDZ4/trunk/libf/phylmd/interface_surf.F90 @ 689

Last change on this file since 689 was 686, checked in by lmdzadmin, 19 years ago

Appel du slab lorsque nisurf=is_sic
IM

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