source: LMDZ.3.3/branches/rel-LF/libf/phylmd/interface_surf.F90 @ 98

Last change on this file since 98 was 98, checked in by lmdzadmin, 24 years ago

Interface avec les differentes surface, version de travail.LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 47.9 KB
Line 
1
2  MODULE interface_surf
3
4! Ce module regroupe toutes les routines gerant l'interface entre le modele
5! atmospherique et les modeles de surface (sols continentaux, oceans, glaces)
6! Les routines sont les suivantes:
7!
8!   interfsurf_*: routines d'aiguillage vers les interfaces avec les
9!                 differents modeles de surface
10!   interfsol\
11!             > routines d'interface proprement dite
12!   interfoce/
13!
14!   interfstart: routine d'initialisation et de lecture de l'etat initial
15!                "interface"
16!   interffin  : routine d'ecriture de l'etat de redemmarage de l'interface
17!
18!
19! L. Fairhead, LMD, 02/2000
20
21  USE ioipsl
22!  USE constantes
23
24  IMPLICIT none
25
26  PRIVATE
27  PUBLIC :: interfsurf,interfsurf_hq
28
29  INTERFACE interfsurf
30    module procedure interfsurf_hq, interfsurf_vent
31  END INTERFACE
32
33  INTERFACE interfoce
34    module procedure interfoce_cpl, interfoce_slab, interfoce_lim
35  END INTERFACE
36
37
38! run_off      ruissellement total
39  real, allocatable, dimension(:),save    :: run_off
40
41
42  CONTAINS
43!
44!############################################################################
45!
46  SUBROUTINE interfsurf_hq(itime, dtime, jour, &
47      & klon, iim, jjm, nisurf, knon, knindex, pctsrf, rlon, rlat, &
48      & debut, lafin, ok_veget, &
49      & zlev,  u1_lay, v1_lay, temp_air, spechum, hum_air, ccanopy, &
50      & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
51      & precip_rain, precip_snow, lwdown, swnet, swdown, &
52      & albedo, snow, qsol, &
53      & tsurf, p1lay, coef1lay, ps, radsol, &
54      & ocean, &
55      & evap, fluxsens, fluxlat, dflux_l, dflux_s, &             
56      & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, pctsrf_new, zmasq)
57
58
59! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
60! (sols continentaux, oceans, glaces) pour les fluxs de chaleur et d'humidite.
61! En pratique l'interface se fait entre la couche limite du modele
62! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
63!
64!
65! L.Fairhead 02/2000
66!
67! input:
68!   itime        numero du pas de temps
69!   klon         nombre total de points de grille
70!   iim, jjm     nbres de pts de grille
71!   dtime        pas de temps de la physique (en s)
72!   jour         jour dans l'annee en cours
73!   nisurf       index de la surface a traiter (1 = sol continental)
74!   knon         nombre de points de la surface a traiter
75!   knindex      index des points de la surface a traiter
76!   pctsrf       tableau des pourcentages de surface de chaque maille
77!   rlon         longitudes
78!   rlat         latitudes
79!   debut        logical: 1er appel a la physique
80!   lafin        logical: dernier appel a la physique
81!   ok_veget     logical: appel ou non au schema de surface continental
82!                     (si false calcul simplifie des fluxs sur les continents)
83!   zlev         hauteur de la premiere couche
84!   u1_lay       vitesse u 1ere couche
85!   v1_lay       vitesse v 1ere couche
86!   temp_air     temperature de l'air 1ere couche
87!   spechum      humidite specifique 1ere couche
88!   hum_air      humidite de l'air
89!   ccanopy      concentration CO2 canopee
90!   tq_cdrag     cdrag
91!   petAcoef     coeff. A de la resolution de la CL pour t
92!   peqAcoef     coeff. A de la resolution de la CL pour q
93!   petBcoef     coeff. B de la resolution de la CL pour t
94!   peqBcoef     coeff. B de la resolution de la CL pour q
95!   precip_rain  precipitation liquide
96!   precip_snow  precipitation solide
97!   lwdown       flux IR entrant a la surface
98!   swnet        flux solaire net
99!   swdown       flux solaire entrant a la surface
100!   albedo       albedo de la surface
101!   tsurf        temperature de surface
102!   p1lay        pression 1er niveau (milieu de couche)
103!   coef1lay     coefficient d'echange
104!   ps           pression au sol
105!   radsol       rayonnement net aus sol (LW + SW)
106!   ocean        type d'ocean utilise (force, slab, couple)
107!
108! output:
109!   evap         evaporation totale
110!   fluxsens     flux de chaleur sensible
111!   fluxlat      flux de chaleur latente
112!   tsol_rad     
113!   tsurf_new    temperature au sol
114!   alb_new      albedo
115!   emis_new     emissivite
116!   z0_new       surface roughness
117!   pctsrf_new   nouvelle repartition des surfaces
118!   zmasq        masque terre/ocean
119
120  include 'indicesol.h'
121
122! Parametres d'entree
123  integer, intent(IN) :: itime
124  integer, intent(IN) :: iim, jjm
125  integer, intent(IN) :: klon
126  real, intent(IN) :: dtime
127  integer, intent(IN) :: jour
128  integer, intent(IN) :: nisurf
129  integer, intent(IN) :: knon
130  integer, dimension(knon), intent(in) :: knindex
131  real, dimension(klon,nbsrf), intent(IN) :: pctsrf
132  logical, intent(IN) :: debut, lafin, ok_veget
133  real, dimension(klon), intent(IN) :: rlon, rlat
134  real, dimension(knon), intent(IN) :: zlev
135  real, dimension(knon), intent(IN) :: u1_lay, v1_lay
136  real, dimension(knon), intent(IN) :: temp_air, spechum
137  real, dimension(knon), intent(IN) :: hum_air, ccanopy
138  real, dimension(knon), intent(IN) :: tq_cdrag, petAcoef, peqAcoef
139  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
140  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
141  real, dimension(knon), intent(IN) :: lwdown, swnet, swdown, ps, albedo
142  real, dimension(knon), intent(IN) :: tsurf, p1lay, coef1lay
143  real, dimension(knon), intent(IN) :: radsol
144  real, dimension(klon), intent(IN) :: zmasq
145  character (len = 6)  :: ocean
146  real, dimension(knon), intent(INOUT) :: evap, snow, qsol
147
148! Parametres de sortie
149  real, dimension(knon), intent(OUT):: fluxsens, fluxlat
150  real, dimension(knon), intent(OUT):: tsol_rad, tsurf_new, alb_new
151  real, dimension(knon), intent(OUT):: emis_new, z0_new
152  real, dimension(knon), intent(OUT):: dflux_l, dflux_s
153  real, dimension(klon,nbsrf), intent(OUT) :: pctsrf_new
154
155! Local
156  character (len = 20) :: modname = 'interfsurf_hq'
157  character (len = 80) :: abort_message
158  logical, save        :: first_call = .true.
159  integer              :: error
160  logical              :: check = .true.
161  real, dimension(knon):: cal, beta, dif_grnd, capsol
162  real, parameter      :: calice=1.0/(5.1444e+06*0.15), tau_gl=1./86400.*5.
163  real, parameter      :: calsno=1./(2.3867e+06*.15)
164
165#include "YOMCST.inc"
166
167  if (check) write(*,*) 'Entree ', modname
168!
169! On doit commencer par appeler les schemas de surfaces continentales
170! car l'ocean a besoin du ruissellement qui est y calcule
171!
172  if (first_call) then
173    if (nisurf /= is_ter .and. klon > 1) then
174      write(*,*)' *** Warning ***'
175      write(*,*)' nisurf = ',nisurf,' /= is_ter = ',is_ter
176      write(*,*)'or on doit commencer par les surfaces continentales'
177      abort_message='voir ci-dessus'
178      call abort_gcm(modname,abort_message,1)
179    endif
180    if (ocean /= 'slab  ' .and. ocean /= 'force ' .and. ocean /= 'couple') then
181      write(*,*)' *** Warning ***'
182      write(*,*)'Option couplage pour l''ocean = ', ocean
183      abort_message='option pour l''ocean non valable'
184      call abort_gcm(modname,abort_message,1)
185    endif
186  endif
187  first_call = .false.
188!
189! Calcul age de la neige
190!
191 
192! Aiguillage vers les differents schemas de surface
193
194  if (nisurf == is_ter) then
195!
196! Surface "terre" appel a l'interface avec les sols continentaux
197!
198! allocation du run-off
199    if (.not. allocated(run_off)) then
200      allocate(run_off(knon), stat = error)
201      if (error /= 0) then
202        abort_message='Pb allocation run_off'
203        call abort_gcm(modname,abort_message,1)
204      endif
205    else if (size(run_off) /= knon) then
206      write(*,*)'Bizarre, le nombre de points continentaux'
207      write(*,*)'a change entre deux appels. Je continue ...'
208      deallocate(run_off, stat = error)
209      allocate(run_off(knon), stat = error)
210      if (error /= 0) then
211        abort_message='Pb allocation run_off'
212        call abort_gcm(modname,abort_message,1)
213      endif
214    endif
215!
216    if (.not. ok_veget) then
217!
218! calcul snow et qsol, hydrol adapté
219!
220      call calbeta(dtime, nisurf, knon, snow, qsol, beta, capsol, dif_grnd)
221      cal = RCPD * capsol
222      call calcul_fluxs( knon, nisurf, dtime, &
223     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
224     &   precip_rain, precip_snow, snow, qsol,  &
225     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
226     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
227     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
228
229!
230! calcul albedo: lecture albedo fichier CL puis ajout albedo neige
231!
232    else
233!
234!  appel a sechiba
235!
236      call interfsol(itime, klon, dtime, nisurf, knon, &
237     &  knindex, rlon, rlat, &
238     &  debut, lafin, ok_veget, &
239     &  zlev,  u1_lay, v1_lay, temp_air, spechum, hum_air, ccanopy, &
240     &  tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
241     &  precip_rain, precip_snow, lwdown, swnet, swdown, &
242     &  tsurf, p1lay, coef1lay, ps, radsol, &
243     &  evap, fluxsens, fluxlat, &             
244     &  tsol_rad, tsurf_new, alb_new, emis_new, z0_new, dflux_l, dflux_s)
245    endif   
246!
247  else if (nisurf == is_oce) then
248
249    if (check) write(*,*)'ocean, nisurf = ',nisurf
250
251
252!
253! Surface "ocean" appel a l'interface avec l'ocean
254!
255!    if (ocean == 'couple') then
256!      call interfoce(nisurf, ocean)
257!    else if (ocean == 'slab  ') then
258!      call interfoce(nisurf)
259!    else                              ! lecture conditions limites
260!      call interfoce(itime, dtime, jour, &
261!     &  klon, nisurf, knon, knindex, &
262!     &  debut, &
263!     &  tsurf_new, alb_new, z0_new, pctsrf_new)
264!    endif
265!
266    cal = 0.
267    beta = 1.
268    dif_grnd = 0.
269
270    call calcul_fluxs( knon, nisurf, dtime, &
271     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
272     &   precip_rain, precip_snow, snow, qsol,  &
273     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
274     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
275     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
276
277
278!
279! calcul albedo
280!
281
282!
283  else if (nisurf == is_sic) then
284
285    if (check) write(*,*)'sea ice, nisurf = ',nisurf
286
287!
288! Surface "glace de mer" appel a l'interface avec l'ocean
289!
290!    call interfoce(nisurf, ocean)
291!
292
293    cal = calice
294    where (snow > 0.0) cal = calsno
295    beta = 1.0
296    dif_grnd = 1.0 / tau_gl
297
298    call calcul_fluxs( knon, nisurf, dtime, &
299     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
300     &   precip_rain, precip_snow, snow, qsol,  &
301     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
302     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
303     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
304
305
306  else if (nisurf == is_lic) then
307
308    if (check) write(*,*)'glacier, nisurf = ',nisurf
309
310!
311! Surface "glacier continentaux" appel a l'interface avec le sol
312!
313!    call interfsol(nisurf)
314
315    cal = calice
316    where (snow > 0.0) cal = calsno
317    beta = 1.0
318    dif_grnd = 0.0
319
320    call calcul_fluxs( knon, nisurf, dtime, &
321     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
322     &   precip_rain, precip_snow, snow, qsol,  &
323     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
324     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
325     &   tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
326
327  else
328    write(*,*)'Index surface = ',nisurf
329    abort_message = 'Index surface non valable'
330    call abort_gcm(modname,abort_message,1)
331  endif
332
333  END SUBROUTINE interfsurf_hq
334
335!
336!#########################################################################
337!
338  SUBROUTINE interfsurf_vent(nisurf, knon   &         
339  &                     )
340!
341! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
342! (sols continentaux, oceans, glaces) pour les tensions de vents.
343! En pratique l'interface se fait entre la couche limite du modele
344! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
345!
346!
347! L.Fairhead 02/2000
348!
349! input:
350!   nisurf       index de la surface a traiter (1 = sol continental)
351!   knon         nombre de points de la surface a traiter
352
353! Parametres d'entree
354  integer, intent(IN) :: nisurf
355  integer, intent(IN) :: knon
356
357
358  return
359  END SUBROUTINE interfsurf_vent
360!
361!#########################################################################
362!
363  SUBROUTINE interfsol(itime, klon, dtime, nisurf, knon, &
364     & knindex, rlon, rlat, &
365     & debut, lafin, ok_veget, &
366     & zlev,  u1_lay, v1_lay, temp_air, spechum, hum_air, ccanopy, &
367     & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
368     & precip_rain, precip_snow, lwdown, swnet, swdown, &
369     & tsurf, p1lay, coef1lay, ps, radsol, &
370     & evap, fluxsens, fluxlat, &             
371     & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, dflux_l, dflux_s)
372
373! Cette routine sert d'interface entre le modele atmospherique et le
374! modele de sol continental. Appel a sechiba
375!
376! L. Fairhead 02/2000
377!
378! input:
379!   itime        numero du pas de temps
380!   klon         nombre total de points de grille
381!   dtime        pas de temps de la physique (en s)
382!   nisurf       index de la surface a traiter (1 = sol continental)
383!   knon         nombre de points de la surface a traiter
384!   knindex      index des points de la surface a traiter
385!   rlon         longitudes de la grille entiere
386!   rlat         latitudes de la grille entiere
387!   debut        logical: 1er appel a la physique (lire les restart)
388!   lafin        logical: dernier appel a la physique (ecrire les restart)
389!   ok_veget     logical: appel ou non au schema de surface continental
390!                     (si false calcul simplifie des fluxs sur les continents)
391!   zlev         hauteur de la premiere couche       
392!   u1_lay       vitesse u 1ere couche
393!   v1_lay       vitesse v 1ere couche
394!   temp_air     temperature de l'air 1ere couche
395!   spechum      humidite specifique 1ere couche
396!   hum_air      humidite de l'air
397!   ccanopy      concentration CO2 canopee
398!   tq_cdrag     cdrag
399!   petAcoef     coeff. A de la resolution de la CL pour t
400!   peqAcoef     coeff. A de la resolution de la CL pour q
401!   petBcoef     coeff. B de la resolution de la CL pour t
402!   peqBcoef     coeff. B de la resolution de la CL pour q
403!   precip_rain  precipitation liquide
404!   precip_snow  precipitation solide
405!   lwdown       flux IR entrant a la surface
406!   swnet        flux solaire net
407!   swdown       flux solaire entrant a la surface
408!   tsurf        temperature de surface
409!   p1lay        pression 1er niveau (milieu de couche)
410!   coef1lay     coefficient d'echange
411!   ps           pression au sol
412!   radsol       rayonnement net aus sol (LW + SW)
413!   
414!
415! input/output
416!   run_off      ruissellement total
417!
418! output:
419!   evap         evaporation totale
420!   fluxsens     flux de chaleur sensible
421!   fluxlat      flux de chaleur latente
422!   tsol_rad     
423!   tsurf_new    temperature au sol
424!   alb_new      albedo
425!   emis_new     emissivite
426!   z0_new       surface roughness
427
428
429! Parametres d'entree
430  integer, intent(IN) :: itime
431  integer, intent(IN) :: klon
432  real, intent(IN)    :: dtime
433  integer, intent(IN) :: nisurf
434  integer, intent(IN) :: knon
435  integer, dimension(knon), intent(IN) :: knindex
436  logical, intent(IN) :: debut, lafin, ok_veget
437  real, dimension(klon), intent(IN) :: rlon, rlat
438  real, dimension(knon), intent(IN) :: zlev
439  real, dimension(knon), intent(IN) :: u1_lay, v1_lay
440  real, dimension(knon), intent(IN) :: temp_air, spechum
441  real, dimension(knon), intent(IN) :: hum_air, ccanopy
442  real, dimension(knon), intent(IN) :: tq_cdrag, petAcoef, peqAcoef
443  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
444  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
445  real, dimension(knon), intent(IN) :: lwdown, swnet, swdown, ps
446  real, dimension(knon), intent(IN) :: tsurf, p1lay, coef1lay
447  real, dimension(knon), intent(IN) :: radsol
448! Parametres de sortie
449  real, dimension(knon), intent(OUT):: evap, fluxsens, fluxlat
450  real, dimension(knon), intent(OUT):: tsol_rad, tsurf_new, alb_new
451  real, dimension(knon), intent(OUT):: emis_new, z0_new
452  real, dimension(knon), intent(OUT):: dflux_s, dflux_l
453
454! Local
455!
456  integer              :: ii
457  integer              :: error
458  character (len = 20) :: modname = 'interfsol'
459  character (len = 80) :: abort_message
460  logical              :: check = .true.
461  real, dimension(knon) :: cal, beta, dif_grnd, capsol
462! type de couplage dans sechiba
463!  character (len=10)   :: coupling = 'implicit'
464! drapeaux controlant les appels dans SECHIBA
465!  type(control_type), save   :: control_in
466! coordonnees geographiques
467  real, allocatable, dimension(:,:), save :: lalo
468! pts voisins
469  integer,allocatable, dimension(:,:), save :: neighbours
470! resolution de la grille
471  real, allocatable, dimension (:,:), save :: resolution
472! Identifieurs des fichiers restart et histoire
473  integer, save          :: rest_id, hist_id
474  integer, save          :: rest_id_stom, hist_id_stom
475
476  real, dimension(knon):: snow, qsol
477
478  if (check) write(*,*)'Entree ', modname
479  if (check) write(*,*)'ok_veget = ',ok_veget
480
481! initialisation
482!    if (debut) then
483!      !
484!      ! Configuration de parametres specifiques a la SSL
485!      !
486!      call intsurf_config(control_in)
487!      !
488!      ! Allouer et initialiser le tableau de coordonnees du sol
489!      !
490!      if (( .not. allocated(lalo))) then
491!        allocate(lalo(knon,2), stat = error)
492!        if (error /= 0) then
493!          abort_message='Pb allocation lalo'
494!          call abort_gcm(modname,abort_message,1)
495!        endif     
496!      endif
497!      do ii = 1, knon
498!        lalo(ii,1) = rlat(knindex(ii))
499!        lalo(ii,2) = rlon(knindex(ii))
500!      enddo
501      !-
502      !- Compute variable to help describe the grid
503      !- once the points are gathered.
504      !-
505!      IF ( (.NOT.ALLOCATED(neighbours))) THEN
506!        ALLOCATE(neighbours(knon,4), stat = error)
507!        if (error /= 0) then
508!          abort_message='Pb allocation neighbours'
509!          call abort_gcm(modname,abort_message,1)
510!        endif
511!      ENDIF
512!      IF ( (.NOT.ALLOCATED(resolution))) THEN
513!        ALLOCATE(resolution(knon,2), stat = error)
514!        if (error /= 0) then
515!          abort_message='Pb allocation resolution'
516!          call abort_gcm(modname,abort_message,1)
517!        endif
518!      ENDIF
519
520! call grid_stuff
521! call sechiba_restart_init
522! call sechiba_history_init
523
524!    endif                          ! (fin debut)
525
526!
527! Appel a la routine sols continentaux
528!
529
530!    call sechiba_main(itime, klon, knon, knindex, dtime, &
531!     & debut, lafin, coupling, control_in, &
532!     & lalo, neighbours, resolution,&
533!     & zlev,  u1_lay, v1_lay, spechum, temp_air,hum_air , ccanopy, &
534!     & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
535!     & precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
536!     & evap, fluxsens, fluxlat, &
537!     & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, &
538!     & rest_id, hist_id, rest_id_stom, hist_id_stom)
539
540!
541! Sauvegarde dans fichiers histoire
542!
543
544  END SUBROUTINE interfsol
545!
546!#########################################################################
547!
548  SUBROUTINE interfoce_cpl(itime, dtime, &
549      & klon, iim, jjm, nisurf, pctsrf, knon, knindex, rlon, rlat, &
550      & ocean, nexca, debut, lafin, &
551      & swdown, lwdown, precip_rain, precip_snow, evap, tsurf, &
552      & fder, albsol, taux, tauy, &
553      & tsurf_new, alb_new, alb_ice, pctsrf_new)
554
555
556
557! Cette routine sert d'interface entre le modele atmospherique et un
558! coupleur avec un modele d'ocean 'complet'
559!
560! L. Fairhead 02/2000
561!
562! input:
563!   itime        numero du pas de temps
564!   iim, jjm     nbres de pts de grille
565!   dtime        pas de tempsde la physique
566!   klon         nombre total de points de grille
567!   nisurf       index de la surface a traiter (1 = sol continental)
568!   pctsrf       tableau des fractions de surface de chaque maille
569!   knon         nombre de points de la surface a traiter
570!   knindex      index des points de la surface a traiter
571!   rlon         longitudes
572!   rlat         latitudes
573!   debut        logical: 1er appel a la physique
574!   lafin        logical: dernier appel a la physique
575!   ocean        type d'ocean
576!   nexca        frequence de couplage
577!   swdown       flux solaire entrant a la surface
578!   lwdown       flux IR entrant a la surface
579!   precip_rain  precipitation liquide
580!   precip_snow  precipitation solide
581!   evap         evaporation
582!   tsurf        temperature de surface
583!   fder         derivee dF/dT
584!   albsol       albedo du sol (coherent avec swdown)
585!   taux         tension de vent en x
586!   tauy         tension de vent en y
587!   nexca        frequence de couplage
588!
589!
590! output:
591!   tsurf_new    temperature au sol
592!   alb_new      albedo
593!   pctsrf_new   nouvelle repartition des surfaces
594!   alb_ice      albedo de la glace
595!
596
597#include 'indicesol.h'
598
599! Parametres d'entree
600  integer, intent(IN) :: itime
601  integer, intent(IN) :: iim, jjm
602  real, intent(IN) :: dtime
603  integer, intent(IN) :: klon
604  integer, intent(IN) :: nisurf
605  integer, intent(IN) :: knon
606  real, dimension(klon,nbsrf), intent(IN) :: pctsrf
607  integer, dimension(knon), intent(in) :: knindex
608  logical, intent(IN) :: debut, lafin
609  real, dimension(klon), intent(IN) :: rlon, rlat
610  character (len = 6)  :: ocean
611  real, dimension(knon), intent(IN) :: lwdown, swdown
612  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
613  real, dimension(knon), intent(IN) :: tsurf, fder, albsol, taux, tauy
614  integer              :: nexca
615
616  real, dimension(knon), intent(INOUT) :: evap
617
618! Parametres de sortie
619  real, dimension(knon), intent(OUT):: tsurf_new, alb_new, alb_ice
620  real, dimension(klon,nbsrf), intent(OUT) :: pctsrf_new
621
622! Variables locales
623  integer                    :: j, error, sum_error, ig
624  integer                    :: npas
625  character (len = 20) :: modname = 'interfoce_cpl'
626  character (len = 80) :: abort_message
627  logical              :: check = .true.
628! variables pour moyenner les variables de couplage
629  real, allocatable, dimension(:),save :: cpl_sols, cpl_nsol, cpl_rain
630  real, allocatable, dimension(:),save :: cpl_snow, cpl_evap, cpl_tsol
631  real, allocatable, dimension(:),save :: cpl_fder, cpl_albe, cpl_taux
632  real, allocatable, dimension(:),save :: cpl_tauy, cpl_ruis
633! variables a passer au coupleur
634  real, dimension(iim, jjm+1) :: wri_sols, wri_nsol, wri_rain
635  real, dimension(iim, jjm+1) :: wri_snow, wri_evap, wri_tsol
636  real, dimension(iim, jjm+1) :: wri_fder, wri_albe, wri_taux
637  real, dimension(iim, jjm+1) :: wri_tauy, wri_ruis
638! variables relues par le coupleur
639  real, dimension(iim, jjm+1) :: read_sst, read_sic
640  real, dimension(iim, jjm+1) :: read_alb_sst, read_alb_sic
641! variable tampon
642  real, dimension(klon)       :: tamp
643  real, dimension(knon)       :: tamp_sic
644
645
646!
647! Initialisation
648!
649  if (debut) then
650    sum_error = 0
651    allocate(cpl_sols(knon), stat = error)
652    sum_error = sum_error + error
653    allocate(cpl_nsol(knon), stat = error)
654    sum_error = sum_error + error
655    allocate(cpl_rain(knon), stat = error)
656    sum_error = sum_error + error
657    allocate(cpl_snow(knon), stat = error)
658    sum_error = sum_error + error
659    allocate(cpl_evap(knon), stat = error)
660    sum_error = sum_error + error
661    allocate(cpl_tsol(knon), stat = error)
662    sum_error = sum_error + error
663    allocate(cpl_fder(knon), stat = error)
664    sum_error = sum_error + error
665    allocate(cpl_albe(knon), stat = error)
666    sum_error = sum_error + error
667    allocate(cpl_taux(knon), stat = error)
668    sum_error = sum_error + error
669    allocate(cpl_tauy(knon), stat = error)
670    sum_error = sum_error + error
671    allocate(cpl_ruis(knon), stat = error)
672    sum_error = sum_error + error
673    if (sum_error /= 0) then
674      abort_message='Pb allocation variables couplees'
675      call abort_gcm(modname,abort_message,1)
676    endif
677    cpl_sols = 0.
678    cpl_nsol = 0.
679    cpl_rain = 0.
680    cpl_snow = 0.
681    cpl_evap = 0.
682    cpl_tsol = 0.
683    cpl_fder = 0.
684    cpl_albe = 0.
685    cpl_taux = 0.
686    cpl_tauy = 0.
687    cpl_ruis = 0.
688!
689! initialisation couplage
690!
691    call inicma(npas, nexca, dtime)
692!
693! 1ere lecture champs ocean
694!
695    call fromcpl(itime,(jjm+1)*iim,                                          &
696     &        read_sst, read_sic, read_alb_sst, read_alb_sic)
697    call cpl2gath(read_sst, tsurf_new, klon, knon,iim,jjm, knindex)
698    call cpl2gath(read_sic, tamp_sic , klon, knon,iim,jjm, knindex)
699    call cpl2gath(read_alb_sst, alb_new, klon, knon,iim,jjm, knindex)
700    call cpl2gath(read_alb_sic, alb_ice, klon, knon,iim,jjm, knindex)
701!
702! transformer tamp_sic en pctsrf_new
703!
704    do ig = 1, klon
705      IF (pctsrf(ig,is_oce) > epsfra .OR.            &
706     &             pctsrf(ig,is_sic) > epsfra) THEN
707            pctsrf_new(ig,is_oce) = pctsrf(ig,is_oce)    &
708     &                        - (tamp_sic(ig)-pctsrf(ig,is_sic))
709            pctsrf_new(ig,is_sic) = tamp_sic(ig)
710      endif
711    enddo
712
713  endif ! fin if (debut)
714
715! fichier restart et fichiers histoires
716
717! calcul des fluxs a passer
718
719  cpl_sols = cpl_sols + swdown      / FLOAT(nexca)
720  cpl_nsol = cpl_nsol + lwdown      / FLOAT(nexca)
721  cpl_rain = cpl_rain + precip_rain / FLOAT(nexca)
722  cpl_snow = cpl_snow + precip_snow / FLOAT(nexca)
723  cpl_evap = cpl_evap + evap        / FLOAT(nexca)
724  cpl_tsol = cpl_tsol + tsurf       / FLOAT(nexca)
725  cpl_fder = cpl_fder + fder        / FLOAT(nexca)
726  cpl_albe = cpl_albe + albsol      / FLOAT(nexca)
727  cpl_taux = cpl_taux + taux        / FLOAT(nexca)
728  cpl_tauy = cpl_tauy + tauy        / FLOAT(nexca)
729  cpl_ruis = cpl_ruis + run_off     / FLOAT(nexca)/dtime
730
731  if (mod(itime, nexca) == 0) then
732!
733! Mise sur la bonne grille des champs a passer au coupleur
734    call gath2cpl(cpl_sols, wri_sols, klon, knon,iim,jjm, knindex)
735    call gath2cpl(cpl_nsol, wri_nsol, klon, knon,iim,jjm, knindex)
736    call gath2cpl(cpl_rain, wri_rain, klon, knon,iim,jjm, knindex)
737    call gath2cpl(cpl_snow, wri_snow, klon, knon,iim,jjm, knindex)
738    call gath2cpl(cpl_evap, wri_evap, klon, knon,iim,jjm, knindex)
739    call gath2cpl(cpl_tsol, wri_tsol, klon, knon,iim,jjm, knindex)
740    call gath2cpl(cpl_fder, wri_fder, klon, knon,iim,jjm, knindex)
741    call gath2cpl(cpl_albe, wri_albe, klon, knon,iim,jjm, knindex)
742    call gath2cpl(cpl_taux, wri_taux, klon, knon,iim,jjm, knindex)
743    call gath2cpl(cpl_tauy, wri_tauy, klon, knon,iim,jjm, knindex)
744    call gath2cpl(cpl_ruis, wri_ruis, klon, knon,iim,jjm, knindex)
745!
746! Passage des champs au coupleur
747!
748    call intocpl(itime, iim, jjm , wri_sols, wri_nsol, wri_rain, wri_snow, &
749      & wri_evap, wri_tsol, wri_fder, wri_albe, wri_taux, wri_tauy,       &
750      & wri_ruis )
751    cpl_sols = 0.
752    cpl_nsol = 0.
753    cpl_rain = 0.
754    cpl_snow = 0.
755    cpl_evap = 0.
756    cpl_tsol = 0.
757    cpl_fder = 0.
758    cpl_albe = 0.
759    cpl_taux = 0.
760    cpl_tauy = 0.
761    cpl_ruis = 0.
762
763    call fromcpl(itime,(jjm+1)*iim,                                          &
764     &        read_sst, read_sic, read_alb_sst, read_alb_sic)
765    call cpl2gath(read_sst, tsurf_new, klon, knon,iim,jjm, knindex)
766    call cpl2gath(read_sic, tamp_sic , klon, knon,iim,jjm, knindex)
767    call cpl2gath(read_alb_sst, alb_new, klon, knon,iim,jjm, knindex)
768    call cpl2gath(read_alb_sic, alb_ice, klon, knon,iim,jjm, knindex)
769! transformer tamp_sic en pctsrf_new
770
771    do ig = 1, klon
772      IF (pctsrf(ig,is_oce) > epsfra .OR.            &
773     &             pctsrf(ig,is_sic) > epsfra) THEN
774            pctsrf_new(ig,is_oce) = pctsrf(ig,is_oce)    &
775     &                        - (tamp_sic(ig)-pctsrf(ig,is_sic))
776            pctsrf_new(ig,is_sic) = tamp_sic(ig)
777      endif
778    enddo
779  endif
780
781  END SUBROUTINE interfoce_cpl
782!
783!#########################################################################
784!
785
786  SUBROUTINE interfoce_slab(nisurf)
787
788! Cette routine sert d'interface entre le modele atmospherique et un
789! modele de 'slab' ocean
790!
791! L. Fairhead 02/2000
792!
793! input:
794!   nisurf       index de la surface a traiter (1 = sol continental)
795!
796! output:
797!
798
799! Parametres d'entree
800  integer, intent(IN) :: nisurf
801
802  END SUBROUTINE interfoce_slab
803!
804!#########################################################################
805!
806  SUBROUTINE interfoce_lim(itime, dtime, jour, &
807     & klon, nisurf, knon, knindex, &
808     & debut,  &
809     & lmt_sst, lmt_alb, lmt_rug, pctsrf_new)
810
811! Cette routine sert d'interface entre le modele atmospherique et un fichier
812! de conditions aux limites
813!
814! L. Fairhead 02/2000
815!
816! input:
817!   itime        numero du pas de temps courant
818!   dtime        pas de temps de la physique (en s)
819!   jour         jour a lire dans l'annee
820!   nisurf       index de la surface a traiter (1 = sol continental)
821!   knon         nombre de points dans le domaine a traiter
822!   knindex      index des points de la surface a traiter
823!   klon         taille de la grille
824!   debut        logical: 1er appel a la physique (initialisation)
825!
826! output:
827!   lmt_sst      SST lues dans le fichier de CL
828!   lmt_alb      Albedo lu
829!   lmt_rug      longueur de rugosité lue
830!   pctsrf_new   sous-maille fractionnelle
831!
832
833#include "indicesol.h"
834
835! Parametres d'entree
836  integer, intent(IN) :: itime
837  real   , intent(IN) :: dtime
838  integer, intent(IN) :: jour
839  integer, intent(IN) :: nisurf
840  integer, intent(IN) :: knon
841  integer, intent(IN) :: klon
842  integer, dimension(knon), intent(in) :: knindex
843  logical, intent(IN) :: debut
844
845! Parametres de sortie
846  real, intent(out), dimension(knon) :: lmt_sst
847  real, intent(out), dimension(knon) :: lmt_alb
848  real, intent(out), dimension(knon) :: lmt_rug
849  real, intent(out), dimension(klon,nbsrf) :: pctsrf_new
850
851! Variables locales
852  integer     :: ii
853  integer     :: lmt_pas     ! frequence de lecture des conditions limites
854                             ! (en pas de physique)
855  logical,save :: deja_lu    ! pour indiquer que le jour a lire a deja
856                             ! lu pour une surface precedente
857  integer,save :: jour_lu
858  integer      :: ierr
859  character (len = 20) :: modname = 'interfoce_lim'
860  character (len = 80) :: abort_message
861  character (len = 20) :: fich ='limit'
862  logical     :: newlmt = .false.
863  logical     :: check = .true.
864! Champs lus dans le fichier de CL
865  real, allocatable , save, dimension(:) :: sst_lu, alb_lu, rug_lu, nat_lu
866  real, allocatable , save, dimension(:,:) :: pct_tmp
867!
868! quelques variables pour netcdf
869!
870#include "netcdf.inc"
871  integer              :: nid, nvarid
872  integer, dimension(2) :: start, epais
873!
874! Fin déclaration
875!
876   
877  if (debut) then
878    lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
879    jour_lu = jour - 1
880    allocate(sst_lu(klon))
881    allocate(alb_lu(klon))
882    allocate(rug_lu(klon))
883    allocate(nat_lu(klon))
884    allocate(pct_tmp(klon,nbsrf))
885  endif
886
887  if ((jour - jour_lu) /= 0) deja_lu = .false.
888 
889  if (check) write(*,*)modname,' :: jour_lu, deja_lu', jour_lu, deja_lu
890
891! Tester d'abord si c'est le moment de lire le fichier
892  if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu) then
893!
894! Ouverture du fichier
895!
896    ierr = NF_OPEN (fich, NF_NOWRITE,nid)
897    if (ierr.NE.NF_NOERR) then
898      abort_message = 'Pb d''ouverture du fichier de conditions aux limites'
899      call abort_gcm(modname,abort_message,1)
900    endif
901!
902! La tranche de donnees a lire:
903!
904    start(1) = 1
905    start(2) = jour + 1
906    epais(1) = klon
907    epais(2) = 1
908!
909    if (newlmt) then
910!
911! Fraction "ocean"
912!
913      ierr = NF_INQ_VARID(nid, 'FOCE', nvarid)
914      if (ierr /= NF_NOERR) then
915        abort_message = 'Le champ <FOCE> est absent'
916        call abort_gcm(modname,abort_message,1)
917      endif
918#ifdef NC_DOUBLE
919      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_oce))
920#else
921      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_oce))
922#endif
923      if (ierr /= NF_NOERR) then
924        abort_message = 'Lecture echouee pour <FOCE>'
925        call abort_gcm(modname,abort_message,1)
926      endif
927!
928! Fraction "glace de mer"
929!
930      ierr = NF_INQ_VARID(nid, 'FSIC', nvarid)
931      if (ierr /= NF_NOERR) then
932        abort_message = 'Le champ <FSIC> est absent'
933        call abort_gcm(modname,abort_message,1)
934      endif
935#ifdef NC_DOUBLE
936      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_sic))
937#else
938      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_sic))
939#endif
940      if (ierr /= NF_NOERR) then
941        abort_message = 'Lecture echouee pour <FSIC>'
942        call abort_gcm(modname,abort_message,1)
943      endif
944!
945! Fraction "terre"
946!
947      ierr = NF_INQ_VARID(nid, 'FTER', nvarid)
948      if (ierr /= NF_NOERR) then
949        abort_message = 'Le champ <FTER> est absent'
950        call abort_gcm(modname,abort_message,1)
951      endif
952#ifdef NC_DOUBLE
953      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_ter))
954#else
955      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_ter))
956#endif
957      if (ierr /= NF_NOERR) then
958        abort_message = 'Lecture echouee pour <FTER>'
959        call abort_gcm(modname,abort_message,1)
960      endif
961!
962! Fraction "glacier terre"
963!
964      ierr = NF_INQ_VARID(nid, 'FLIC', nvarid)
965      if (ierr /= NF_NOERR) then
966        abort_message = 'Le champ <FLIC> est absent'
967        call abort_gcm(modname,abort_message,1)
968      endif
969#ifdef NC_DOUBLE
970      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_lic))
971#else
972      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_lic))
973#endif
974      if (ierr /= NF_NOERR) then
975        abort_message = 'Lecture echouee pour <FLIC>'
976        call abort_gcm(modname,abort_message,1)
977      endif
978!
979    else  ! on en est toujours a rnatur
980!
981      ierr = NF_INQ_VARID(nid, 'NAT', nvarid)
982      if (ierr /= NF_NOERR) then
983        abort_message = 'Le champ <NAT> est absent'
984        call abort_gcm(modname,abort_message,1)
985      endif
986#ifdef NC_DOUBLE
987      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, nat_lu)
988#else
989      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, nat_lu)
990#endif
991      if (ierr /= NF_NOERR) then
992        abort_message = 'Lecture echouee pour <NAT>'
993        call abort_gcm(modname,abort_message,1)
994      endif
995!
996! Remplissage des fractions de surface
997! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice
998!
999      pct_tmp = 0.0
1000      do ii = 1, klon
1001        pct_tmp(ii,nint(nat_lu(ii)) + 1) = 1.
1002      enddo
1003
1004!
1005!  On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire
1006!
1007      pctsrf_new = pct_tmp
1008      pctsrf_new (:,2)= pct_tmp (:,1)
1009      pctsrf_new (:,1)= pct_tmp (:,2)
1010      pct_tmp = pctsrf_new
1011    endif ! fin test sur newlmt
1012!
1013! Lecture SST
1014!
1015    ierr = NF_INQ_VARID(nid, 'SST', nvarid)
1016    if (ierr /= NF_NOERR) then
1017      abort_message = 'Le champ <SST> est absent'
1018      call abort_gcm(modname,abort_message,1)
1019    endif
1020#ifdef NC_DOUBLE
1021    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, sst_lu)
1022#else
1023    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, sst_lu)
1024#endif
1025    if (ierr /= NF_NOERR) then
1026      abort_message = 'Lecture echouee pour <SST>'
1027      call abort_gcm(modname,abort_message,1)
1028    endif   
1029!
1030! Lecture Albedo
1031!
1032    ierr = NF_INQ_VARID(nid, 'ALB', nvarid)
1033    if (ierr /= NF_NOERR) then
1034      abort_message = 'Le champ <ALB> est absent'
1035      call abort_gcm(modname,abort_message,1)
1036    endif
1037#ifdef NC_DOUBLE
1038    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, alb_lu)
1039#else
1040    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, alb_lu)
1041#endif
1042    if (ierr /= NF_NOERR) then
1043      abort_message = 'Lecture echouee pour <ALB>'
1044      call abort_gcm(modname,abort_message,1)
1045    endif
1046!
1047! Lecture rugosité
1048!
1049    ierr = NF_INQ_VARID(nid, 'RUG', nvarid)
1050    if (ierr /= NF_NOERR) then
1051      abort_message = 'Le champ <RUG> est absent'
1052      call abort_gcm(modname,abort_message,1)
1053    endif
1054#ifdef NC_DOUBLE
1055    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, rug_lu)
1056#else
1057    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, rug_lu)
1058#endif
1059    if (ierr /= NF_NOERR) then
1060      abort_message = 'Lecture echouee pour <RUG>'
1061      call abort_gcm(modname,abort_message,1)
1062    endif
1063
1064!
1065! Fin de lecture
1066!
1067    ierr = NF_CLOSE(nid)
1068    deja_lu = .true.
1069    jour_lu = jour
1070  endif
1071!
1072! Recopie des variables dans les champs de sortie
1073!
1074  do ii = 1, knon
1075    lmt_sst(ii) = sst_lu(knindex(ii))
1076    lmt_alb(ii) = alb_lu(knindex(ii))
1077    lmt_rug(ii) = rug_lu(knindex(ii))
1078  enddo
1079  pctsrf_new = pct_tmp
1080
1081  END SUBROUTINE interfoce_lim
1082
1083!
1084!#########################################################################
1085!
1086
1087  SUBROUTINE calcul_fluxs( knon, nisurf, dtime, &
1088     & tsurf, p1lay, cal, beta, coef1lay, ps, &
1089     & precip_rain, precip_snow, snow, qsol, &
1090     & radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, &
1091     & petAcoef, peqAcoef, petBcoef, peqBcoef, &
1092     & tsurf_new, evap, fluxlat, fluxsens, dflux_s, dflux_l)
1093
1094! Cette routine calcule les fluxs en h et q a l'interface et eventuellement
1095! une temperature de surface (au cas ou ok_veget = false)
1096!
1097! L. Fairhead 4/2000
1098!
1099! input:
1100!   knon         nombre de points a traiter
1101!   nisurf       surface a traiter
1102!   tsurf        temperature de surface
1103!   p1lay        pression 1er niveau (milieu de couche)
1104!   cal          capacite calorifique du sol
1105!   beta         evap reelle
1106!   coef1lay     coefficient d'echange
1107!   ps           pression au sol
1108!   precip_rain  precipitations liquides
1109!   precip_snow  precipitations solides
1110!   snow         champs hauteur de neige
1111!   qsol         humidite du sol
1112!   runoff       runoff en cas de trop plein
1113!   petAcoef     coeff. A de la resolution de la CL pour t
1114!   peqAcoef     coeff. A de la resolution de la CL pour q
1115!   petBcoef     coeff. B de la resolution de la CL pour t
1116!   peqBcoef     coeff. B de la resolution de la CL pour q
1117!   radsol       rayonnement net aus sol (LW + SW)
1118!   dif_grnd     coeff. diffusion vers le sol profond
1119!
1120! output:
1121!   tsurf_new    temperature au sol
1122!   fluxsens     flux de chaleur sensible
1123!   fluxlat      flux de chaleur latente
1124!   dflux_s      derivee du flux de chaleur sensible / Ts
1125!   dflux_l      derivee du flux de chaleur latente  / Ts
1126!
1127
1128#include "YOMCST.inc"
1129#include "YOETHF.inc"
1130#include "FCTTRE.inc"
1131#include 'indicesol.h'
1132
1133! Parametres d'entree
1134  integer, intent(IN) :: knon, nisurf
1135  real   , intent(IN) :: dtime
1136  real, dimension(knon), intent(IN) :: petAcoef, peqAcoef
1137  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
1138  real, dimension(knon), intent(IN) :: ps, q1lay
1139  real, dimension(knon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
1140  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
1141  real, dimension(knon), intent(IN) :: radsol, dif_grnd
1142  real, dimension(knon), intent(IN) :: t1lay, u1lay, v1lay
1143  real, dimension(knon), intent(INOUT) :: snow, qsol
1144
1145! Parametres sorties
1146  real, dimension(knon), intent(OUT):: tsurf_new, evap, fluxsens, fluxlat
1147  real, dimension(knon), intent(OUT):: dflux_s, dflux_l
1148
1149! Variables locales
1150  integer :: i
1151  real, dimension(knon) :: zx_mh, zx_nh, zx_oh
1152  real, dimension(knon) :: zx_mq, zx_nq, zx_oq
1153  real, dimension(knon) :: zx_pkh, zx_dq_s_dt, zx_qsat, zx_coef
1154  real, dimension(knon) :: zx_sl, zx_k1,  zx_dq,  zx_cq,  zx_dh, zx_ch
1155  real, dimension(knon) :: zx_h_ts, zx_q_0 , d_ts
1156  real                  :: zdelta, zcvm5, zx_qs, zcor, zx_dq_s_dh
1157  real                  :: bilan_f, fq_fonte
1158  real, parameter :: t_grnd = 271.35, t_coup = 273.15
1159  logical         :: check = .true.
1160  character (len = 20)  :: modname = 'calcul_fluxs'
1161  logical         :: fonte_neige = .false.
1162  real            :: max_eau_sol = 150.0
1163  character (len = 80) :: abort_message
1164
1165  if (check) write(*,*)'Entree ', modname
1166
1167  if (size(run_off) /= knon .AND. nisurf == is_ter) then
1168    write(*,*)'Bizarre, le nombre de points continentaux'
1169    write(*,*)'a change entre deux appels. J''arrete ...'
1170    abort_message='Pb run_off'
1171    call abort_gcm(modname,abort_message,1)
1172  endif
1173!
1174! Traitement neige et humidite du sol
1175!
1176    if (nisurf == is_oce) then
1177      snow = 0.
1178      qsol = max_eau_sol
1179    else
1180      snow = max(0.0, snow + (precip_snow - evap) * dtime)
1181      qsol = qsol + (precip_rain - evap) * dtime
1182    endif
1183
1184
1185!
1186! Initialisation
1187!
1188!
1189! zx_qs = qsat en kg/kg
1190!
1191  DO i = 1, knon
1192    zx_pkh(i) = (ps(i)/ps(i))**RKAPPA
1193  IF (thermcep) THEN
1194      zdelta=MAX(0.,SIGN(1.,rtt-tsurf(i)))
1195      zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta
1196      zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q1lay(i))
1197      zx_qs= r2es * FOEEW(tsurf(i),zdelta)/ps(i)
1198      zx_qs=MIN(0.5,zx_qs)
1199      zcor=1./(1.-retv*zx_qs)
1200      zx_qs=zx_qs*zcor
1201      zx_dq_s_dh = FOEDE(tsurf(i),zdelta,zcvm5,zx_qs,zcor) &
1202     &                 /RLVTT / zx_pkh(i)
1203    ELSE
1204      IF (tsurf(i).LT.t_coup) THEN
1205        zx_qs = qsats(tsurf(i)) / ps(i)
1206        zx_dq_s_dh = dqsats(tsurf(i),zx_qs)/RLVTT &
1207     &                    / zx_pkh(i)
1208      ELSE
1209        zx_qs = qsatl(tsurf(i)) / ps(i)
1210        zx_dq_s_dh = dqsatl(tsurf(i),zx_qs)/RLVTT &
1211     &               / zx_pkh(i)
1212      ENDIF
1213    ENDIF
1214    zx_dq_s_dt(i) = RCPD * zx_pkh(i) * zx_dq_s_dh
1215    zx_qsat(i) = zx_qs
1216    zx_coef(i) = coef1lay(i) &
1217     & * (1.0+SQRT(u1lay(i)**2+v1lay(i)**2)) &
1218     & * p1lay(i)/(RD*t1lay(i))
1219   
1220  ENDDO
1221
1222
1223! === Calcul de la temperature de surface ===
1224!
1225! zx_sl = chaleur latente d'evaporation ou de sublimation
1226!
1227  do i = 1, knon
1228    zx_sl(i) = RLVTT
1229    if (tsurf(i) .LT. RTT) zx_sl(i) = RLSTT
1230    zx_k1(i) = zx_coef(i)
1231  enddo
1232
1233
1234  do i = 1, knon
1235! Q
1236    zx_oq(i) = 1. - (beta(i) * zx_k1(i) * peqBcoef(i) * dtime)
1237    zx_mq(i) = beta(i) * zx_k1(i) * &
1238     &             (peqAcoef(i) - zx_qsat(i) &
1239     &                          + zx_dq_s_dt(i) * tsurf(i)) &
1240     &             / zx_oq(i)
1241    zx_nq(i) = beta(i) * zx_k1(i) * (-1. * zx_dq_s_dt(i)) &
1242     &                              / zx_oq(i)
1243
1244! H
1245    zx_oh(i) = 1. - (zx_k1(i) * petBcoef(i) * dtime)
1246    zx_mh(i) = zx_k1(i) * petAcoef(i) / zx_oh(i)
1247    zx_nh(i) = - (zx_k1(i) * RCPD * zx_pkh(i))/ zx_oh(i)
1248
1249! Tsurface
1250    tsurf_new(i) = (tsurf(i) + cal(i)/(RCPD * zx_pkh(i)) * dtime * &
1251     &             (radsol(i) + zx_mh(i) + zx_sl(i) * zx_mq(i)) &
1252     &                 + dif_grnd(i) * t_grnd * dtime)/ &
1253     &          ( 1. - dtime * cal(i)/(RCPD * zx_pkh(i)) * ( &
1254     &                       zx_nh(i) + zx_sl(i) * zx_nq(i)) & 
1255     &                     + dtime * dif_grnd(i))
1256
1257!
1258! Y'a-t-il fonte de neige?
1259!
1260    fonte_neige = (nisurf /= is_oce) .AND. &
1261     & (snow(i) > 0. .OR. nisurf == is_sic .OR. nisurf == is_lic) &
1262     & .AND. (tsurf_new(i) >= RTT)
1263    if (fonte_neige) tsurf_new(i) = RTT 
1264    zx_h_ts(i) = tsurf_new(i) * RCPD * zx_pkh(i)
1265    d_ts(i) = tsurf_new(i) - tsurf(i)
1266    zx_q_0(i) = zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)
1267!== flux_q est le flux de vapeur d'eau: kg/(m**2 s)  positive vers bas
1268!== flux_t est le flux de cpt (energie sensible): j/(m**2 s)
1269    evap(i) = - zx_mq(i) - zx_nq(i) * tsurf_new(i)
1270    fluxlat(i) = - evap(i) * zx_sl(i)
1271    fluxsens(i) = zx_mh(i) + zx_nh(i) * tsurf_new(i)
1272! Derives des flux dF/dTs (W m-2 K-1):
1273    dflux_s(i) = zx_nh(i)
1274    dflux_l(i) = (zx_sl(i) * zx_nq(i))
1275!
1276! en cas de fonte de neige
1277!
1278    if (fonte_neige) then
1279      bilan_f = radsol(i) + fluxsens(i) - (zx_sl(i) * evap (i)) - &
1280     &          dif_grnd(i) * (tsurf_new(i) - t_grnd) - &
1281     &          RCPD * (zx_pkh(i))/cal(i)/dtime * (tsurf_new(i) - tsurf(i))
1282      bilan_f = max(0., bilan_f)
1283      fq_fonte = bilan_f / zx_sl(i)
1284      snow(i) = max(0., snow(i) - fq_fonte * dtime)
1285      qsol(i) = qsol(i) + (fq_fonte * dtime)
1286    endif
1287    if (nisurf == is_ter)  &
1288     &  run_off(i) = run_off(i) + max(qsol(i) - max_eau_sol, 0.0)
1289    qsol(i) = min(qsol(i), max_eau_sol)
1290  ENDDO
1291
1292  END SUBROUTINE calcul_fluxs
1293!
1294!#########################################################################
1295!
1296
1297  SUBROUTINE sol_dem_write(itime, klon, rlon, rlat, &
1298  &                     pctsrf_new,tsurf_new,alb_new)
1299
1300! Routine d'ecriture de l'etat de redemarrage pour le sol
1301!
1302! L.Fairhead
1303!
1304! input:
1305!   itime        numero du pas de temps
1306!   klon         nombre total de points de grille
1307!   rlon         longitudes
1308!   rlat         latitudes
1309!   tsurf_new    temperature au sol
1310!   alb_new      albedo
1311!   pctsrf_new   repartition des surfaces
1312
1313  include 'indicesol.h'
1314#include 'temps.inc'
1315  include 'netcdf.inc'
1316
1317! Parametres d'entree
1318  integer, intent(IN) :: itime
1319  integer, intent(IN) :: klon
1320  real, dimension(klon), intent(IN) :: rlon, rlat
1321  real, dimension(klon,nbsrf), intent(IN)  :: tsurf_new, alb_new
1322  real, dimension(klon,nbsrf), intent(IN) :: pctsrf_new
1323
1324! Variables locales
1325  integer             :: ierr, nid
1326  integer             :: idim1, idim2, idim3
1327  integer,parameter   :: length = 100
1328  character (len = 20) :: modname = 'sol_dem_write'
1329  character (len = 80) :: abort_message
1330  real, dimension(length) :: tab_cntrl = 0.
1331  integer                 :: nvarid
1332
1333  ierr = NF_CREATE('restartsol', NF_CLOBBER, nid)
1334  IF (ierr.NE.NF_NOERR) THEN
1335    abort_message=' Pb d''ouverture du fichier restartsol'
1336    CALL abort_gcm(modname,abort_message,ierr)
1337  ENDIF
1338
1339  ierr = NF_PUT_ATT_TEXT (nid, NF_GLOBAL, "title", 23,    &
1340     &                    "Fichier redemmarage sol")
1341  ierr = NF_DEF_DIM (nid, "index", length, idim1)
1342  ierr = NF_DEF_DIM (nid, "points_physiques", klon, idim2)
1343  ierr = NF_DEF_DIM (nid, "nombre_surfaces", nbsrf, idim3)
1344  ierr = NF_ENDDEF(nid)
1345
1346  tab_cntrl(13) = day_end
1347  tab_cntrl(14) = anne_ini
1348
1349  ierr = NF_REDEF (nid)
1350#ifdef NC_DOUBLE
1351  ierr = NF_DEF_VAR (nid, "controle", NF_DOUBLE, 1, idim1,nvarid)
1352#else
1353  ierr = NF_DEF_VAR (nid, "controle", NF_FLOAT, 1, idim1,nvarid)
1354#endif
1355  ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 22,  &
1356     &                        "Parametres de controle")
1357  ierr = NF_ENDDEF(nid)
1358#ifdef NC_DOUBLE
1359  ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,tab_cntrl)
1360#else
1361  ierr = NF_PUT_VAR_REAL (nid,nvarid,tab_cntrl)
1362#endif
1363
1364  ierr = NF_REDEF (nid)
1365#ifdef NC_DOUBLE
1366  ierr = NF_DEF_VAR (nid, "longitude", NF_DOUBLE, 1, idim2,nvarid)
1367#else
1368  ierr = NF_DEF_VAR (nid, "longitude", NF_FLOAT, 1, idim2,nvarid)
1369#endif
1370  ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 32,  &
1371     &               "Longitudes de la grille physique")
1372  ierr = NF_ENDDEF(nid)
1373#ifdef NC_DOUBLE
1374  ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlon)
1375#else
1376  ierr = NF_PUT_VAR_REAL (nid,nvarid,rlon)
1377#endif
1378!
1379  ierr = NF_REDEF (nid)
1380#ifdef NC_DOUBLE
1381  ierr = NF_DEF_VAR (nid, "latitude", NF_DOUBLE, 1, idim2,nvarid)
1382#else
1383  ierr = NF_DEF_VAR (nid, "latitude", NF_FLOAT, 1, idim2,nvarid)
1384#endif
1385  ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 31,    &
1386     &                        "Latitudes de la grille physique")
1387  ierr = NF_ENDDEF(nid)
1388#ifdef NC_DOUBLE
1389  ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlat)
1390#else
1391  ierr = NF_PUT_VAR_REAL (nid,nvarid,rlat)
1392#endif
1393        ierr = NF_REDEF (nid)
1394#ifdef NC_DOUBLE
1395        ierr = NF_DEF_VAR (nid, "TS", NF_DOUBLE, 1, idim2,nvarid)
1396#else
1397        ierr = NF_DEF_VAR (nid, "TS", NF_FLOAT, 1, idim2,nvarid)
1398#endif
1399        ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 22, &
1400     &                        "Temperature de surface")
1401        ierr = NF_ENDDEF(nid)
1402
1403
1404
1405
1406  END SUBROUTINE sol_dem_write
1407!
1408!#########################################################################
1409!
1410  SUBROUTINE gath2cpl(champ_in, champ_out, klon, knon, iim, jjm, knindex)
1411
1412! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
1413! au coupleur.
1414!
1415!
1416! input:         
1417!   champ_in     champ sur la grille gathere       
1418!   knon         nombre de points dans le domaine a traiter
1419!   knindex      index des points de la surface a traiter
1420!   klon         taille de la grille
1421!   iim,jjm      dimension de la grille 2D
1422!
1423! output:
1424!   champ_out    champ sur la grille 2D
1425!
1426! input
1427  integer                   :: klon, knon, iim, jjm
1428  real, dimension(knon)     :: champ_in
1429  integer, dimension(knon)  :: knindex
1430! output
1431  real, dimension(iim,jjm+1)  :: champ_out
1432! local
1433  integer                   :: i, ig, j
1434  real, dimension(klon)     :: tamp
1435
1436  do i = 1, knon
1437    ig = knindex(i)
1438    tamp(ig) = champ_in(i)
1439  enddo   
1440  champ_out(:,1) = tamp(1)
1441  do j = 2, jjm
1442    do i = 1, iim
1443      champ_out(i,j) = tamp((j-2)*jjm + i + 1)
1444    enddo
1445  enddo
1446  champ_out(:,jjm+1) = tamp(klon)
1447
1448  END SUBROUTINE gath2cpl
1449!
1450!#########################################################################
1451!
1452  SUBROUTINE cpl2gath(champ_in, champ_out, klon, knon, iim, jjm, knindex)
1453
1454! Cette routine ecrit un champ 'gathered' sur la grille 2D pour le passer
1455! au coupleur.
1456!
1457!
1458! input:         
1459!   champ_in     champ sur la grille gathere       
1460!   knon         nombre de points dans le domaine a traiter
1461!   knindex      index des points de la surface a traiter
1462!   klon         taille de la grille
1463!   iim,jjm      dimension de la grille 2D
1464!
1465! output:
1466!   champ_out    champ sur la grille 2D
1467!
1468! input
1469  integer                   :: klon, knon, iim, jjm
1470  real, dimension(iim,jjm+1)     :: champ_in
1471  integer, dimension(knon)  :: knindex
1472! output
1473  real, dimension(knon)  :: champ_out
1474! local
1475  integer                   :: i, ig, j
1476  real, dimension(klon)     :: tamp
1477
1478  tamp(1) = champ_in(1,1)
1479  do j = 2, jjm
1480    do i = 1, iim
1481      tamp((j-2)*jjm + i + 1) = champ_in(i,j)
1482    enddo
1483  enddo
1484  tamp(klon) = champ_in(1,jjm+1)
1485
1486  do i = 1, knon
1487    ig = knindex(i)
1488    champ_out(i) = tamp(ig)
1489  enddo   
1490
1491  END SUBROUTINE cpl2gath
1492!
1493!#########################################################################
1494!
1495  END MODULE interface_surf
Note: See TracBrowser for help on using the repository browser.