/[lmdze]/trunk/Sources/dyn3d/Guide/Read_reanalyse/read_reanalyse.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/Guide/Read_reanalyse/read_reanalyse.f

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

revision 171 by guez, Tue Sep 29 19:48:59 2015 UTC revision 172 by guez, Wed Sep 30 15:59:14 2015 UTC
# Line 4  module read_reanalyse_m Line 4  module read_reanalyse_m
4    
5  contains  contains
6    
7    subroutine read_reanalyse(timestep, psi, u, v, t, q, masse, nlevnc)    subroutine read_reanalyse(psi, u, v, t, q)
8    
9      ! From LMDZ4/libf/dyn3d/read_reanalyse.F, version 1.3, 2005/04/15 12:31:21      ! From LMDZ4/libf/dyn3d/read_reanalyse.F, version 1.3, 2005/04/15 12:31:21
10    
11      USE conf_guide_m, ONLY: guide_q, guide_t, guide_u, guide_v, ncep      USE conf_guide_m, ONLY: guide_q, guide_t, guide_u, guide_v, ncep
12      USE dimens_m, ONLY: iim, jjm, llm      USE dimens_m, ONLY: iim, jjm, llm
13      USE netcdf, ONLY: nf90_get_var, nf90_inq_varid, nf90_nowrite, nf90_open      USE netcdf, ONLY: nf90_nowrite
14        USE netcdf95, ONLY: nf95_get_var, nf95_inq_dimid, nf95_inq_varid, &
15             nf95_inquire_dimension, nf95_open
16      USE paramet_m, ONLY: iip1, jjp1      USE paramet_m, ONLY: iip1, jjp1
17      use reanalyse2nat_m, only: reanalyse2nat      use reanalyse2nat_m, only: reanalyse2nat
18    
19      integer, intent(in):: timestep      real, intent(in):: psi(:, :) ! (iip1, jjp1)
20      real, intent(in):: psi(iip1, jjp1)      real, intent(out):: u(:, :, :) ! (iip1, jjp1, llm)
21      real u(iip1, jjp1, llm), v(iip1, jjm, llm)      real, intent(out):: v(:, :, :) ! (iip1, jjm, llm)
22      real t(iip1, jjp1, llm), q(iip1, jjp1, llm)      real, intent(out):: t(:, :, :), q(:, :, :) ! (iip1, jjp1, llm)
     real masse(iip1, jjp1, llm)  
     integer nlevnc  
23    
24      ! Local:      ! Local:
25        integer, save:: nlevnc
26      integer l      integer:: timestep = 0
27      real pk(iip1, jjp1, llm)      real pk(iip1, jjp1, llm)
28      integer, save:: ncidu, varidu, ncidv, varidv, ncidt, varidt      integer, save:: ncidu, varidu, ncidv, varidv, ncidt, varidt, ncidQ, varidQ
29      integer, save:: ncidpl      integer ncid, varid, dimid
30      integer, save:: varidpl, ncidQ, varidQ      real, allocatable, save:: unc(:, :, :) ! (iip1, jjp1, nlevnc)
31      real unc(iip1, jjp1, nlevnc), vnc(iip1, jjm, nlevnc)      real, allocatable, save:: vnc(:, :, :) ! (iip1, jjm, nlevnc)
32      real tnc(iip1, jjp1, nlevnc)      real, allocatable, save:: tnc(:, :, :), Qnc(:, :, :) ! (iip1, jjp1, nlevnc)
33      real Qnc(iip1, jjp1, nlevnc)      real, allocatable, save:: pl(:) ! (nlevnc)
     real pl(nlevnc)  
     integer start(4), count(4), status  
     real rcode  
34      logical:: first = .true.      logical:: first = .true.
35        character(len = 8) name
36    
37      ! -----------------------------------------------------------------      ! -----------------------------------------------------------------
38    
39      !   Initialisation de la lecture des fichiers      ! Initialisation de la lecture des fichiers
40    
41      if (first) then      if (first) then
        ncidpl=-99  
42         print *, 'Intitialisation de read reanalsye'         print *, 'Intitialisation de read reanalsye'
43    
44         ! Vent zonal         ! Vent zonal
45         if (guide_u) then         if (guide_u) then
46            rcode=nf90_open('u.nc', nf90_nowrite, ncidu)            call nf95_open('u.nc', nf90_nowrite, ncidu)
47            rcode = nf90_inq_varid(ncidu, 'UWND', varidu)            call nf95_inq_varid(ncidu, 'UWND', varidu)
           if (ncidpl.eq.-99) ncidpl=ncidu  
48         endif         endif
49    
50         ! Vent meridien         ! Vent meridien
51         if (guide_v) then         if (guide_v) then
52            rcode=nf90_open('v.nc', nf90_nowrite, ncidv)            call nf95_open('v.nc', nf90_nowrite, ncidv)
53            rcode = nf90_inq_varid(ncidv, 'VWND', varidv)            call nf95_inq_varid(ncidv, 'VWND', varidv)
           if (ncidpl.eq.-99) ncidpl=ncidv  
54         endif         endif
55    
56         ! Temperature         ! Temperature
57         if (guide_T) then         if (guide_T) then
58            rcode=nf90_open('T.nc', nf90_nowrite, ncidt)            call nf95_open('T.nc', nf90_nowrite, ncidt)
59            rcode = nf90_inq_varid(ncidt, 'AIR', varidt)            call nf95_inq_varid(ncidt, 'AIR', varidt)
           if (ncidpl.eq.-99) ncidpl=ncidt  
