/[lmdze]/trunk/libf/IOIPSL/flinget.f90
ViewVC logotype

Annotation of /trunk/libf/IOIPSL/flinget.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 48 - (hide annotations)
Tue Jul 19 12:54:20 2011 UTC (12 years, 11 months ago) by guez
File size: 15413 byte(s)
Replaced calls to "flinget" by calls to "NetCDF95".
1 guez 32 MODULE flinget_m
2    
3     ! From flincom.f90, version 2.2 2006/03/07 09:21:51
4    
5     IMPLICIT NONE
6    
7     PRIVATE
8     PUBLIC flinget
9    
10     INTERFACE flinget
11     MODULE PROCEDURE flinget_r3d, flinget_r2d
12 guez 48 ! The difference between the procedures is the rank of argument "var".
13 guez 32 END INTERFACE
14    
15     CONTAINS
16    
17 guez 48 SUBROUTINE flinget_r2d(fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
18     itau_fin, var)
19 guez 32
20 guez 48 INTEGER, intent(in):: fid_in
21     CHARACTER(LEN=*), intent(in):: varname
22     INTEGER, intent(in):: iim, jjm, llm, ttm, itau_dep, itau_fin
23     REAL, intent(out):: var(:, :)
24 guez 32
25 guez 48 ! Local:
26 guez 32 INTEGER :: jl, jj, ji
27 guez 48 REAL, DIMENSION(:), ALLOCATABLE, SAVE :: buff_tmp
28 guez 32 LOGICAL :: check = .FALSE.
29 guez 48
30 guez 32 !---------------------------------------------------------------------
31 guez 48
32 guez 32 IF (.NOT.ALLOCATED(buff_tmp)) THEN
33 guez 48 IF (check) WRITE(*, *) &
34     "flinget_r2d : allocate buff_tmp for buff_sz = ", SIZE(var)
35 guez 32 ALLOCATE (buff_tmp(SIZE(var)))
36     ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
37 guez 48 IF (check) WRITE(*, *) &
38     "flinget_r2d : re-allocate buff_tmp for buff_sz = ", SIZE(var)
39 guez 32 DEALLOCATE (buff_tmp)
40     ALLOCATE (buff_tmp(SIZE(var)))
41     ENDIF
42    
43 guez 48 CALL flinget_mat(fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
44     itau_fin, 1, iim, 1, jjm, buff_tmp)
45 guez 32
46     jl=0
47 guez 48 DO jj=1, SIZE(var, 2)
48     DO ji=1, SIZE(var, 1)
49 guez 32 jl=jl+1
50 guez 48 var(ji, jj) = buff_tmp(jl)
51 guez 32 ENDDO
52     ENDDO
53 guez 48
54 guez 32 END SUBROUTINE flinget_r2d
55    
56 guez 48 !****************************************************************
57 guez 32
58 guez 48 SUBROUTINE flinget_r3d(fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
59     itau_fin, var)
60 guez 32
61     INTEGER, intent(in):: fid_in
62     CHARACTER(LEN=*), intent(in):: varname
63     INTEGER, intent(in):: iim, jjm, llm, ttm, itau_dep, itau_fin
64 guez 48 REAL, intent(out):: var(:, :, :)
65 guez 32
66 guez 48 ! Local:
67 guez 32 INTEGER :: jl, jk, jj, ji
68 guez 48 REAL, DIMENSION(:), ALLOCATABLE, SAVE :: buff_tmp
69 guez 32 LOGICAL :: check = .FALSE.
70 guez 48
71 guez 32 !---------------------------------------------------------------------
72 guez 48
73 guez 32 IF (.NOT.ALLOCATED(buff_tmp)) THEN
74 guez 48 IF (check) WRITE(*, *) &
75     "flinget_r3d : allocate buff_tmp for buff_sz = ", SIZE(var)
76 guez 32 ALLOCATE (buff_tmp(SIZE(var)))
77     ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
78 guez 48 IF (check) WRITE(*, *) &
79     "flinget_r3d : re-allocate buff_tmp for buff_sz = ", SIZE(var)
80 guez 32 DEALLOCATE (buff_tmp)
81     ALLOCATE (buff_tmp(SIZE(var)))
82     ENDIF
83    
84 guez 48 CALL flinget_mat (fid_in, varname, iim, jjm, llm, ttm, &
85     itau_dep, itau_fin, 1, iim, 1, jjm, buff_tmp)
86 guez 32
87     jl=0
88 guez 48 DO jk=1, SIZE(var, 3)
89     DO jj=1, SIZE(var, 2)
90     DO ji=1, SIZE(var, 1)
91 guez 32 jl=jl+1
92 guez 48 var(ji, jj, jk) = buff_tmp(jl)
93 guez 32 ENDDO
94     ENDDO
95     ENDDO
96 guez 48
97 guez 32 END SUBROUTINE flinget_r3d
98    
99 guez 48 !****************************************************************
100 guez 32
101 guez 48 SUBROUTINE flinget_mat(fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
102 guez 32 itau_fin, iideb, iilen, jjdeb, jjlen, var)
103 guez 48
104 guez 32 !- This subroutine will read the variable named varname from
105     !- the file previously opened by flinopen and identified by fid
106    
107     !- It is checked that the dimensions of the variable to be read
108     !- correspond to what the user requested when he specified
109     !- iim, jjm and llm. The only exception which is allowed is
110     !- for compressed data where the horizontal grid is not expected
111     !- to be iim x jjm.
112    
113     !- If variable is of size zero a global attribute is read.
114     !- This global attribute will be of type real
115    
116     !- INPUT
117    
118     !- fid : File ID returned by flinopen
119     !- varname : Name of the variable to be read from the file
120     !- iim : | These three variables give the size of the variables
121     !- jjm : | to be read. It will be verified that the variables
122     !- llm : | fits in there.
123     !- ttm : |
124     !- itau_dep : Time step at which we will start to read
125     !- itau_fin : Time step until which we are going to read
126     !- For the moment this is done on indexes
127     !- but it should be in the physical space.
128     !- If there is no time-axis in the file then use a
129     !- itau_fin < itau_dep, this will tell flinget not to
130     !- expect a time-axis in the file.
131     !- iideb : index i for zoom
132     !- iilen : length of zoom
133     !- jjdeb : index j for zoom
134     !- jjlen : length of zoom
135    
136     !- OUTPUT
137    
138     !- var : array that will contain the data
139 guez 48
140 guez 32 USE strlowercase_m, ONLY : strlowercase
141     USE errioipsl, ONLY : histerr
142     USE netcdf, ONLY : nf90_byte, nf90_double, nf90_float, nf90_get_att, &
143     nf90_get_var, nf90_inquire_attribute, nf90_inquire_dimension, &
144     nf90_inquire_variable, nf90_inq_attname, nf90_inq_varid, nf90_int, &
145     nf90_max_var_dims, nf90_noerr, nf90_short, nf90_strerror
146     use flincom, only: ncids
147    
148     ! ARGUMENTS
149    
150     INTEGER, intent(in):: fid_in
151     CHARACTER(LEN=*), intent(in):: varname
152     INTEGER, intent(in):: iim, jjm, llm, ttm, itau_dep, itau_fin
153     INTEGER :: iideb
154     integer, intent(in):: iilen
155     integer jjdeb
156     integer, intent(in):: jjlen
157     REAL :: var(:)
158    
159     ! LOCAL
160    
161     INTEGER, SAVE :: cind_vid
162     INTEGER, SAVE :: cind_fid
163     INTEGER, SAVE :: cind_len
164 guez 48 INTEGER, DIMENSION(:), ALLOCATABLE, SAVE :: cindex
165     INTEGER, DIMENSION(4) :: w_sta, w_len, w_dim
166 guez 32 INTEGER :: iret, fid
167     INTEGER :: vid, cvid, clen
168     CHARACTER(LEN=70) :: str1
169     CHARACTER(LEN=250) :: att_n, tmp_n
170     INTEGER :: tmp_i
171 guez 48 REAL, SAVE :: mis_v=0.
172 guez 32 REAL :: tmp_r
173     INTEGER :: ndims, x_typ, nb_atts
174 guez 48 INTEGER, DIMENSION(NF90_MAX_VAR_DIMS) :: dimids
175 guez 32 INTEGER :: i, i2d, cnd
176 guez 48 REAL, DIMENSION(:), ALLOCATABLE, SAVE :: var_tmp
177 guez 32 LOGICAL :: uncompress = .FALSE.
178     LOGICAL :: check = .FALSE.
179 guez 48
180 guez 32 !---------------------------------------------------------------------
181 guez 48
182 guez 32 fid = ncids(fid_in)
183    
184     IF (check) THEN
185 guez 48 WRITE(*, *) &
186 guez 32 'flinget_mat : fid_in, fid, varname :', fid_in, fid, TRIM(varname)
187 guez 48 WRITE(*, *) &
188 guez 32 'flinget_mat : iim, jjm, llm, ttm, itau_dep, itau_fin :', &
189     iim, jjm, llm, ttm, itau_dep, itau_fin
190 guez 48 WRITE(*, *) &
191 guez 32 'flinget_mat : iideb, iilen, jjdeb, jjlen :', &
192     iideb, iilen, jjdeb, jjlen
193     ENDIF
194    
195     uncompress = .FALSE.
196    
197     ! 1.0 We get first all the details on this variable from the file
198    
199     vid = -1
200     iret = NF90_INQ_VARID (fid, varname, vid)
201    
202     IF (vid < 0 .OR. iret /= NF90_NOERR) THEN
203 guez 48 CALL histerr (3, 'flinget', &
204     'Variable '//TRIM(varname)//' not found in file', ' ', ' ')
205 guez 32 ENDIF
206    
207     iret = NF90_INQUIRE_VARIABLE (fid, vid, &
208     ndims=ndims, dimids=dimids, nAtts=nb_atts)
209     IF (check) THEN
210 guez 48 WRITE(*, *) &
211 guez 32 'flinget_mat : fid, vid :', fid, vid
212 guez 48 WRITE(*, *) &
213 guez 32 'flinget_mat : ndims, dimids(1:ndims), nb_atts :', &
214     ndims, dimids(1:ndims), nb_atts
215     ENDIF
216    
217     w_dim(:) = 0
218 guez 48 DO i=1, ndims
219 guez 32 iret = NF90_INQUIRE_DIMENSION (fid, dimids(i), len=w_dim(i))
220     ENDDO
221 guez 48 IF (check) WRITE(*, *) &
222 guez 32 'flinget_mat : w_dim :', w_dim(1:ndims)
223    
224     mis_v = 0.0
225    
226     IF (nb_atts > 0) THEN
227     IF (check) THEN
228 guez 48 WRITE(*, *) 'flinget_mat : attributes for variable :'
229 guez 32 ENDIF
230     ENDIF
231 guez 48 DO i=1, nb_atts
232 guez 32 iret = NF90_INQ_ATTNAME (fid, vid, i, att_n)
233     iret = NF90_INQUIRE_ATTRIBUTE (fid, vid, att_n, xtype=x_typ)
234     CALL strlowercase (att_n)
235     IF ( (x_typ == NF90_INT).OR.(x_typ == NF90_SHORT) &
236     .OR.(x_typ == NF90_BYTE) ) THEN
237     iret = NF90_GET_ATT (fid, vid, att_n, tmp_i)
238     IF (check) THEN
239 guez 48 WRITE(*, *) ' ', TRIM(att_n), ' : ', tmp_i
240 guez 32 ENDIF
241     ELSE IF ( (x_typ == NF90_FLOAT).OR.(x_typ == NF90_DOUBLE) ) THEN
242     iret = NF90_GET_ATT (fid, vid, att_n, tmp_r)
243     IF (check) THEN
244 guez 48 WRITE(*, *) ' ', TRIM(att_n), ' : ', tmp_r
245 guez 32 ENDIF
246 guez 48 IF (index(att_n, 'missing_value') > 0) THEN
247 guez 32 mis_v = tmp_r
248     ENDIF
249     ELSE
250     tmp_n = ''
251     iret = NF90_GET_ATT (fid, vid, att_n, tmp_n)
252     IF (check) THEN
253 guez 48 WRITE(*, *) ' ', TRIM(att_n), ' : ', TRIM(tmp_n)
254 guez 32 ENDIF
255     ENDIF
256     ENDDO
257     !?
258     !!!!!!!!!! We will need a verification on the type of the variable
259     !?
260    
261     ! 2.0 The dimensions are analysed to determine what is to be read
262    
263     ! 2.1 the longitudes
264    
265     IF ( w_dim(1) /= iim .OR. w_dim(2) /= jjm) THEN
266     !---
267     !-- There is a possibility that we have to deal with a compressed axis !
268     !---
269     iret = NF90_INQUIRE_DIMENSION (fid, dimids(1), &
270     name=tmp_n, len=clen)
271     iret = NF90_INQ_VARID (fid, tmp_n, cvid)
272     !---
273 guez 48 IF (check) WRITE(*, *) &
274     'Dimname, iret , NF90_NOERR : ', TRIM(tmp_n), iret, NF90_NOERR
275 guez 32 !---
276     !-- If we have an axis which has the same name
277     !-- as the dimension we can see if it is compressed
278     !---
279     !-- TODO TODO for zoom2d
280     !---
281     IF (iret == NF90_NOERR) THEN
282     iret = NF90_GET_ATT (fid, cvid, 'compress', str1)
283     !-----
284     IF (iret == NF90_NOERR) THEN
285 guez 48 iret = NF90_INQUIRE_VARIABLE (fid, cvid, xtype=x_typ, ndims=cnd)
286 guez 32 !-------
287     IF ( cnd /= 1 .AND. x_typ /= NF90_INT) THEN
288 guez 48 CALL histerr (3, 'flinget', &
289 guez 32 'Variable '//TRIM(tmp_n)//' can not be a compressed axis', &
290     'Either it has too many dimensions'// &
291     ' or it is not of type integer', ' ')
292     ELSE
293     !---------
294     !-------- Let us see if we already have that index table
295     !---------
296     IF ( (cind_len /= clen).OR.(cind_vid /= cvid) &
297     .OR.(cind_fid /= fid) ) THEN
298     IF (ALLOCATED(cindex)) DEALLOCATE(cindex)
299     ALLOCATE(cindex(clen))
300     cind_len = clen
301     cind_vid = cvid
302     cind_fid = fid
303     iret = NF90_GET_VAR (fid, cvid, cindex)
304     ENDIF
305     !---------
306     !-------- In any case we need to set the slab of data to be read
307     !---------
308     uncompress = .TRUE.
309     w_sta(1) = 1
310     w_len(1) = clen
311     i2d = 1
312     ENDIF
313     ELSE
314     str1 = 'The horizontal dimensions of '//varname
315 guez 48 CALL histerr (3, 'flinget', str1, &
316 guez 32 'is not compressed and does not'// &
317 guez 48 ' correspond to the requested size', ' ')
318 guez 32 ENDIF
319     ELSE
320     IF (w_dim(1) /= iim) THEN
321     str1 = 'The longitude dimension of '//varname
322 guez 48 CALL histerr (3, 'flinget', str1, &
323 guez 32 'in the file is not equal to the dimension', &
324     'that should be read')
325     ENDIF
326     IF (w_dim(2) /= jjm) THEN
327     str1 = 'The latitude dimension of '//varname
328 guez 48 CALL histerr (3, 'flinget', str1, &
329 guez 32 'in the file is not equal to the dimension', &
330     'that should be read')
331     ENDIF
332     ENDIF
333     ELSE
334     w_sta(1:2) = (/ iideb, jjdeb /)
335     w_len(1:2) = (/ iilen, jjlen /)
336     i2d = 2
337     ENDIF
338    
339     ! 2.3 Now the difficult part, the 3rd dimension which can be
340     ! time or levels.
341    
342     ! Priority is given to the time axis if only three axes are present.
343    
344     IF (ndims > i2d) THEN
345     !---
346     !-- 2.3.1 We have a vertical axis
347     !---
348     IF (llm == 1 .AND. ndims == i2d+2 .OR. llm == w_dim(i2d+1)) THEN
349     !-----
350     IF (w_dim(i2d+1) /= llm) THEN
351 guez 48 CALL histerr (3, 'flinget', &
352 guez 32 'The vertical dimension of '//varname, &
353     'in the file is not equal to the dimension', &
354     'that should be read')
355     ELSE
356     w_sta(i2d+1) = 1
357     IF (llm > 0) THEN
358     w_len(i2d+1) = llm
359     ELSE
360     w_len(i2d+1) = w_sta(i2d+1)
361     ENDIF
362     ENDIF
363     !-----
364     IF ((itau_fin-itau_dep) >= 0) THEN
365     IF (ndims /= i2d+2) THEN
366 guez 48 CALL histerr (3, 'flinget', &
367 guez 32 'You attempt to read a time slab', &
368     'but there is no time axis on this variable', varname)
369     ELSE IF ((itau_fin - itau_dep) <= w_dim(i2d+2)) THEN
370     w_sta(i2d+2) = itau_dep
371     w_len(i2d+2) = itau_fin-itau_dep+1
372     ELSE
373 guez 48 CALL histerr (3, 'flinget', &
374 guez 32 'The time step you try to read is not', &
375     'in the file (1)', varname)
376     ENDIF
377     ELSE IF (ndims == i2d+2 .AND. w_dim(i2d+2) > 1) THEN
378 guez 48 CALL histerr (3, 'flinget', &
379 guez 32 'There is a time axis in the file but no', &
380     'time step give in the call', varname)
381     ELSE
382     w_sta(i2d+2) = 1
383     w_len(i2d+2) = 1
384     ENDIF
385     ELSE
386     !-----
387     !---- 2.3.2 We do not have any vertical axis
388     !-----
389     IF (ndims == i2d+2) THEN
390 guez 48 CALL histerr (3, 'flinget', &
391 guez 32 'The file contains 4 dimensions', &
392     'but only 3 are requestes for variable ', varname)
393     ENDIF
394     IF ((itau_fin-itau_dep) >= 0) THEN
395     IF (ndims == i2d+1) THEN
396     IF ((itau_fin-itau_dep) < w_dim(i2d+1) ) THEN
397     w_sta(i2d+1) = itau_dep
398     w_len(i2d+1) = itau_fin-itau_dep+1
399     ELSE
400 guez 48 CALL histerr (3, 'flinget', &
401 guez 32 'The time step you try to read is not', &
402     'in the file (2)', varname)
403     ENDIF
404     ELSE
405 guez 48 CALL histerr (3, 'flinget', &
406 guez 32 'From your input you sould have 3 dimensions', &
407     'in the file but there are 4', varname)
408     ENDIF
409     ELSE
410     IF (ndims == i2d+1 .AND. w_dim(i2d+1) > 1) THEN
411 guez 48 CALL histerr (3, 'flinget', &
412 guez 32 'There is a time axis in the file but no', &
413     'time step given in the call', varname)
414     ELSE
415     w_sta(i2d+1) = 1
416     w_len(i2d+1) = 1
417     ENDIF
418     ENDIF
419     ENDIF
420     ELSE
421     !---
422     !-- 2.3.3 We do not have any vertical axis
423     !---
424     w_sta(i2d+1:i2d+2) = (/ 0, 0 /)
425     w_len(i2d+1:i2d+2) = (/ 0, 0 /)
426     ENDIF
427    
428     ! 3.0 Reading the data
429    
430 guez 48 IF (check) WRITE(*, *) &
431 guez 32 'flinget_mat 3.0 : ', uncompress, w_sta, w_len
432     !---
433     IF (uncompress) THEN
434     !---
435     IF (ALLOCATED(var_tmp)) THEN
436     IF (SIZE(var_tmp) < clen) THEN
437     DEALLOCATE(var_tmp)
438     ALLOCATE(var_tmp(clen))
439     ENDIF
440     ELSE
441     ALLOCATE(var_tmp(clen))
442     ENDIF
443     !---
444     iret = NF90_GET_VAR (fid, vid, var_tmp, &
445     start=w_sta(:), count=w_len(:))
446     !---
447     var(:) = mis_v
448     var(cindex(:)) = var_tmp(:)
449     !---
450     ELSE
451     iret = NF90_GET_VAR (fid, vid, var, &
452     start=w_sta(:), count=w_len(:))
453     ENDIF
454    
455 guez 48 IF (check) WRITE(*, *) 'flinget_mat 3.1 : ', NF90_STRERROR (iret)
456    
457 guez 32 END SUBROUTINE flinget_mat
458    
459     END MODULE flinget_m

  ViewVC Help
Powered by ViewVC 1.1.21