/[lmdze]/trunk/Sources/phylmd/Interface_surf/interfoce_lim.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/Interface_surf/interfoce_lim.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 82 by guez, Wed Mar 5 14:57:53 2014 UTC revision 104 by guez, Thu Sep 4 10:05:52 2014 UTC
# Line 4  module interfoce_lim_m Line 4  module interfoce_lim_m
4    
5  contains  contains
6    
7    SUBROUTINE interfoce_lim(itime, dtime, jour,  &    SUBROUTINE interfoce_lim(itime, dtime, jour, klon, knon, knindex, &
8         klon, nisurf, knon, knindex,  &         debut, lmt_sst, pctsrf_new)
        debut,  &  
        lmt_sst, pctsrf_new)  
9    
10      ! Cette routine sert d'interface entre le modele atmospherique et      ! Cette routine sert d'interface entre le modèle atmosphérique et
11      ! un fichier de conditions aux limites      ! un fichier de conditions aux limites.
12    
13      ! L. Fairhead 02/2000      ! Laurent FAIRHEAD, February 2000
14    
15      use abort_gcm_m, only: abort_gcm      use abort_gcm_m, only: abort_gcm
16      use indicesol      USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
17      use netcdf      USE netcdf, ONLY: nf90_noerr, nf90_nowrite
18        use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open
19    
20      integer, intent(IN) :: itime ! numero du pas de temps courant      integer, intent(IN):: itime ! numero du pas de temps courant
21      real , intent(IN) :: dtime ! pas de temps de la physique (en s)      real, intent(IN):: dtime ! pas de temps de la physique (en s)
22      integer, intent(IN) :: jour ! jour a lire dans l'annee      integer, intent(IN):: jour ! jour a lire dans l'annee
23      integer, intent(IN) :: nisurf ! index de la surface a traiter (1 = sol continental)      integer, intent(IN):: klon ! taille de la grille
24      integer, intent(IN) :: knon ! nombre de points dans le domaine a traiter      integer, intent(IN):: knon ! nombre de points dans le domaine a traiter
25      integer, intent(IN) :: klon ! taille de la grille  
26      integer, dimension(klon), intent(in) :: knindex ! index des points de la surface a traiter      integer, intent(in):: knindex(klon)
27      logical, intent(IN) :: debut ! logical: 1er appel a la physique (initialisation)      ! index des points de la surface a traiter
28    
29      ! Parametres de sortie      logical, intent(IN):: debut ! 1er appel a la physique (initialisation)
30      ! output:  
31      ! lmt_sst SST lues dans le fichier de CL      real, intent(out), dimension(knon):: lmt_sst
32      ! pctsrf_new sous-maille fractionnelle      ! SST lues dans le fichier de CL
33      real, intent(out), dimension(klon) :: lmt_sst  
34      real, intent(out), dimension(klon, nbsrf) :: pctsrf_new      real, intent(out), dimension(klon, nbsrf):: pctsrf_new
35        ! sous-maille fractionnelle
36      ! Variables locales  
37      integer :: ii      ! Local:
38      INTEGER, save :: lmt_pas ! frequence de lecture des conditions limites  
39        integer ii
40    
41        INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites
42      ! (en pas de physique)      ! (en pas de physique)
43      logical, save :: deja_lu ! pour indiquer que le jour a lire a deja  
44        logical, save:: deja_lu ! pour indiquer que le jour a lire a deja
45      ! lu pour une surface precedente      ! lu pour une surface precedente
     integer, save :: jour_lu  
     integer :: ierr  
     character (len = 20) :: modname = 'interfoce_lim'  
     character (len = 80) :: abort_message  
     logical, save :: newlmt = .TRUE.  
     logical, save :: check = .FALSE.  
     ! Champs lus dans le fichier de CL  
     real, allocatable , save, dimension(:) :: sst_lu, rug_lu, nat_lu  
     real, allocatable , save, dimension(:, :) :: pct_tmp  
46    
47      ! quelques variables pour netcdf      integer, save:: jour_lu
48    
49      integer :: nid, nvarid      ! Champs lus dans le fichier de conditions aux limites :
50      integer, dimension(2) :: start, epais      real, allocatable, save, dimension(:):: sst_lu
51        real, allocatable, save, dimension(:, :):: pct_tmp
52    
53        ! Pour Netcdf :
54        integer nid, nvarid
55        integer, dimension(2):: start, epais
56    
57      ! --------------------------------------------------      ! --------------------------------------------------
58    
# Line 61  contains Line 60  contains
60         lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour         lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
61         jour_lu = jour - 1         jour_lu = jour - 1
62         allocate(sst_lu(klon))         allocate(sst_lu(klon))
        allocate(nat_lu(klon))  
63         allocate(pct_tmp(klon, nbsrf))         allocate(pct_tmp(klon, nbsrf))
64      endif      endif
65    
66      if ((jour - jour_lu) /= 0) deja_lu = .false.      if ((jour - jour_lu) /= 0) deja_lu = .false.
67    
     if (check) write(*, *)modname, ' :: jour, jour_lu, deja_lu', jour, jour_lu, &  
          deja_lu  
     if (check) write(*, *)modname, ' :: itime, lmt_pas ', itime, lmt_pas, dtime  
   