60         endif         endif
61    
62         ! Humidite         ! Humidite
63         if (guide_Q) then         if (guide_Q) then
64            rcode=nf90_open('hur.nc', nf90_nowrite, ncidQ)            call nf95_open('hur.nc', nf90_nowrite, ncidQ)
65            rcode = nf90_inq_varid(ncidQ, 'RH', varidQ)            call nf95_inq_varid(ncidQ, 'RH', varidQ)
           if (ncidpl.eq.-99) ncidpl=ncidQ  
        endif  
   
        ! Coordonnee verticale  
        if (ncep) then  
           print *, 'Vous etes entrain de lire des donnees NCEP'  
           rcode = nf90_inq_varid(ncidpl, 'LEVEL', varidpl)  
        else  
           print *, 'Vous etes entrain de lire des donnees ECMWF'  
           rcode = nf90_inq_varid(ncidpl, 'PRESSURE', varidpl)  
66         endif         endif
     endif  
   
     ! Niveaux de pression  
   
     ! Warning: il n y a pas de test de coherence sur le nombre de  
     ! niveaux verticaux dans le fichier nc'  
     status=NF90_GET_VAR(ncidpl, varidpl, pl)  
     !  passage en pascal  
     pl(:)=100.*pl(:)  
     if (first) then  
        do l=1, nlevnc  
           print *, 'PL(', l, ')=', pl(l)  
        enddo  
     endif  
   
     !   lecture des champs u, v, T  
   
     !  dimensions pour les champs scalaires et le vent zonal  
67    
68      start(1)=1         ! Coordonn\'ee verticale :
     start(2)=1  
     start(3)=1  
     start(4)=timestep  
69    
70      count(1)=iip1         if (guide_u) then
71      count(2)=jjp1            ncid = ncidu
72      count(3)=nlevnc         else if (guide_v) then
73      count(4)=1            ncid = ncidv
74           else if (guide_T) then
75      ! mise a zero des tableaux            ncid = ncidt
76           else
77      unc(:, :, :)=0.            ncid = ncidq
78      vnc(:, :, :)=0.         end if
     tnc(:, :, :)=0.  
     Qnc(:, :, :)=0.  
79    
80      !  Vent zonal         name = merge('LEVEL   ', 'PRESSURE', ncep)
81           call nf95_inq_dimid(ncid, name, dimid)
82           call nf95_inquire_dimension(ncid, dimid, nclen = nlevnc)
83           call nf95_inq_varid(ncid, name, varid)
84           PRINT *, 'nlevnc = ', nlevnc
85           allocate(unc(iip1, jjp1, nlevnc), vnc(iip1, jjm, nlevnc))
86           allocate(tnc(iip1, jjp1, nlevnc), Qnc(iip1, jjp1, nlevnc), pl(nlevnc))
87           call NF95_GET_VAR(ncid, varid, pl)
88           pl = 100. * pl ! passage en pascal
89           first = .false.
90        endif
91    
92        ! lecture des champs u, v, T
93    
94        timestep = timestep + 1
95        unc = 0.
96        vnc = 0.
97        tnc = 0.
98        Qnc = 0.
99    
100        ! Vent zonal
101      if (guide_u) then      if (guide_u) then
102         status=NF90_GET_VAR(ncidu, varidu, unc, start, count)         call NF95_GET_VAR(ncidu, varidu, unc, start = (/1, 1, 1, timestep/))
103         ! Warning Correction bidon pour palier a un probleme dans la         call correctbid(iim, jjp1 * nlevnc, unc)
        ! creation des fichiers nc  
        call correctbid(iim, jjp1*nlevnc, unc)  
104      endif      endif
105    
106      !  Temperature      ! Temperature
   
107      if (guide_T) then      if (guide_T) then
108         status=NF90_GET_VAR(ncidt, varidt, tnc, start, count)         call NF95_GET_VAR(ncidt, varidt, tnc, start = (/1, 1, 1, timestep/))
109         call correctbid(iim, jjp1*nlevnc, tnc)         call correctbid(iim, jjp1 * nlevnc, tnc)
110      endif      endif
111    
112      !  Humidite      ! Humidite
   
113      if (guide_Q) then      if (guide_Q) then
114         status=NF90_GET_VAR(ncidQ, varidQ, Qnc, start, count)         call NF95_GET_VAR(ncidQ, varidQ, Qnc, start = (/1, 1, 1, timestep/))
115         call correctbid(iim, jjp1*nlevnc, Qnc)         call correctbid(iim, jjp1 * nlevnc, Qnc)
116      endif      endif
117    
118      count(2)=jjm      ! Vent meridien
     !  Vent meridien  
   
119      if (guide_v) then      if (guide_v) then
120         status=NF90_GET_VAR(ncidv, varidv, vnc, start, count)         call NF95_GET_VAR(ncidv, varidv, vnc, start = (/1, 1, 1, timestep/))
121         call correctbid(iim, jjm*nlevnc, vnc)         call correctbid(iim, jjm * nlevnc, vnc)
122      endif      endif
123    
124      start(3)=timestep      call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, pk)
     start(4)=0  
     count(2)=jjp1  
     count(3)=1  
     count(4)=0  
   
     !  Interpolation verticale sur les niveaux modele  
   
     call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, &  
          masse, pk)  
   
     !  Passage aux variables du modele (vents covariants, temperature  
     !  potentielle et humidite specifique)  
   
125      call nat2gcm(u, v, t, Q, pk, u, v, t, Q)      call nat2gcm(u, v, t, Q, pk, u, v, t, Q)
     first=.false.  
126    
127    end subroutine read_reanalyse    end subroutine read_reanalyse
128    

Legend:
Removed from v.171  
changed lines
  Added in v.172

  ViewVC Help
Powered by ViewVC 1.1.21