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

Last change on this file since 524 was 524, checked in by lmdzadmin, 20 years ago

Initial revision

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