68      ! Tester d'abord si c'est le moment de lire le fichier      ! Tester d'abord si c'est le moment de lire le fichier
69      if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu) then      if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu) then
70           call NF95_OPEN ('limit.nc', NF90_NOWRITE, nid)
        ! Ouverture du fichier  
   
        ierr = NF90_OPEN ('limit.nc', NF90_NOWRITE, nid)  
        if (ierr.NE.NF90_NOERR) then  
           abort_message &  
                = 'Pb d''ouverture du fichier de conditions aux limites'  
           call abort_gcm(modname, abort_message, 1)  
        endif  
71    
72         ! La tranche de donnees a lire:         ! La tranche de donnees a lire:
   
73         start(1) = 1         start(1) = 1
74         start(2) = jour         start(2) = jour
75         epais(1) = klon         epais(1) = klon
76         epais(2) = 1         epais(2) = 1
77    
78         if (newlmt) then         ! Fraction "ocean"
79           call NF95_INQ_VARID(nid, 'FOCE', nvarid)
80            ! Fraction "ocean"         call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_oce), start, epais)
81    
82            ierr = NF90_INQ_VARID(nid, 'FOCE', nvarid)         ! Fraction "glace de mer"
83            if (ierr /= NF90_NOERR) then         call NF95_INQ_VARID(nid, 'FSIC', nvarid)
84               abort_message = 'Le champ <FOCE> est absent'         call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_sic), start, epais)
85               call abort_gcm(modname, abort_message, 1)  
86            endif         ! Fraction "terre"
87            ierr = NF90_GET_VAR(nid, nvarid, pct_tmp(:, is_oce), start, epais)         call NF95_INQ_VARID(nid, 'FTER', nvarid)
88            if (ierr /= NF90_NOERR) then         call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_ter), start, epais)
89               abort_message = 'Lecture echouee pour <FOCE>'  
90               call abort_gcm(modname, abort_message, 1)         ! Fraction "glacier terre"
91            endif         call NF95_INQ_VARID(nid, 'FLIC', nvarid)
92           call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_lic), start, epais)
           ! Fraction "glace de mer"  
   
           ierr = NF90_INQ_VARID(nid, 'FSIC', nvarid)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Le champ <FSIC> est absent'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
           ierr = NF90_GET_VAR(nid, nvarid, pct_tmp(:, is_sic), start, epais)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Lecture echouee pour <FSIC>'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
   
           ! Fraction "terre"  
   
           ierr = NF90_INQ_VARID(nid, 'FTER', nvarid)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Le champ <FTER> est absent'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
           ierr = NF90_GET_VAR(nid, nvarid, pct_tmp(:, is_ter), start, epais)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Lecture echouee pour <FTER>'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
   
           ! Fraction "glacier terre"  
   
           ierr = NF90_INQ_VARID(nid, 'FLIC', nvarid)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Le champ <FLIC> est absent'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
           ierr = NF90_GET_VAR(nid, nvarid, pct_tmp(:, is_lic), start, epais)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Lecture echouee pour <FLIC>'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
   
        else ! on en est toujours a rnatur  
   
           ierr = NF90_INQ_VARID(nid, 'NAT', nvarid)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Le champ <NAT> est absent'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
           ierr = NF90_GET_VAR(nid, nvarid, nat_lu, start, epais)  
           if (ierr /= NF90_NOERR) then  
              abort_message = 'Lecture echouee pour <NAT>'  
              call abort_gcm(modname, abort_message, 1)  
           endif  
   
           ! Remplissage des fractions de surface  
           ! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice  
   
           pct_tmp = 0.0  
           do ii = 1, klon  
              pct_tmp(ii, nint(nat_lu(ii)) + 1) = 1.  
           enddo  
   
   
           ! On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire  
   
           pctsrf_new = pct_tmp  
           pctsrf_new (:, 2)= pct_tmp (:, 1)  
           pctsrf_new (:, 1)= pct_tmp (:, 2)  
           pct_tmp = pctsrf_new  
        endif ! fin test sur newlmt  
   
        ! Lecture SST  
   
        ierr = NF90_INQ_VARID(nid, 'SST', nvarid)  
        if (ierr /= NF90_NOERR) then  
           abort_message = 'Le champ <SST> est absent'  
           call abort_gcm(modname, abort_message, 1)  
        endif  
        ierr = NF90_GET_VAR(nid, nvarid, sst_lu, start, epais)  
        if (ierr /= NF90_NOERR) then  
           abort_message = 'Lecture echouee pour <SST>'  
           call abort_gcm(modname, abort_message, 1)  
        endif  
   
93    
94         ! Fin de lecture         call NF95_INQ_VARID(nid, 'SST', nvarid)
95           call NF95_GET_VAR(nid, nvarid, sst_lu, start, epais)
96    
97         ierr = NF90_CLOSE(nid)         call NF95_CLOSE(nid)
98         deja_lu = .true.         deja_lu = .true.
99         jour_lu = jour         jour_lu = jour
100      endif      endif
101    
102      ! Recopie des variables dans les champs de sortie      ! Recopie des variables dans les champs de sortie
103    
     lmt_sst = 999999999.  
104      do ii = 1, knon      do ii = 1, knon
105         lmt_sst(ii) = sst_lu(knindex(ii))         lmt_sst(ii) = sst_lu(knindex(ii))
106      enddo      enddo

Legend:
Removed from v.82  
changed lines
  Added in v.104

  ViewVC Help
Powered by ViewVC 1.1.21