source: BOL/IPCC_AR4/tostdlev.F90 @ 5442

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

Un filtre permettant de relire a posteriori un fichier netcdf contenant
la serie temporelle d'une variable 3d et de l'interpoler sur des niveaux standards
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 7.5 KB
Line 
1
2  PROGRAM tostdlev
3
4!
5! Pour passer a posteriori des champs 3d du modele sur des niveaux standards
6!
7! L. Fairhead 2004/12
8!
9! Ce programme est appelé avec un argument, le nom du fichier à traiter
10! Il nécessite aussi un fichier config.def contenant diverses informations
11! (voir plus bas) et l'accès à un tableau faisant la correspondance entre
12! les noms de variables du modèle et les noms imposés par l'IPCC
13! Pour l'instant on ne traite que les fichiers contenant la serie
14! temporelle d'une seule variable
15
16  implicit none
17
18#include "netcdf.inc"
19
20  INTEGER :: lunout, ierr, icount, ivar, nvars, nlev, len
21  integer :: press_id, data_id, newdata_id
22  INTEGER :: lonid, ilon, newlonid
23  INTEGER :: latid, ilat, newlatid
24  INTEGER :: levid, ilev, newlevid
25  INTEGER :: vartype, ndims
26  INTEGER :: timid, itime, newtimid
27  INTEGER :: varlonid, newvarlonid
28  INTEGER :: varlatid, newvarlatid
29  INTEGER :: vartimid, newvartimid
30  INTEGER :: varlevid, newvarlevid
31  INTEGER :: varid, newvarid
32  integer :: presvar_id
33  INTEGER, dimension(4) :: dimids
34  integer :: natts
35  CHARACTER (len=80) :: varname, attname
36  INTEGER, DIMENSION(4) :: start, count
37 
38  REAL, DIMENSION(:), ALLOCATABLE :: lon, lat, time
39  REAL, DIMENSION(:,:,:), ALLOCATABLE :: pression, champ, newchamp
40
41
42! niveaux standards:
43  INTEGER,parameter :: nlevstd = 17
44  real       :: rlevstd(nlevstd)
45  DATA rlevstd /100000., 92500., 85000., 70000.,             &
46    & 60000., 50000., 40000., 30000., 25000., 20000.,        &
47    & 15000., 10000., 7000., 5000., 3000., 2000., 1000./
48!
49! quelques initialisations
50  lunout = 6
51
52!
53! Ouverture du fichier contenant les pressions modele
54  ierr = nf_open('pression.nc', NF_NOWRITE, press_id)
55  IF (ierr /= NF_NOERR) then
56    WRITE(lunout,*)NF_STRERROR(ierr)
57    stop
58  endif
59
60!
61! Ouverture du fichier contenant le champ a interpoler
62  ierr = nf_open('data.nc', NF_NOWRITE, data_id)
63  IF (ierr /= NF_NOERR) then
64    WRITE(lunout,*)NF_STRERROR(ierr)
65    stop
66  ENDIF
67
68!
69! Ouverture du fichier contenant le champ interpolé
70  ierr = nf_create('newdata.nc', NF_CLOBBER, newdata_id)
71  IF (ierr /= NF_NOERR) then
72    WRITE(lunout,*)NF_STRERROR(ierr)
73    stop
74  endif
75
76!
77! Début lecture fichier origine/ecriture fichier sortie
78!
79! Definition des dimensions:
80  ierr = nf_inq_ndims(data_id,ndims)
81  DO icount = 1, ndims
82    ierr = nf_inq_dim(data_id, icount, varname, len)
83    select case (trim(varname))
84    case ('lon')
85      ilon = len
86!     Definition de la longitude
87      ierr = nf_def_dim(newdata_id,varname,ilon, newlonid)
88    case ('lat')
89      ilat = len
90!     Definition de la latitude
91      ierr = nf_def_dim(newdata_id,varname,ilat, newlatid)
92    case('presnivs')
93      ilev = len
94!     Definition niveaux verticaux
95      ierr = nf_def_dim(newdata_id,'presnivs',nlevstd,newlevid)
96    case ('time_counter')
97      itime = len
98!     Definition du temps
99      ierr = nf_def_dim(newdata_id,varname,itime, newtimid)
100    case default
101      WRITE(lunout,*)'je ne reconnais pas cette dimension: ',varname
102      stop
103    end select
104  enddo
105!
106! Definition des variables
107!
108!
109!
110  ierr = nf_inq_nvars(data_id, nvars)
111  DO ivar = 1, nvars
112    ierr = nf_inq_var(data_id, ivar, varname, vartype, ndims, dimids, natts)
113    if (ierr /= 0) call handle_err(ierr)
114    selectcase (trim(varname))
115    case('lon')
116!     definition de la longitude
117      ierr = nf_inq_varid(data_id,varname, varlonid)
118      dimids(1)=newlonid
119      ierr = nf_def_var(newdata_id, varname, vartype , ndims, dimids,newvarlonid)
120!     recopie des attributs de la variable:
121      DO icount = 1, natts 
122        ierr = nf_inq_attname(data_id, varlonid, icount, attname)
123        ierr = nf_copy_att(data_id, varlonid, attname, newdata_id, newvarlonid)
124      enddo
125    case('lat')
126!     definition de la latitude
127      ierr = nf_inq_varid(data_id, varname, varlatid)
128      dimids(1)=newlatid
129      ierr = nf_def_var(newdata_id, varname, vartype , ndims, dimids,newvarlatid)
130!     recopie des attributs de la variable:
131      DO icount = 1, natts 
132        ierr = nf_inq_attname(data_id, varlatid, icount, attname)
133        ierr = nf_copy_att(data_id, varlatid, attname, newdata_id, newvarlatid)
134      enddo
135    case('time_counter')
136!     definition du temps
137      ierr = nf_inq_varid(data_id, varname, vartimid)
138      dimids(1)=newtimid
139      ierr = nf_def_var(newdata_id, varname, vartype, ndims, dimids,newvartimid)
140!     recopie des attributs de la variable:
141      DO icount = 1, natts 
142        ierr = nf_inq_attname(data_id, vartimid, icount, attname)
143        ierr = nf_copy_att(data_id, vartimid, attname, newdata_id, newvartimid)
144      enddo
145    case('presnivs')
146!     definition des niveaux de pression
147      ierr = nf_inq_varid(data_id,varname,varlevid)
148      dimids(1)=newlevid
149      ierr = nf_def_var(newdata_id,  varname, vartype, ndims, dimids,newvarlevid)
150!     recopie des attributs de la variable:
151      ierr = nf_copy_att(data_id, varlevid, 'units', newdata_id, newvarlevid)
152      ierr = nf_copy_att(data_id, varlevid, 'title', newdata_id, newvarlevid)
153      ierr = nf_copy_att(data_id, varlevid, 'long_name', newdata_id, newvarlevid)
154    case default
155!     normalement il ne reste que la variable à interpoler 4d
156      IF (ndims /= 4) then
157        WRITE(lunout,*)'La variable principale du fichier n''est pas 4D'
158        stop
159      endif
160      ierr = nf_inq_varid(data_id,varname,varid)
161      if (ierr /= 0) call handle_err(ierr)
162      dimids(1) = newlonid
163      dimids(2) = newlatid
164      dimids(3) = newlevid
165      dimids(4) = newtimid
166      ierr = nf_def_var(newdata_id,  varname, vartype, ndims, dimids,newvarid)
167      if (ierr /= 0) call handle_err(ierr)
168      DO icount = 1, natts 
169        ierr = nf_inq_attname(data_id, varid, icount, attname)
170        ierr = nf_copy_att(data_id, varid, attname, newdata_id, newvarid)
171      enddo
172    end select
173  enddo
174!
175
176!
177! fermeture du mode definition du fichier
178  ierr = nf_enddef(newdata_id)
179!
180! ecriture des variables:
181! longitude
182  allocate(lon(ilon))
183  ierr = nf_get_var_real(data_id, varlonid, lon)
184  ierr = nf_put_var_real(newdata_id, newvarlonid,lon)
185!
186! latitude
187  allocate(lat(ilat))
188  ierr = nf_get_var_real(data_id, varlatid, lat)
189  ierr = nf_put_var_real(newdata_id, newvarlatid, lat)
190!
191! niveaux de pression
192  ierr = nf_put_var_real(newdata_id, newvarlevid, rlevstd)
193!
194! temps
195  allocate(time(itime))
196  ierr = nf_get_var_real(data_id, vartimid, time)
197  ierr = nf_put_var_real(newdata_id, newvartimid, time)
198 
199! Traitement de la variable
200  ALLOCATE(pression(ilon,ilat,ilev))
201  ierr = nf_inq_varid(press_id, 'pres', presvar_id)
202  ALLOCATE(champ(ilon,ilat,ilev))
203  ALLOCATE(newchamp(ilon,ilat,nlevstd))
204  start(1) = 1
205  start(2) = 1
206  start(3) = 1
207  COUNT(1) = ilon
208  COUNT(2) = ilat
209  COUNT(3) = ilev
210  count(4) = 1
211  DO icount = 1, itime
212    COUNT(3) = ilev
213    start(4) = icount
214    ierr = nf_get_vara_real(press_id, presvar_id, start, count, pression)
215    if (ierr /= 0) call handle_err(ierr)
216    ierr = nf_get_vara_real(data_id, varid, start, count, champ)
217    if (ierr /= 0) call handle_err(ierr)
218    DO nlev = 1, nlevstd
219      CALL plevel(ilon*ilat, ilev,.TRUE., pression, rlevstd(nlev),champ, &
220 &                newchamp(:,:,nlev))
221    enddo
222    count(3) = nlevstd
223    ierr = nf_put_vara_real(newdata_id, newvarid, start, count, newchamp)
224    if (ierr /= 0) call handle_err(ierr)
225  enddo
226
227! On ferme!
228  ierr = nf_close(press_id)
229  ierr = nf_close(data_id)
230  ierr = nf_close(newdata_id)
231
232contains
233
234subroutine handle_err(status)
235  use netcdf
236
237  implicit none
238  integer, intent(in) :: status
239
240  write(lunout,*)nf90_strerror(status)
241  stop
242
243end subroutine handle_err
244
245
246END PROGRAM tostdlev
247
Note: See TracBrowser for help on using the repository browser.