source: LMDZ.3.3/trunk/libf/phylmd/interface_surf.F90 @ 90

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

Version de base de la routine d'interface avec la surface. Les
aiguillages
sont faits, il suffit de rajouter les routines decrivant les differents
types de sol. Pour l'instant, seul un calcul de flux a la surface, se
trouvant auparavant dans clmain, est effectue

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 31.7 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, nisurf, knon, knindex, rlon, rlat, &
48      & debut, lafin, ok_veget, &
49      & zlev, zlflu, 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      & tsurf, p1lay, cal, beta, coef1lay, ps, radsol, dif_grnd, &
53      & ocean, &
54      & evap, fluxsens, fluxlat, dflux_l, dflux_s, &             
55      & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, pctsrf_new)
56
57
58! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
59! (sols continentaux, oceans, glaces) pour les fluxs de chaleur et d'humidite.
60! En pratique l'interface se fait entre la couche limite du modele
61! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
62!
63!
64! L.Fairhead 02/2000
65!
66! input:
67!   itime        numero du pas de temps
68!   klon         nombre total de points de grille
69!   dtime        pas de temps de la physique (en s)
70!   jour         jour dans l'annee en cours
71!   nisurf       index de la surface a traiter (1 = sol continental)
72!   knon         nombre de points de la surface a traiter
73!   knindex      index des points de la surface a traiter
74!   rlon         longitudes
75!   rlat         latitudes
76!   debut        logical: 1er appel a la physique
77!   lafin        logical: dernier appel a la physique
78!   ok_veget     logical: appel ou non au schema de surface continental
79!                     (si false calcul simplifie des fluxs sur les continents)
80!   zlev         hauteur de la premiere couche
81!   zlflu        epaisseur de la premier couche
82!   u1_lay       vitesse u 1ere couche
83!   v1_lay       vitesse v 1ere couche
84!   temp_air     temperature de l'air 1ere couche
85!   spechum      humidite specifique 1ere couche
86!   hum_air      humidite de l'air
87!   ccanopy      concentration CO2 canopee
88!   tq_cdrag     cdrag
89!   petAcoef     coeff. A de la resolution de la CL pour t
90!   peqAcoef     coeff. A de la resolution de la CL pour q
91!   petBcoef     coeff. B de la resolution de la CL pour t
92!   peqBcoef     coeff. B de la resolution de la CL pour q
93!   precip_rain  precipitation liquide
94!   precip_snow  precipitation solide
95!   lwdown       flux IR entrant a la surface
96!   swnet        flux solaire net
97!   swdown       flux solaire entrant a la surface
98!   tsurf        temperature de surface
99!   p1lay        pression 1er niveau (milieu de couche)
100!   cal          capacite calorifique du sol
101!   beta         evap reelle
102!   coef1lay     coefficient d'echange
103!   ps           pression au sol
104!   radsol       rayonnement net aus sol (LW + SW)
105!   dif_grnd     coeff. diffusion vers le sol profond
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
118  include 'indicesol.h'
119
120! Parametres d'entree
121  integer, intent(IN) :: itime
122  integer, intent(IN) :: klon
123  real, intent(IN) :: dtime
124  integer, intent(IN) :: jour
125  integer, intent(IN) :: nisurf
126  integer, intent(IN) :: knon
127  integer, dimension(knon), intent(in) :: knindex
128  logical, intent(IN) :: debut, lafin, ok_veget
129  real, dimension(klon), intent(IN) :: rlon, rlat
130  real, dimension(knon), intent(IN) :: zlev, zlflu
131  real, dimension(knon), intent(IN) :: u1_lay, v1_lay
132  real, dimension(knon), intent(IN) :: temp_air, spechum
133  real, dimension(knon), intent(IN) :: hum_air, ccanopy
134  real, dimension(knon), intent(IN) :: tq_cdrag, petAcoef, peqAcoef
135  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
136  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
137  real, dimension(knon), intent(IN) :: lwdown, swnet, swdown, ps
138  real, dimension(knon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
139  real, dimension(knon), intent(IN) :: radsol, dif_grnd
140  character (len = 6)  :: ocean
141
142! Parametres de sortie
143  real, dimension(knon), intent(OUT):: evap, fluxsens, fluxlat
144  real, dimension(knon), intent(OUT):: tsol_rad, tsurf_new, alb_new
145  real, dimension(knon), intent(OUT):: emis_new, z0_new
146  real, dimension(knon), intent(OUT):: dflux_l, dflux_s
147  real, dimension(klon,nbsrf), intent(OUT) :: pctsrf_new
148
149! Local
150  character (len = 20) :: modname = 'interfsurf_hq'
151  character (len = 80) :: abort_message
152  logical, save        :: first_call = .true.
153  integer              :: error
154  logical              :: check = .true.
155
156  if (check) write(*,*) 'Entree ', modname
157!
158! On doit commencer par appeler les schemas de surfaces continentales
159! car l'ocean a besoin du ruissellement qui est y calcule
160!
161  if (first_call) then
162    if (nisurf /= is_ter .and. klon > 1) then
163      write(*,*)' *** Warning ***'
164      write(*,*)' nisurf = ',nisurf,' /= is_ter = ',is_ter
165      write(*,*)'or on doit commencer par les surfaces continentales'
166      abort_message='voir ci-dessus'
167      call abort_gcm(modname,abort_message,1)
168    endif
169    if (ocean /= 'slab  ' .and. ocean /= 'force ' .and. ocean /= 'couple') then
170      write(*,*)' *** Warning ***'
171      write(*,*)'Option couplage pour l''ocean = ', ocean
172      abort_message='option pour l''ocean non valable'
173      call abort_gcm(modname,abort_message,1)
174    endif
175  endif
176  first_call = .false.
177
178! Aiguillage vers les differents schemas de surface
179
180  if (nisurf == is_ter) then
181!
182! Surface "terre" appel a l'interface avec les sols continentaux
183!
184! allocation du run-off
185    if (.not. allocated(run_off)) then
186      allocate(run_off(knon), stat = error)
187      if (error /= 0) then
188        abort_message='Pb allocation run_off'
189        call abort_gcm(modname,abort_message,1)
190      endif
191    else if (size(run_off) /= knon) then
192      write(*,*)'Bizarre, le nombre de points continentaux'
193      write(*,*)'a change entre deux appels. Je continue ...'
194      deallocate(run_off, stat = error)
195      allocate(run_off(knon), stat = error)
196      if (error /= 0) then
197        abort_message='Pb allocation run_off'
198        call abort_gcm(modname,abort_message,1)
199      endif
200    endif
201!
202    call interfsol(itime, klon, dtime, nisurf, knon, &
203     &  knindex, rlon, rlat, &
204     &  debut, lafin, ok_veget, &
205     &  zlev, zlflu, u1_lay, v1_lay, temp_air, spechum, hum_air, ccanopy, &
206     &  tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
207     &  precip_rain, precip_snow, lwdown, swnet, swdown, &
208     &  tsurf, p1lay, cal, beta, coef1lay, ps, radsol, dif_grnd, &
209     &  evap, fluxsens, fluxlat, &             
210     &  tsol_rad, tsurf_new, alb_new, emis_new, z0_new, dflux_l, dflux_s)
211     
212!
213  else if (nisurf == is_oce) then
214
215    if (check) write(*,*)'ocean, nisurf = ',nisurf
216!
217! Surface "ocean" appel a l'interface avec l'ocean
218!
219!    if (ocean == 'couple') then
220!      call interfoce(nisurf, ocean)
221!    else if (ocean == 'slab  ') then
222!      call interfoce(nisurf)
223!    else                              ! lecture conditions limites
224!      call interfoce(itime, dtime, jour, &
225!     &  klon, nisurf, knon, knindex, &
226!     &  debut, &
227!     &  tsurf_new, alb_new, z0_new, pctsrf_new)
228!    endif
229!
230    call calcul_fluxs( knon, dtime, &
231     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
232     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
233     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
234     &   tsurf_new, fluxlat, fluxsens, dflux_s, dflux_l)
235!
236  else if (nisurf == is_sic) then
237
238    if (check) write(*,*)'sea ice, nisurf = ',nisurf
239
240!
241! Surface "glace de mer" appel a l'interface avec l'ocean
242!
243!    call interfoce(nisurf, ocean)
244!
245    call calcul_fluxs( knon, dtime, &
246     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
247     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
248     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
249     &   tsurf_new, fluxlat, fluxsens, dflux_s, dflux_l)
250
251  else if (nisurf == is_lic) then
252
253    if (check) write(*,*)'glacier, nisurf = ',nisurf
254
255!
256! Surface "glacier continentaux" appel a l'interface avec le sol
257!
258!    call interfsol(nisurf)
259    call calcul_fluxs( knon, dtime, &
260     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
261     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
262     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
263     &   tsurf_new, fluxlat, fluxsens, dflux_s, dflux_l)
264
265  else
266    write(*,*)'Index surface = ',nisurf
267    abort_message = 'Index surface non valable'
268    call abort_gcm(modname,abort_message,1)
269  endif
270
271  END SUBROUTINE interfsurf_hq
272
273!
274!#########################################################################
275!
276  SUBROUTINE interfsurf_vent(nisurf, knon   &         
277  &                     )
278!
279! Cette routine sert d'aiguillage entre l'atmosphere et la surface en general
280! (sols continentaux, oceans, glaces) pour les tensions de vents.
281! En pratique l'interface se fait entre la couche limite du modele
282! atmospherique (clmain.F) et les routines de surface (sechiba, oasis, ...)
283!
284!
285! L.Fairhead 02/2000
286!
287! input:
288!   nisurf       index de la surface a traiter (1 = sol continental)
289!   knon         nombre de points de la surface a traiter
290
291! Parametres d'entree
292  integer, intent(IN) :: nisurf
293  integer, intent(IN) :: knon
294
295
296  return
297  END SUBROUTINE interfsurf_vent
298!
299!#########################################################################
300!
301  SUBROUTINE interfsol(itime, klon, dtime, nisurf, knon, &
302     & knindex, rlon, rlat, &
303     & debut, lafin, ok_veget, &
304     & zlev, zlflu, u1_lay, v1_lay, temp_air, spechum, hum_air, ccanopy, &
305     & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
306     & precip_rain, precip_snow, lwdown, swnet, swdown, &
307     & tsurf, p1lay, cal, beta, coef1lay, ps, radsol, dif_grnd, &
308     & evap, fluxsens, fluxlat, &             
309     & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, dflux_l, dflux_s)
310
311! Cette routine sert d'interface entre le modele atmospherique et le
312! modele de sol continental. Appel a sechiba
313!
314! L. Fairhead 02/2000
315!
316! input:
317!   itime        numero du pas de temps
318!   klon         nombre total de points de grille
319!   dtime        pas de temps de la physique (en s)
320!   nisurf       index de la surface a traiter (1 = sol continental)
321!   knon         nombre de points de la surface a traiter
322!   knindex      index des points de la surface a traiter
323!   rlon         longitudes de la grille entiere
324!   rlat         latitudes de la grille entiere
325!   debut        logical: 1er appel a la physique (lire les restart)
326!   lafin        logical: dernier appel a la physique (ecrire les restart)
327!   ok_veget     logical: appel ou non au schema de surface continental
328!                     (si false calcul simplifie des fluxs sur les continents)
329!   zlev         hauteur de la premiere couche
330!   zlflu       
331!   u1_lay       vitesse u 1ere couche
332!   v1_lay       vitesse v 1ere couche
333!   temp_air     temperature de l'air 1ere couche
334!   spechum      humidite specifique 1ere couche
335!   hum_air      humidite de l'air
336!   ccanopy      concentration CO2 canopee
337!   tq_cdrag     cdrag
338!   petAcoef     coeff. A de la resolution de la CL pour t
339!   peqAcoef     coeff. A de la resolution de la CL pour q
340!   petBcoef     coeff. B de la resolution de la CL pour t
341!   peqBcoef     coeff. B de la resolution de la CL pour q
342!   precip_rain  precipitation liquide
343!   precip_snow  precipitation solide
344!   lwdown       flux IR entrant a la surface
345!   swnet        flux solaire net
346!   swdown       flux solaire entrant a la surface
347!   tsurf        temperature de surface
348!   p1lay        pression 1er niveau (milieu de couche)
349!   cal          capacite calorifique du sol
350!   beta         evap reelle
351!   coef1lay     coefficient d'echange
352!   ps           pression au sol
353!   radsol       rayonnement net aus sol (LW + SW)
354!   dif_grnd     coeff. diffusion vers le sol profond     
355!   
356!
357! input/output
358!   run_off      ruissellement total
359!
360! output:
361!   evap         evaporation totale
362!   fluxsens     flux de chaleur sensible
363!   fluxlat      flux de chaleur latente
364!   tsol_rad     
365!   tsurf_new    temperature au sol
366!   alb_new      albedo
367!   emis_new     emissivite
368!   z0_new       surface roughness
369
370! Parametres d'entree
371  integer, intent(IN) :: itime
372  integer, intent(IN) :: klon
373  real, intent(IN)    :: dtime
374  integer, intent(IN) :: nisurf
375  integer, intent(IN) :: knon
376  integer, dimension(knon), intent(IN) :: knindex
377  logical, intent(IN) :: debut, lafin, ok_veget
378  real, dimension(klon), intent(IN) :: rlon, rlat
379  real, dimension(knon), intent(IN) :: zlev, zlflu
380  real, dimension(knon), intent(IN) :: u1_lay, v1_lay
381  real, dimension(knon), intent(IN) :: temp_air, spechum
382  real, dimension(knon), intent(IN) :: hum_air, ccanopy
383  real, dimension(knon), intent(IN) :: tq_cdrag, petAcoef, peqAcoef
384  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
385  real, dimension(knon), intent(IN) :: precip_rain, precip_snow
386  real, dimension(knon), intent(IN) :: lwdown, swnet, swdown, ps
387  real, dimension(knon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
388  real, dimension(knon), intent(IN) :: radsol, dif_grnd
389! Parametres de sortie
390  real, dimension(knon), intent(OUT):: evap, fluxsens, fluxlat
391  real, dimension(knon), intent(OUT):: tsol_rad, tsurf_new, alb_new
392  real, dimension(knon), intent(OUT):: emis_new, z0_new
393  real, dimension(knon), intent(OUT):: dflux_s, dflux_l
394
395! Local
396!
397  integer              :: ii
398  integer              :: error
399  character (len = 20) :: modname = 'interfsol'
400  character (len = 80) :: abort_message
401  logical              :: check = .true.
402! type de couplage dans sechiba
403!  character (len=10)   :: coupling = 'implicit'
404! drapeaux controlant les appels dans SECHIBA
405!  type(control_type), save   :: control_in
406! coordonnees geographiques
407  real, allocatable, dimension(:,:), save :: lalo
408! pts voisins
409  integer,allocatable, dimension(:,:), save :: neighbours
410! resolution de la grille
411  real, allocatable, dimension (:,:), save :: resolution
412! Identifieurs des fichiers restart et histoire
413  integer, save          :: rest_id, hist_id
414  integer, save          :: rest_id_stom, hist_id_stom
415
416  if (check) write(*,*)'Entree ', modname
417  if (check) write(*,*)'ok_veget = ',ok_veget
418
419! initialisation
420  if (.not. ok_veget) then
421    call calcul_fluxs( knon, dtime, &
422     &   tsurf, p1lay, cal, beta, coef1lay, ps, &
423     &   radsol, dif_grnd, temp_air, spechum, u1_lay, v1_lay, &
424     &   petAcoef, peqAcoef, petBcoef, peqBcoef, &
425     &   tsurf_new, fluxlat, fluxsens, dflux_s, dflux_l)
426  else
427!    if (debut) then
428!      !
429!      ! Configuration de parametres specifiques a la SSL
430!      !
431!      call intsurf_config(control_in)
432!      !
433!      ! Allouer et initialiser le tableau de coordonnees du sol
434!      !
435!      if (( .not. allocated(lalo))) then
436!        allocate(lalo(knon,2), stat = error)
437!        if (error /= 0) then
438!          abort_message='Pb allocation lalo'
439!          call abort_gcm(modname,abort_message,1)
440!        endif     
441!      endif
442!      do ii = 1, knon
443!        lalo(ii,1) = rlat(knindex(ii))
444!        lalo(ii,2) = rlon(knindex(ii))
445!      enddo
446      !-
447      !- Compute variable to help describe the grid
448      !- once the points are gathered.
449      !-
450!      IF ( (.NOT.ALLOCATED(neighbours))) THEN
451!        ALLOCATE(neighbours(knon,4), stat = error)
452!        if (error /= 0) then
453!          abort_message='Pb allocation neighbours'
454!          call abort_gcm(modname,abort_message,1)
455!        endif
456!      ENDIF
457!      IF ( (.NOT.ALLOCATED(resolution))) THEN
458!        ALLOCATE(resolution(knon,2), stat = error)
459!        if (error /= 0) then
460!          abort_message='Pb allocation resolution'
461!          call abort_gcm(modname,abort_message,1)
462!        endif
463!      ENDIF
464
465! call grid_stuff
466! call sechiba_restart_init
467! call sechiba_history_init
468
469!    endif                          ! (fin debut)
470
471!
472! Appel a la routine sols continentaux
473!
474
475!    call sechiba_main(itime, klon, knon, knindex, dtime, &
476!     & debut, lafin, coupling, control_in, &
477!     & lalo, neighbours, resolution,&
478!     & zlev, zlflu, u1_lay, v1_lay, spechum, temp_air,hum_air , ccanopy, &
479!     & tq_cdrag, petAcoef, peqAcoef, petBcoef, peqBcoef, &
480!     & precip_rain, precip_snow, lwdown, swnet, swdown, ps, &
481!     & evap, fluxsens, fluxlat, &
482!     & tsol_rad, tsurf_new, alb_new, emis_new, z0_new, &
483!     & rest_id, hist_id, rest_id_stom, hist_id_stom)
484
485!
486! Sauvegarde dans fichiers histoire
487!
488  endif  ! fin ok_veget
489
490  END SUBROUTINE interfsol
491!
492!#########################################################################
493!
494  SUBROUTINE interfoce_cpl(nisurf, ocean)
495
496! Cette routine sert d'interface entre le modele atmospherique et un
497! coupleur avec un modele d'ocean 'complet'
498!
499! L. Fairhead 02/2000
500!
501! input:
502!   nisurf       index de la surface a traiter (1 = sol continental)
503!   ocean        type d'ocean
504!
505! output:
506!
507
508! Parametres d'entree
509  integer, intent(IN) :: nisurf
510  character (len = 6)  :: ocean
511
512! Parametres de sortie
513
514! Variables locales
515
516
517! Initialisation
518! fichier restart et fichiers histoires
519
520! calcul des fluxs a passer
521
522  END SUBROUTINE interfoce_cpl
523!
524!#########################################################################
525!
526  SUBROUTINE interfoce_slab(nisurf)
527
528! Cette routine sert d'interface entre le modele atmospherique et un
529! modele de 'slab' ocean
530!
531! L. Fairhead 02/2000
532!
533! input:
534!   nisurf       index de la surface a traiter (1 = sol continental)
535!
536! output:
537!
538
539! Parametres d'entree
540  integer, intent(IN) :: nisurf
541
542  END SUBROUTINE interfoce_slab
543!
544!#########################################################################
545!
546  SUBROUTINE interfoce_lim(itime, dtime, jour, &
547     & klon, nisurf, knon, knindex, &
548     & debut,  &
549     & lmt_sst, lmt_alb, lmt_rug, pctsrf_new)
550
551! Cette routine sert d'interface entre le modele atmospherique et un fichier
552! de conditions aux limites
553!
554! L. Fairhead 02/2000
555!
556! input:
557!   itime        numero du pas de temps courant
558!   dtime        pas de temps de la physique (en s)
559!   jour         jour a lire dans l'annee
560!   nisurf       index de la surface a traiter (1 = sol continental)
561!   knon         nombre de points dans le domaine a traiter
562!   knindex      index des points de la surface a traiter
563!   klon         taille de la grille
564!   debut        logical: 1er appel a la physique (initialisation)
565!
566! output:
567!   lmt_sst      SST lues dans le fichier de CL
568!   lmt_alb      Albedo lu
569!   lmt_rug      longueur de rugosité lue
570!   pctsrf_new   sous-maille fractionnelle
571!
572
573#include "indicesol.h"
574
575! Parametres d'entree
576  integer, intent(IN) :: itime
577  real   , intent(IN) :: dtime
578  integer, intent(IN) :: jour
579  integer, intent(IN) :: nisurf
580  integer, intent(IN) :: knon
581  integer, intent(IN) :: klon
582  integer, dimension(knon), intent(in) :: knindex
583  logical, intent(IN) :: debut
584
585! Parametres de sortie
586  real, intent(out), dimension(knon) :: lmt_sst
587  real, intent(out), dimension(knon) :: lmt_alb
588  real, intent(out), dimension(knon) :: lmt_rug
589  real, intent(out), dimension(klon,nbsrf) :: pctsrf_new
590
591! Variables locales
592  integer     :: ii
593  integer     :: lmt_pas     ! frequence de lecture des conditions limites
594                             ! (en pas de physique)
595  logical,save :: deja_lu    ! pour indiquer que le jour a lire a deja
596                             ! lu pour une surface precedente
597  integer,save :: jour_lu
598  integer      :: ierr
599  character (len = 20) :: modname = 'interfoce_lim'
600  character (len = 80) :: abort_message
601  character (len = 20) :: fich ='limit'
602  logical     :: newlmt = .false.
603  logical     :: check = .true.
604! Champs lus dans le fichier de CL
605  real, allocatable , save, dimension(:) :: sst_lu, alb_lu, rug_lu, nat_lu
606  real, allocatable , save, dimension(:,:) :: pct_tmp
607!
608! quelques variables pour netcdf
609!
610#include "netcdf.inc"
611  integer              :: nid, nvarid
612  integer, dimension(2) :: start, epais
613!
614! Fin déclaration
615!
616   
617  if (debut) then
618    lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
619    jour_lu = jour - 1
620    allocate(sst_lu(klon))
621    allocate(alb_lu(klon))
622    allocate(rug_lu(klon))
623    allocate(nat_lu(klon))
624    allocate(pct_tmp(klon,nbsrf))
625  endif
626
627  if ((jour - jour_lu) /= 0) deja_lu = .false.
628 
629  if (check) write(*,*)modname,' :: jour_lu, deja_lu', jour_lu, deja_lu
630
631! Tester d'abord si c'est le moment de lire le fichier
632  if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu) then
633!
634! Ouverture du fichier
635!
636    ierr = NF_OPEN (fich, NF_NOWRITE,nid)
637    if (ierr.NE.NF_NOERR) then
638      abort_message = 'Pb d''ouverture du fichier de conditions aux limites'
639      call abort_gcm(modname,abort_message,1)
640    endif
641!
642! La tranche de donnees a lire:
643!
644    start(1) = 1
645    start(2) = jour + 1
646    epais(1) = klon
647    epais(2) = 1
648!
649    if (newlmt) then
650!
651! Fraction "ocean"
652!
653      ierr = NF_INQ_VARID(nid, 'FOCE', nvarid)
654      if (ierr /= NF_NOERR) then
655        abort_message = 'Le champ <FOCE> est absent'
656        call abort_gcm(modname,abort_message,1)
657      endif
658#ifdef NC_DOUBLE
659      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_oce))
660#else
661      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_oce))
662#endif
663      if (ierr /= NF_NOERR) then
664        abort_message = 'Lecture echouee pour <FOCE>'
665        call abort_gcm(modname,abort_message,1)
666      endif
667!
668! Fraction "glace de mer"
669!
670      ierr = NF_INQ_VARID(nid, 'FSIC', nvarid)
671      if (ierr /= NF_NOERR) then
672        abort_message = 'Le champ <FSIC> est absent'
673        call abort_gcm(modname,abort_message,1)
674      endif
675#ifdef NC_DOUBLE
676      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_sic))
677#else
678      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_sic))
679#endif
680      if (ierr /= NF_NOERR) then
681        abort_message = 'Lecture echouee pour <FSIC>'
682        call abort_gcm(modname,abort_message,1)
683      endif
684!
685! Fraction "terre"
686!
687      ierr = NF_INQ_VARID(nid, 'FTER', nvarid)
688      if (ierr /= NF_NOERR) then
689        abort_message = 'Le champ <FTER> est absent'
690        call abort_gcm(modname,abort_message,1)
691      endif
692#ifdef NC_DOUBLE
693      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_ter))
694#else
695      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_ter))
696#endif
697      if (ierr /= NF_NOERR) then
698        abort_message = 'Lecture echouee pour <FTER>'
699        call abort_gcm(modname,abort_message,1)
700      endif
701!
702! Fraction "glacier terre"
703!
704      ierr = NF_INQ_VARID(nid, 'FLIC', nvarid)
705      if (ierr /= NF_NOERR) then
706        abort_message = 'Le champ <FLIC> est absent'
707        call abort_gcm(modname,abort_message,1)
708      endif
709#ifdef NC_DOUBLE
710      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pctsrf_new(1,is_lic))
711#else
712      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pctsrf_new(1,is_lic))
713#endif
714      if (ierr /= NF_NOERR) then
715        abort_message = 'Lecture echouee pour <FLIC>'
716        call abort_gcm(modname,abort_message,1)
717      endif
718!
719    else  ! on en est toujours a rnatur
720!
721      ierr = NF_INQ_VARID(nid, 'NAT', nvarid)
722      if (ierr /= NF_NOERR) then
723        abort_message = 'Le champ <NAT> est absent'
724        call abort_gcm(modname,abort_message,1)
725      endif
726#ifdef NC_DOUBLE
727      ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, nat_lu)
728#else
729      ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, nat_lu)
730#endif
731      if (ierr /= NF_NOERR) then
732        abort_message = 'Lecture echouee pour <NAT>'
733        call abort_gcm(modname,abort_message,1)
734      endif
735!
736! Remplissage des fractions de surface
737! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice
738!
739      pct_tmp = 0.0
740      do ii = 1, klon
741        pct_tmp(ii,nint(nat_lu(ii)) + 1) = 1.
742      enddo
743
744!
745!  On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire
746!
747      pctsrf_new = pct_tmp
748      pctsrf_new (:,2)= pct_tmp (:,1)
749      pctsrf_new (:,1)= pct_tmp (:,2)
750      pct_tmp = pctsrf_new
751    endif ! fin test sur newlmt
752!
753! Lecture SST
754!
755    ierr = NF_INQ_VARID(nid, 'SST', nvarid)
756    if (ierr /= NF_NOERR) then
757      abort_message = 'Le champ <SST> est absent'
758      call abort_gcm(modname,abort_message,1)
759    endif
760#ifdef NC_DOUBLE
761    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, sst_lu)
762#else
763    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, sst_lu)
764#endif
765    if (ierr /= NF_NOERR) then
766      abort_message = 'Lecture echouee pour <SST>'
767      call abort_gcm(modname,abort_message,1)
768    endif   
769!
770! Lecture Albedo
771!
772    ierr = NF_INQ_VARID(nid, 'ALB', nvarid)
773    if (ierr /= NF_NOERR) then
774      abort_message = 'Le champ <ALB> est absent'
775      call abort_gcm(modname,abort_message,1)
776    endif
777#ifdef NC_DOUBLE
778    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, alb_lu)
779#else
780    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, alb_lu)
781#endif
782    if (ierr /= NF_NOERR) then
783      abort_message = 'Lecture echouee pour <ALB>'
784      call abort_gcm(modname,abort_message,1)
785    endif
786!
787! Lecture rugosité
788!
789    ierr = NF_INQ_VARID(nid, 'RUG', nvarid)
790    if (ierr /= NF_NOERR) then
791      abort_message = 'Le champ <RUG> est absent'
792      call abort_gcm(modname,abort_message,1)
793    endif
794#ifdef NC_DOUBLE
795    ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, rug_lu)
796#else
797    ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, rug_lu)
798#endif
799    if (ierr /= NF_NOERR) then
800      abort_message = 'Lecture echouee pour <RUG>'
801      call abort_gcm(modname,abort_message,1)
802    endif
803
804!
805! Fin de lecture
806!
807    ierr = NF_CLOSE(nid)
808    deja_lu = .true.
809    jour_lu = jour
810  endif
811!
812! Recopie des variables dans les champs de sortie
813!
814  do ii = 1, knon
815    lmt_sst(ii) = sst_lu(knindex(ii))
816    lmt_alb(ii) = alb_lu(knindex(ii))
817    lmt_rug(ii) = rug_lu(knindex(ii))
818  enddo
819  pctsrf_new = pct_tmp
820
821  END SUBROUTINE interfoce_lim
822
823!
824!#########################################################################
825!
826
827  SUBROUTINE calcul_fluxs( knon, dtime, &
828     & tsurf, p1lay, cal, beta, coef1lay, ps, &
829     & radsol, dif_grnd, t1lay, q1lay, u1lay, v1lay, &
830     & petAcoef, peqAcoef, petBcoef, peqBcoef, &
831     & tsurf_new, fluxlat, fluxsens, dflux_s, dflux_l)
832
833! Cette routine calcule les fluxs en h et q a l'interface et eventuellement
834! une temperature de surface (au cas ou ok_veget = false)
835!
836! L. Fairhead 4/2000
837!
838! input:
839!   knon         nombre de points a traiter
840!   tsurf        temperature de surface
841!   p1lay        pression 1er niveau (milieu de couche)
842!   cal          capacite calorifique du sol
843!   beta         evap reelle
844!   coef1lay     coefficient d'echange
845!   ps           pression au sol
846!   petAcoef     coeff. A de la resolution de la CL pour t
847!   peqAcoef     coeff. A de la resolution de la CL pour q
848!   petBcoef     coeff. B de la resolution de la CL pour t
849!   peqBcoef     coeff. B de la resolution de la CL pour q
850!   radsol       rayonnement net aus sol (LW + SW)
851!   dif_grnd     coeff. diffusion vers le sol profond
852!
853! output:
854!   tsurf_new    temperature au sol
855!   fluxsens     flux de chaleur sensible
856!   fluxlat      flux de chaleur latente
857!   dflux_s      derivee du flux de chaleur sensible / Ts
858!   dflux_l      derivee du flux de chaleur latente  / Ts
859!
860
861#include "YOMCST.inc"
862#include "YOETHF.inc"
863#include "FCTTRE.inc"
864
865! Parametres d'entree
866  integer, intent(IN) :: knon
867  real   , intent(IN) :: dtime
868  real, dimension(knon), intent(IN) :: petAcoef, peqAcoef
869  real, dimension(knon), intent(IN) :: petBcoef, peqBcoef
870  real, dimension(knon), intent(IN) :: ps, q1lay
871  real, dimension(knon), intent(IN) :: tsurf, p1lay, cal, beta, coef1lay
872  real, dimension(knon), intent(IN) :: radsol, dif_grnd
873  real, dimension(knon), intent(IN) :: t1lay, u1lay, v1lay
874
875! Parametres sorties
876  real, dimension(knon), intent(OUT):: tsurf_new, fluxlat, fluxsens
877  real, dimension(knon), intent(OUT):: dflux_s, dflux_l
878
879! Variables locales
880  integer :: i
881  real, dimension(knon) :: zx_mh, zx_nh, zx_oh
882  real, dimension(knon) :: zx_mq, zx_nq, zx_oq
883  real, dimension(knon) :: zx_pkh, zx_dq_s_dt, zx_qsat, zx_coef
884  real, dimension(knon) :: zx_sl, zx_k1,  zx_dq,  zx_cq,  zx_dh, zx_ch
885  real, dimension(knon) :: zx_h_ts, zx_q_0 , d_ts
886  real                  :: zdelta, zcvm5, zx_qs, zcor, zx_dq_s_dh
887  real, parameter :: t_grnd = 271.35, t_coup = 273.15
888  logical         :: check = .true.
889  character (len = 20)  :: modname = 'calcul_fluxs'
890
891  if (check) write(*,*)'Entree ', modname
892!
893! Initialisation
894!
895!
896! zx_qs = qsat en kg/kg
897!
898  DO i = 1, knon
899    zx_pkh(i) = (ps(i)/ps(i))**RKAPPA
900  IF (thermcep) THEN
901      zdelta=MAX(0.,SIGN(1.,rtt-tsurf(i)))
902      zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta
903      zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q1lay(i))
904      zx_qs= r2es * FOEEW(tsurf(i),zdelta)/ps(i)
905      zx_qs=MIN(0.5,zx_qs)
906      zcor=1./(1.-retv*zx_qs)
907      zx_qs=zx_qs*zcor
908      zx_dq_s_dh = FOEDE(tsurf(i),zdelta,zcvm5,zx_qs,zcor) &
909     &                 /RLVTT / zx_pkh(i)
910    ELSE
911      IF (tsurf(i).LT.t_coup) THEN
912        zx_qs = qsats(tsurf(i)) / ps(i)
913        zx_dq_s_dh = dqsats(tsurf(i),zx_qs)/RLVTT &
914     &                    / zx_pkh(i)
915      ELSE
916        zx_qs = qsatl(tsurf(i)) / ps(i)
917        zx_dq_s_dh = dqsatl(tsurf(i),zx_qs)/RLVTT &
918     &               / zx_pkh(i)
919      ENDIF
920    ENDIF
921    zx_dq_s_dt(i) = RCPD * zx_pkh(i) * zx_dq_s_dh
922    zx_qsat(i) = zx_qs
923    zx_coef(i) = coef1lay(i) &
924     & * (1.0+SQRT(u1lay(i)**2+v1lay(i)**2)) &
925     & * p1lay(i)/(RD*t1lay(i))
926   
927  ENDDO
928
929
930! === Calcul de la temperature de surface ===
931!
932! zx_sl = chaleur latente d'evaporation ou de sublimation
933!
934  do i = 1, knon
935    zx_sl(i) = RLVTT
936    if (tsurf(i) .LT. RTT) zx_sl(i) = RLSTT
937    zx_k1(i) = zx_coef(i)
938  enddo
939
940
941  do i = 1, knon
942! Q
943    zx_oq(i) = 1. - (beta(i) * zx_k1(i) * peqBcoef(i) * dtime)
944    zx_mq(i) = beta(i) * zx_k1(i) * &
945     &             (peqAcoef(i) - zx_qsat(i) &
946     &                          + zx_dq_s_dt(i) * tsurf(i)) &
947     &             / zx_oq(i)
948    zx_nq(i) = beta(i) * zx_k1(i) * (-1. * zx_dq_s_dt(i)) &
949     &                              / zx_oq(i)
950
951! H
952    zx_oh(i) = 1. - (zx_k1(i) * petBcoef(i) * dtime)
953    zx_mh(i) = zx_k1(i) * petAcoef(i) / zx_oh(i)
954    zx_nh(i) = - (zx_k1(i) * RCPD * zx_pkh(i))/ zx_oh(i)
955
956! Tsurface
957    tsurf_new(i) = (tsurf(i) + cal(i)/(RCPD * zx_pkh(i)) * dtime * &
958     &             (radsol(i) + zx_mh(i) + zx_sl(i) * zx_mq(i)) &
959     &                 + dif_grnd(i) * t_grnd * dtime)/ &
960     &          ( 1. - dtime * cal(i)/(RCPD * zx_pkh(i)) * ( &
961     &                       zx_nh(i) + zx_sl(i) * zx_nq(i)) & 
962     &                     + dtime * dif_grnd(i))
963    zx_h_ts(i) = tsurf_new(i) * RCPD * zx_pkh(i)
964    d_ts(i) = tsurf_new(i) - tsurf(i)
965    zx_q_0(i) = zx_qsat(i) + zx_dq_s_dt(i) * d_ts(i)
966!== flux_q est le flux de vapeur d'eau: kg/(m**2 s)  positive vers bas
967!== flux_t est le flux de cpt (energie sensible): j/(m**2 s)
968    fluxlat(i) = zx_mq(i) + zx_nq(i) * tsurf_new(i)
969    fluxsens(i) = zx_mh(i) + zx_nh(i) * tsurf_new(i)
970! Derives des flux dF/dTs (W m-2 K-1):
971    dflux_s(i) = zx_nh(i)
972    dflux_l(i) = (zx_sl(i) * zx_nq(i))
973  ENDDO
974
975  END SUBROUTINE calcul_fluxs
976!
977!#########################################################################
978!
979
980  END MODULE interface_surf
Note: See TracBrowser for help on using the repository browser.