source: IOIPSL/trunk/src/flincom.f90 @ 4

Last change on this file since 4 was 4, checked in by rblod, 19 years ago

First import of IOIPSL sources

File size: 57.0 KB
Line 
1!$Header: /home/ioipsl/CVSROOT/IOIPSL/src/flincom.f90,v 2.2 2006/03/07 09:21:51 adm Exp $
2!-
3MODULE flincom
4!---------------------------------------------------------------------
5  USE netcdf
6!-
7  USE calendar,  ONLY : ju2ymds, ymds2ju, ioconf_calendar
8  USE errioipsl, ONLY : histerr
9  USE stringop,  ONLY : strlowercase
10!-
11  IMPLICIT NONE
12!-
13  PRIVATE
14  PUBLIC :: flinput, flincre, flinget, flinclo, &
15            flinopen, flininfo, flininspect, flinquery_var
16!-
17  INTERFACE flinopen
18!---------------------------------------------------------------------
19!- The "flinopen" routines will open an input file
20!-
21!- INPUT
22!-
23!- filename  : Name of the netCDF file to be opened
24!-
25!- iideb     : index i for zoom     !
26!- iilen     : length  of  zoom     !   for
27!- jjdeb     : index j for zoom     !   zoom
28!- jjlen     : length  of  zoom     !
29!-
30!- do_test   : A flag that enables the testing of the content
31!-             of the file against the input from the model
32!-
33!- INPUT if do_test=TRUE      OUTPUT else
34!-
35!- iim       : size in the x direction in the file (longitude)
36!- jjm       : size in the y direction
37!- llm       : number of levels
38!-             (llm = 0 means no axis to be expected)
39!- lon       : array of (iilen,jjlen) (zoom), or (iim,jjm) (no zoom),
40!-             that contains the longitude of each point
41!- lat       : same for latitude
42!- lev       : An array of llm for the latitude
43!-
44!- WARNING :
45!- In the case of do_test=FALSE it is for the user to check
46!- that the dimensions of lon lat and lev are correct when passed to
47!- flinopen. This can be done after the call when iim and jjm have
48!- been retrieved from the netCDF file. In F90 this problem will
49!- be solved with an internal assign
50!- IF iim, jjm, llm or ttm are parameters in the calling program and
51!- you use the option do_test=FALSE it will create a segmentation fault
52!-
53!-  OUTPUT
54!-
55!-  ttm       : size of time axis
56!-  itaus     : Time steps within this file
57!-  date0     : Julian date at which itau = 0
58!-  dt        : length of the time steps of the data
59!-  fid       : returned file ID which is later used to read the data
60!---------------------------------------------------------------------
61    MODULE PROCEDURE flinopen_zoom2d, flinopen_nozoom
62  END INTERFACE
63!-
64  INTERFACE flinput
65!---------------------------------------------------------------------
66!- The "flinput" routines will put a variable
67!- on the netCDF file created by flincre.
68!- If the sizes of the axis do not match the one of the IDs
69!- then a new axis is created.
70!- That is we loose the possibility of writting hyperslabs of data.
71!-
72!- Again here if iim = jjm = llm = ttm = 0
73!- then a global attribute is added to the file.
74!-
75!- INPUT
76!-
77!- fid      : Identification of the file in which we will write
78!- varname  : Name of variable to be written
79!- iim      : size in x of variable
80!- nlonid   : ID of x axis which could fit for this axis
81!- jjm      : size in y of variable
82!- nlatid   : ID of y axis which could fit for this axis
83!- llm      : size in z of variable
84!- zdimid   : ID of z axis which could fit for this axis
85!- ttm      : size in t of variable
86!- tdimid   : ID of t axis which could fit for this axis
87!-
88!- OUTPUT
89!-
90!- NONE
91!---------------------------------------------------------------------
92    MODULE PROCEDURE flinput_r4d, flinput_r3d, flinput_r2d, &
93                     flinput_r1d, flinput_scal
94  END INTERFACE
95!-
96  INTERFACE flinget
97    MODULE PROCEDURE flinget_r4d, flinget_r3d, flinget_r2d, &
98                     flinget_r1d, flinget_scal, &
99                     flinget_r4d_zoom2d, flinget_r3d_zoom2d, &
100                     flinget_r2d_zoom2d
101  END INTERFACE
102!-
103! This is the data we keep on each file we open
104!-
105  INTEGER, PARAMETER :: nbfile_max = 200
106  INTEGER, SAVE :: nbfiles = 0
107  INTEGER, SAVE :: ncids(nbfile_max), ncnbd(nbfile_max), &
108                   ncfunli(nbfile_max), ncnba(nbfile_max)
109  INTEGER, SAVE :: ncnbva(nbfile_max), ncdims(nbfile_max,4)
110  LOGICAL, SAVE :: ncfileopen(nbfile_max)=.FALSE.
111!-
112  INTEGER, SAVE :: cind_vid, cind_fid, cind_len
113  INTEGER,DIMENSION(:),ALLOCATABLE,SAVE :: cindex
114!-
115  INTEGER,DIMENSION(4) :: w_sta, w_len, w_dim
116!-
117CONTAINS
118!-
119!===
120!-
121SUBROUTINE flincre &
122  (filename, iim1, jjm1, lon1, lat1, llm1, lev1, ttm1, itaus, &
123   time0, dt, fid_out, nlonid1, nlatid1, zdimid1, tdimid1)
124!---------------------------------------------------------------------
125!- This is a "low level" subroutine for opening netCDF files wich
126!- contain the major coordinate system of the model.
127!- Other coordinates needed for other variables
128!- will be added as they are needed.
129!-
130!- INPUT
131!-
132!- filename    : Name of the file to be created
133!- iim1, jjm1  : Horizontal size of the grid
134!-               which will be stored in the file
135!- lon1, lat1  : Horizontal grids
136!- llm1        : Size of the vertical grid
137!- lev1        : Vertical grid
138!- ttm1        : Size of time axis
139!- itaus       : time steps on the time axis
140!- time0       : Time in julian days at which itau = 0
141!- dt          : time step in seconds between itaus
142!-               (one step of itau)
143!-
144!- OUTPUT
145!-
146!- fid         : File identification
147!- nlonid1     : Identification of longitudinal axis
148!- nlatid1     : Identification of latitudinal axis
149!- zdimid1     : ID of vertical axis
150!- tdimid1     : ID of time axis
151!---------------------------------------------------------------------
152  IMPLICIT NONE
153!-
154! ARGUMENTS
155!-
156  CHARACTER(LEN=*) :: filename
157  INTEGER :: iim1, jjm1, llm1, ttm1
158  REAL :: lon1(iim1,jjm1)
159  REAL :: lat1(iim1,jjm1)
160  REAL :: lev1(llm1)
161  INTEGER :: itaus(ttm1)
162  REAL :: time0
163  REAL :: dt
164  INTEGER :: fid_out, zdimid1, nlonid1, nlatid1, tdimid1
165!-
166! LOCAL
167!-
168  INTEGER :: iret, lll, fid
169  INTEGER :: lonid, latid, levid, timeid
170  INTEGER :: year, month, day
171  REAL :: sec
172  CHARACTER(LEN=250):: name
173!-
174  LOGICAL :: check = .FALSE.
175!---------------------------------------------------------------------
176  lll = LEN_TRIM(filename)
177  IF (filename(lll-2:lll) /= '.nc') THEN
178    name=filename(1:lll)//'.nc'
179  ELSE
180    name=filename(1:lll)
181  ENDIF
182!-
183  iret = NF90_CREATE (name, NF90_CLOBBER, fid)
184!-
185  iret = NF90_DEF_DIM (fid, 'x',     iim1, nlonid1)
186  iret = NF90_DEF_DIM (fid, 'y',     jjm1, nlatid1)
187  iret = NF90_DEF_DIM (fid, 'lev',   llm1, zdimid1)
188  iret = NF90_DEF_DIM (fid, 'tstep', ttm1, tdimid1)
189!-
190! Vertical axis
191!-
192  IF (check) WRITE(*,*) 'flincre Vertical axis'
193!-
194  iret = NF90_DEF_VAR (fid, 'lev', NF90_FLOAT, zdimid1, levid)
195  iret = NF90_PUT_ATT (fid, levid, 'units',     '-')
196  iret = NF90_PUT_ATT (fid, levid, 'title',     'levels')
197  iret = NF90_PUT_ATT (fid, levid, 'long_name', 'Sigma Levels')
198!-
199! Time axis
200!-
201  IF (check) WRITE(*,*) 'flincre time axis'
202!-
203  iret = NF90_DEF_VAR (fid, 'tstep', NF90_FLOAT, tdimid1, timeid)
204  iret = NF90_PUT_ATT (fid, timeid, 'units',     '-')
205  iret = NF90_PUT_ATT (fid, timeid, 'title',     'time')
206  iret = NF90_PUT_ATT (fid, timeid, 'long_name', 'time steps')
207!-
208! The longitude
209!-
210  IF (check) WRITE(*,*) 'flincre Longitude axis'
211!-
212  iret = NF90_DEF_VAR (fid, "nav_lon", NF90_FLOAT, &
213                       (/ nlonid1, nlatid1 /), lonid)
214  iret = NF90_PUT_ATT (fid, lonid, 'units', "degrees_east")
215  iret = NF90_PUT_ATT (fid, lonid, 'title', "Longitude")
216  iret = NF90_PUT_ATT (fid, lonid, 'nav_model', &
217                       "Lambert projection of PROMES")
218  iret = NF90_PUT_ATT (fid, lonid, 'valid_min', &
219                       REAL(MINVAL(lon1),KIND=4))
220  iret = NF90_PUT_ATT (fid, lonid, 'valid_max', &
221                       REAL(MAXVAL(lon1),KIND=4))
222!-
223! The Latitude
224!-
225  IF (check) WRITE(*,*) 'flincre Latitude axis'
226!-
227  iret = NF90_DEF_VAR (fid, "nav_lat", NF90_FLOAT, &
228                     (/ nlonid1, nlatid1 /), latid)
229  iret = NF90_PUT_ATT (fid, latid, 'units', "degrees_north")
230  iret = NF90_PUT_ATT (fid, latid, 'title', "Latitude")
231  iret = NF90_PUT_ATT (fid, latid, 'nav_model', &
232                       "Lambert projection of PROMES")
233  iret = NF90_PUT_ATT (fid, latid, 'valid_min', &
234                       REAL(MINVAL(lat1),KIND=4))
235  iret = NF90_PUT_ATT (fid, latid, 'valid_max', &
236                       REAL(MAXVAL(lat1),KIND=4))
237!-
238! The time coordinates
239!-
240  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, 'delta_tstep_sec', &
241                       REAL(dt,KIND=4))
242!-
243  CALL ju2ymds (time0, year, month, day, sec)
244!-
245  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, 'year0',  REAL(year,KIND=4))
246  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, 'month0', REAL(month,KIND=4))
247  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, 'day0',   REAL(day,KIND=4))
248  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, 'sec0',   REAL(sec,KIND=4))
249!-
250  iret = NF90_ENDDEF (fid)
251!-
252  IF (check) WRITE(*,*) 'flincre Variable'
253!-
254  iret = NF90_PUT_VAR (fid, levid, lev1(1:llm1))
255!-
256  IF (check) WRITE(*,*) 'flincre Time Variable'
257!-
258  iret = NF90_PUT_VAR (fid, timeid, REAL(itaus(1:ttm1)))
259!-
260  IF (check) WRITE(*,*) 'flincre Longitude'
261!-
262  iret = NF90_PUT_VAR (fid, lonid, lon1(1:iim1,1:jjm1))
263!-
264  IF (check) WRITE(*,*) 'flincre Latitude'
265!-
266  iret = NF90_PUT_VAR (fid, latid, lat1(1:iim1,1:jjm1))
267!-
268! Keep all this information
269!-
270  nbfiles = nbfiles+1
271!-
272  IF (nbfiles > nbfile_max) THEN
273    CALL histerr (3,'flincre', &
274     'Too many files. Please increase nbfil_max', &
275     'in program flincom.F90.',' ')
276  ENDIF
277!-
278  ncids(nbfiles) = fid
279  ncnbd(nbfiles) = 4
280!-
281  ncdims(nbfiles,1:4) = (/ iim1, jjm1, llm1, ttm1 /)
282!-
283  ncfunli(nbfiles) = -1
284  ncnba(nbfiles)   =  4
285  ncnbva(nbfiles)  =  0
286  ncfileopen(nbfiles) = .TRUE.
287!-
288  fid_out = nbfiles
289!---------------------
290END SUBROUTINE flincre
291!-
292!===
293!-
294SUBROUTINE flinopen_zoom2d &
295  (filename, iideb, iilen, jjdeb, jjlen, do_test, &
296  iim, jjm, llm, lon, lat, lev, ttm, itaus, date0, dt, fid_out)
297!---------------------------------------------------------------------
298  IMPLICIT NONE
299!-
300! ARGUMENTS
301!-
302  CHARACTER(LEN=*) :: filename
303  LOGICAL :: do_test
304  INTEGER :: iim, jjm, llm, ttm, iideb, iilen, jjdeb, jjlen
305  REAL :: lon(iilen,jjlen), lat(iilen,jjlen), lev(llm)
306  INTEGER :: itaus(ttm)
307  REAL :: date0, dt
308  INTEGER :: fid_out
309!-
310  LOGICAL :: check = .FALSE.
311!---------------------------------------------------------------------
312  IF (check) WRITE (*,*) ' iideb, iilen, jjdeb, jjlen, iim, jjm ', &
313                           iideb, iilen, jjdeb, jjlen, iim, jjm
314  IF (check) WRITE (*,*) ' lon ', lon(1,1), lon(iilen,jjlen)
315  IF (check) WRITE (*,*) ' lat ', lat(1,1), lat(iilen,jjlen)
316!-
317  CALL flinopen_work &
318    (filename, iideb, iilen, jjdeb, jjlen, do_test, &
319     iim, jjm, llm, lon, lat, lev, ttm, itaus, date0, dt, fid_out)
320!-----------------------------
321END SUBROUTINE flinopen_zoom2d
322!-
323!===
324!-
325SUBROUTINE flinopen_nozoom &
326  (filename, do_test, iim, jjm, llm, lon, lat, lev, ttm, &
327   itaus, date0, dt, fid_out)
328!---------------------------------------------------------------------
329  IMPLICIT NONE
330!-
331! ARGUMENTS
332!-
333  CHARACTER(LEN=*) :: filename
334  LOGICAL :: do_test
335  INTEGER :: iim, jjm, llm, ttm
336  REAL :: lon(iim,jjm), lat(iim,jjm), lev(llm)
337  INTEGER :: itaus(ttm)
338  REAL :: date0, dt
339  INTEGER :: fid_out
340!---------------------------------------------------------------------
341  CALL flinopen_work &
342    (filename, 1, iim, 1, jjm, do_test, &
343     iim, jjm, llm, lon, lat, lev, ttm, itaus, date0, dt, fid_out)
344!-------------------------
345END SUBROUTINE flinopen_nozoom
346!-
347!===
348!-
349SUBROUTINE flinopen_work &
350  (filename, iideb, iilen, jjdeb, jjlen, do_test, &
351   iim, jjm, llm, lon, lat, lev, ttm, itaus, date0, dt, fid_out)
352!---------------------------------------------------------------------
353  IMPLICIT NONE
354!-
355! ARGUMENTS
356!-
357  CHARACTER(LEN=*) :: filename
358  LOGICAL :: do_test
359  INTEGER :: iim, jjm, llm, ttm, iideb, iilen, jjdeb, jjlen
360  REAL :: lon(iilen,jjlen), lat(iilen,jjlen), lev(llm)
361  INTEGER :: itaus(ttm)
362  REAL :: date0, dt
363  INTEGER :: fid_out
364!-
365! LOCAL
366!-
367  REAL, PARAMETER :: eps = 1.e-4
368!-
369  INTEGER :: iret, vid, fid, nbdim, i, iilast, jjlast
370  INTEGER :: it, len, gdtt_id, old_id, iv, gdtmaf_id
371  CHARACTER(LEN=250) :: name
372  CHARACTER(LEN=80) :: units, calendar
373  INTEGER :: tmp_iim, tmp_jjm, tmp_llm, tmp_ttm
374  REAL :: x_first, x_last
375  INTEGER :: year, month, day
376  REAL :: r_year, r_month, r_day
377  INTEGER :: year0, month0, day0, hours0, minutes0, seci
378  REAL :: sec, sec0
379  CHARACTER :: strc
380!-
381  REAL,DIMENSION(:),ALLOCATABLE :: vec_tmp
382!-
383  LOGICAL :: open_file
384  LOGICAL :: check = .FALSE.
385!---------------------------------------------------------------------
386  iilast = iideb+iilen-1
387  jjlast = jjdeb+jjlen-1
388  IF (check) WRITE (*,*) &
389    ' flinopen_work zoom 2D information '// &
390    ' iideb, iilen, iilast, jjdeb, jjlen, jjlast ', &
391      iideb, iilen, iilast, jjdeb, jjlen, jjlast
392!-
393! 1.0 get all infos on the file  we if it is not yet opened
394!-
395! Either the fid_out has not been initialized (0 or very large)
396! then we have to open anyway. Else we only need to open the file
397! if it has not been opened before.
398!-
399  IF ( (fid_out < 1).OR.(fid_out > nbfile_max) ) THEN
400    open_file = .TRUE.
401  ELSE IF (.NOT.ncfileopen(fid_out)) THEN
402    open_file = .TRUE.
403  ELSE
404    open_file = .FALSE.
405  ENDIF
406!-
407  IF (open_file) THEN
408    CALL flininfo (filename,tmp_iim,tmp_jjm,tmp_llm,tmp_ttm,fid_out)
409  ELSE
410!-- The user has already opened the file
411!-- and we trust that he knows the dimensions
412    tmp_iim = iim
413    tmp_jjm = jjm
414    tmp_llm = llm
415    tmp_ttm = ttm
416  ENDIF
417!-
418  IF (check) &
419    WRITE(*,*) 'OUT OF flininfo :',tmp_iim,tmp_jjm,tmp_llm,tmp_ttm
420!-
421  fid = ncids(fid_out)
422!-
423! 2.0 get the sizes and names of the different coordinates
424!     and do a first set of verification.
425!-
426! 2.2 We test the axis if we have to.
427!-
428  IF (check) &
429    WRITE(*,*) 'flininfo 2.2 We test if we have to test : ',do_test
430!-
431  IF (do_test) THEN
432    IF      (iim /= tmp_iim) THEN
433      CALL histerr (3,'flinopen', &
434        'file '//filename//' does not have the ', &
435        'required dimension in x direction (longitude)',' ')
436    ELSE IF (jjm /= tmp_jjm) THEN
437      CALL histerr (3,'flinopen', &
438        'file '//filename//' does not have the ', &
439        'required dimension in y direction (latitude)',' ')
440    ELSE IF ( llm /= tmp_llm .AND. llm > 0 ) THEN
441      CALL histerr (3,'flinopen', &
442        'file '//filename//' does not have the ', &
443        'required dimension in the vertical',' ')
444    ENDIF
445  ELSE
446!---
447!-- 2.3 Else the sizes of the axes are returned to the user
448!---
449    IF (check) WRITE(*,*) 'flinopen 2.3 Else sizes are returned'
450!---
451    iim = tmp_iim
452    jjm = tmp_jjm
453    llm = tmp_llm
454  ENDIF
455!-
456  ttm = tmp_ttm
457!-
458! 3.0 Check if we are realy talking about the same coodinate system
459!     if not then we get the lon, lat and lev variables from the file
460!-
461  IF (check) WRITE(*,*) 'flinopen 3.0 we are realy talking'
462!-
463  IF (do_test) THEN
464!---
465    CALL flinfindcood (fid_out, 'lon', vid, nbdim)
466    iret = NF90_GET_VAR (fid, vid, x_first, start=(/ iideb, jjdeb /))
467    iret = NF90_GET_VAR (fid, vid, x_last, start=(/ iilast, jjlast /))
468!---
469    IF (check) &
470      WRITE(*,*) 'from file lon first and last, modulo 360. ', &
471        x_first, x_last, MODULO(x_first,360.), MODULO(x_last,360.)
472    IF (check) &
473      WRITE(*,*) 'from model lon first and last, modulo 360. ', &
474        lon(1,1),lon(iilen,jjlen), &
475        MODULO(lon(1,1),360.), MODULO(lon(iilen,jjlen),360.)
476    IF (    (ABS( MODULO(x_first,360.) &
477                 -MODULO(lon(1,1),360.)) > eps) &
478        .OR.(ABS( MODULO(x_last,360.) &
479                 -MODULO(lon(iilen ,jjlen),360.)) > eps ) ) THEN
480      CALL histerr (3,'flinopen', &
481        'file '//filename//' and the model do not', &
482        'share the same longitude coordinate', &
483        'Obtained by comparing the first and last values ')
484    ENDIF
485!---
486    CALL flinfindcood (fid_out, 'lat', vid, nbdim)
487    iret = NF90_GET_VAR (fid, vid, x_first, start=(/ iideb, jjdeb /))
488    iret = NF90_GET_VAR (fid, vid, x_last, start=(/ iilast, jjlast /))
489!---
490    IF (check) WRITE(*,*) &
491      'from file lat first and last ',x_first,x_last
492    IF (check) WRITE(*,*) &
493      'from model lat first and last ',lat(1,1),lat(iilen,jjlen)
494!---
495    IF (    (ABS(x_first-lat(1,1)) > eps) &
496        .OR.(ABS(x_last-lat(iilen,jjlen)) > eps) ) THEN
497      CALL histerr (3,'flinopen', &
498        'file '//filename//' and the model do not', &
499        'share the same latitude coordinate', &
500        'Obtained by comparing the first and last values ')
501    ENDIF
502!---
503    IF (llm > 0) THEN
504      CALL flinfindcood (fid_out, 'lev', vid, nbdim)
505      iret = NF90_GET_VAR (fid, vid, x_first, start=(/ 1 /))
506      iret = NF90_GET_VAR (fid, vid, x_last, start=(/ llm /))
507!-----
508      IF (check) WRITE(*,*) &
509        'from file lev first and last ',x_first ,x_last
510      IF (check) WRITE(*,*) &
511        'from model lev first and last ',lev(1),lev(llm)
512!-----
513      IF (    (ABS(x_first-lev(1)) > eps) &
514          .OR.(ABS(x_last-lev(llm)) > eps) ) THEN
515        CALL histerr (3,'flinopen', &
516          'file '//filename//' and the model do not', &
517          'share the same vertical coordinate', &
518          'Obtained by comparing the first and last values')
519      ENDIF
520    ENDIF
521!---
522  ELSE
523!---
524!-- 4.0 extracting the coordinates if we do not check
525!---
526    IF (check) WRITE(*,*) 'flinopen 4.0 extracting the coordinates'
527!---
528    CALL flinfindcood (fid_out, 'lon', vid, nbdim)
529    IF (nbdim == 2) THEN
530      iret = NF90_GET_VAR (fid, vid, lon, &
531               start=(/ iideb, jjdeb /), count=(/ iilen, jjlen /))
532    ELSE
533      ALLOCATE(vec_tmp(iilen))
534      iret = NF90_GET_VAR (fid, vid, vec_tmp, &
535               start=(/ iideb /), count=(/ iilen /))
536      DO i=1,jjlen
537        lon(:,i) = vec_tmp(:)
538      ENDDO
539      DEALLOCATE(vec_tmp)
540    ENDIF
541!---
542    CALL flinfindcood (fid_out, 'lat', vid, nbdim)
543    IF (nbdim == 2) THEN
544      iret = NF90_GET_VAR (fid, vid, lat, &
545               start=(/ iideb, jjdeb /), count=(/ iilen, jjlen /))
546    ELSE
547      ALLOCATE(vec_tmp(jjlen))
548      iret = NF90_GET_VAR (fid, vid, vec_tmp, &
549               start=(/ jjdeb /), count=(/ jjlen /))
550      DO i=1,iilen
551        lat(i,:) = vec_tmp(:)
552      ENDDO
553      DEALLOCATE(vec_tmp)
554    ENDIF
555!---
556    IF (llm > 0) THEN
557      CALL flinfindcood (fid_out, 'lev', vid, nbdim)
558      IF (nbdim == 1) THEN
559        iret = NF90_GET_VAR (fid, vid, lev, &
560                 start=(/ 1 /), count=(/ llm /))
561      ELSE
562        CALL histerr (3,'flinopen', &
563          'Can not handle vertical coordinates that have more',&
564          'than 1 dimension',' ')
565      ENDIF
566    ENDIF
567  ENDIF
568!-
569! 5.0 Get all the details for the time if possible needed
570!-
571  IF (check) WRITE(*,*) 'flinopen 5.0 Get time'
572!-
573  IF (ttm > 0) THEN
574!---
575!-- 5.1 Find the time axis. Prefered method is the 'timestep since'
576!---
577    gdtmaf_id = -1
578    gdtt_id = -1
579    old_id = -1
580    DO iv=1,ncnbva(fid_out)
581      name=''
582      iret = NF90_INQUIRE_VARIABLE (fid, iv, name=name)
583      units=''
584      iret = NF90_GET_ATT (fid, iv, 'units', units)
585      IF (INDEX(units,'seconds since') > 0) gdtmaf_id = iv
586      IF (INDEX(units,'timesteps since') > 0) gdtt_id = iv
587      IF (INDEX(name, 'tstep') > 0) old_id = iv
588    ENDDO
589!---
590    IF (gdtt_id > 0) THEN
591      vid = gdtt_id
592    ELSE IF (gdtmaf_id > 0) THEN
593      vid = gdtmaf_id
594    ELSE IF (old_id > 0) THEN
595      vid = old_id
596    ELSE
597      CALL histerr (3, 'flinopen', 'No time axis found',' ',' ')
598    ENDIF
599!---
600    ALLOCATE(vec_tmp(ttm))
601    iret = NF90_GET_VAR (fid, vid, vec_tmp, &
602             start=(/ 1 /), count=(/ ttm /))
603    itaus(1:ttm) = NINT(vec_tmp(1:ttm))
604    DEALLOCATE(vec_tmp)
605!---
606    IF (check) WRITE(*,*) 'flinopen 5.1 Times ',itaus
607!---
608!-- Getting all the details for the time axis
609!---
610!-- Find the calendar
611    calendar='XXXX'
612    iret = NF90_GET_ATT (fid, gdtmaf_id, 'calendar', calendar)
613    IF ( INDEX(calendar,'XXXX') < 1 ) THEN
614       CALL ioconf_calendar(calendar)
615    ENDIF
616!--
617    units = ''
618    iret = NF90_GET_ATT (fid, vid, 'units', units)
619    IF (gdtt_id > 0) THEN
620      units = units(INDEX(units,'since')+6:LEN_TRIM(units))
621      READ (units,'(I4.4,5(a,I2.2))') &
622        year0, strc, month0, strc, day0, &
623               strc, hours0, strc, minutes0, strc, seci
624      sec0 = hours0*3600. + minutes0*60. + seci
625      CALL ymds2ju (year0, month0, day0, sec0, date0)
626      IF (check) &
627        WRITE(*,*) 'flinopen 5.1 gdtt_id year0 ... date0 ', &
628                   year0, month0, day0, sec0, date0
629!-----
630      iret = NF90_GET_ATT (fid, gdtt_id, 'tstep_sec', dt)
631    ELSE IF (gdtmaf_id > 0) THEN
632      units = units(INDEX(units,'since')+6:LEN_TRIM(units))
633      READ (units,'(I4.4,5(a,I2.2))') &
634        year0, strc, month0, strc, day0, &
635               strc, hours0, strc, minutes0, strc, seci
636      sec0 = hours0*3600. + minutes0*60. + seci
637      CALL ymds2ju (year0, month0, day0, sec0, date0)
638!-----
639      IF (check) &
640        WRITE(*,*) 'flinopen 5.1 gdtmaf_id year0 ... date0 ', &
641                   year0, month0, day0, sec0, date0
642    ELSE IF (old_id > 0) THEN
643      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'delta_tstep_sec', dt)
644      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'day0', r_day)
645      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'sec0', sec)
646      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'year0', r_year)
647      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'month0', r_month)
648!-----
649      day = INT(r_day)
650      month = INT(r_month)
651      year = INT(r_year)
652!-----
653      CALL ymds2ju (year, month, day, sec, date0)
654    ENDIF
655  ENDIF
656!-
657  IF (check) WRITE(*,*) 'flinopen 6.0 File opened', date0, dt
658!---------------------------
659END SUBROUTINE flinopen_work
660!-
661!===
662!-
663SUBROUTINE flininfo (filename, iim, jjm, llm, ttm, fid_out)
664!---------------------------------------------------------------------
665!- This subroutine allows to get some information.
666!- It is usualy done within flinopen but the user may want to call
667!- it before in order to allocate the space needed to extract the
668!- data from the file.
669!---------------------------------------------------------------------
670  IMPLICIT NONE
671!-
672! ARGUMENTS
673!-
674  CHARACTER(LEN=*) :: filename
675  INTEGER :: iim, jjm, llm, ttm, fid_out
676!-
677! LOCAL
678!-
679  INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim
680  INTEGER :: iv, lll
681  INTEGER :: xid, yid, zid, tid
682  CHARACTER(LEN=80) :: name
683  CHARACTER(LEN=30) :: axname
684!-
685  LOGICAL :: check = .FALSE.
686!---------------------------------------------------------------------
687  lll = LEN_TRIM(filename)
688  IF (filename(lll-2:lll) /= '.nc') THEN
689    name = filename(1:lll)//'.nc'
690  ELSE
691    name = filename(1:lll)
692  ENDIF
693!-
694  iret = NF90_OPEN (name, NF90_NOWRITE, fid)
695  IF (iret /= NF90_NOERR) THEN
696    CALL histerr(3, 'flininfo','Could not open file :',TRIM(name),' ')
697  ENDIF
698!-
699  iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, &
700                      nAttributes=nb_atts, unlimitedDimId=id_unlim)
701!-
702  xid = -1; iim = 0;
703  yid = -1; jjm = 0;
704  zid = -1; llm = 0;
705  tid = -1; ttm = 0;
706!-
707  DO iv=1,ndims
708!---
709    iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll)
710    CALL strlowercase (axname)
711    axname = ADJUSTL(axname)
712!---
713    IF (check) WRITE(*,*) &
714      'flininfo - getting axname',iv,axname,lll
715!---
716    IF      (    (INDEX(axname,'x') == 1) &
717             .OR.(INDEX(axname,'lon') == 1) ) THEN
718      xid = iv; iim = lll;
719    ELSE IF (    (INDEX(axname,'y') == 1) &
720             .OR.(INDEX(axname,'lat') == 1) ) THEN
721      yid = iv; jjm = lll;
722    ELSE IF (    (INDEX(axname,'lev') == 1) &
723             .OR.(INDEX(axname,'plev') == 1) &
724             .OR.(INDEX(axname,'z') == 1) &
725             .OR.(INDEX(axname,'depth') == 1) ) THEN
726      zid = iv; llm = lll;
727    ELSE IF (    (INDEX(axname,'tstep') == 1) &
728             .OR.(INDEX(axname,'time_counter') == 1) ) THEN
729!---- For the time we certainly need to allow for other names
730      tid = iv; ttm = lll;
731    ELSE IF (ndims == 1) THEN
732!---- Nothing was found and ndims=1 then we have a vector of data
733      xid = 1; iim = lll;
734    ENDIF
735!---
736  ENDDO
737!-
738! Keep all this information
739!-
740  nbfiles = nbfiles+1
741!-
742  IF (nbfiles > nbfile_max) THEN
743    CALL histerr (3,'flininfo', &
744      'Too many files. Please increase nbfil_max', &
745      'in program flincom.F90.',' ')
746  ENDIF
747!-
748  ncids(nbfiles) = fid
749  ncnbd(nbfiles) = ndims
750!-
751  ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /)
752!-
753  ncfunli(nbfiles) = id_unlim
754  ncnba(nbfiles)   = nb_atts
755  ncnbva(nbfiles)  = nvars
756  ncfileopen(nbfiles) = .TRUE.
757!-
758  fid_out = nbfiles
759!----------------------
760END SUBROUTINE flininfo
761!-
762!===
763!-
764SUBROUTINE flinput_r1d &
765  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
766!---------------------------------------------------------------------
767  IMPLICIT NONE
768!-
769  INTEGER :: fid_in
770  CHARACTER(LEN=*) :: varname
771  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
772  REAL :: var(:)
773!-
774  INTEGER :: fid, ncvarid, ndim, iret
775  LOGICAL :: check = .FALSE.
776!---------------------------------------------------------------------
777  IF (check) WRITE(*,*) &
778     "flinput_r1d : SIZE(var) = ",SIZE(var)
779!-
780  CALL flinput_mat &
781    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
782     fid,ncvarid,ndim)
783!-
784  iret = NF90_PUT_VAR (fid, ncvarid, var, &
785           start=w_sta(1:ndim), count=w_len(1:ndim))
786!-------------------------
787END SUBROUTINE flinput_r1d
788!-
789!===
790!-
791SUBROUTINE flinput_r2d &
792  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
793!---------------------------------------------------------------------
794  IMPLICIT NONE
795!-
796  INTEGER :: fid_in
797  CHARACTER(LEN=*) :: varname
798  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
799  REAL :: var(:,:)
800!-
801  INTEGER :: fid, ncvarid, ndim, iret
802  LOGICAL :: check = .FALSE.
803!---------------------------------------------------------------------
804  IF (check) WRITE(*,*) &
805     "flinput_r2d : SIZE(var) = ",SIZE(var)
806!-
807  CALL flinput_mat &
808    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
809     fid,ncvarid,ndim)
810!-
811  iret = NF90_PUT_VAR (fid, ncvarid, var, &
812           start=w_sta(1:ndim), count=w_len(1:ndim))
813!-------------------------
814END SUBROUTINE flinput_r2d
815!-
816!===
817!-
818SUBROUTINE flinput_r3d &
819  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
820!---------------------------------------------------------------------
821  IMPLICIT NONE
822!-
823  INTEGER :: fid_in
824  CHARACTER(LEN=*) :: varname
825  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
826  REAL :: var(:,:,:)
827!-
828  INTEGER :: fid, ncvarid, ndim, iret
829  LOGICAL :: check = .FALSE.
830!---------------------------------------------------------------------
831  IF (check) WRITE(*,*) &
832     "flinput_r3d : SIZE(var) = ",SIZE(var)
833!-
834  CALL flinput_mat &
835    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
836     fid,ncvarid,ndim)
837!-
838  iret = NF90_PUT_VAR (fid, ncvarid, var, &
839           start=w_sta(1:ndim), count=w_len(1:ndim))
840!-------------------------
841END SUBROUTINE flinput_r3d
842!-
843!===
844!-
845SUBROUTINE flinput_r4d &
846  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
847!---------------------------------------------------------------------
848  IMPLICIT NONE
849!-
850  INTEGER :: fid_in
851  CHARACTER(LEN=*) :: varname
852  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
853  REAL :: var(:,:,:,:)
854!-
855  INTEGER :: fid, ncvarid, ndim, iret
856  LOGICAL :: check = .FALSE.
857!---------------------------------------------------------------------
858  IF (check) WRITE(*,*) &
859     "flinput_r4d : SIZE(var) = ",SIZE(var)
860!-
861  CALL flinput_mat &
862    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
863     fid,ncvarid,ndim)
864!-
865  iret = NF90_PUT_VAR (fid, ncvarid, var, &
866           start=w_sta(1:ndim), count=w_len(1:ndim))
867!-------------------------
868END SUBROUTINE flinput_r4d
869!-
870!===
871!-
872SUBROUTINE flinput_mat &
873  (fid_in,varname,iim,nlonid,jjm,nlatid, &
874                  llm,zdimid,ttm,tdimid,fid,ncvarid,ndim)
875!---------------------------------------------------------------------
876  IMPLICIT NONE
877!-
878  INTEGER :: fid_in
879  CHARACTER(LEN=*) :: varname
880  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
881  INTEGER :: fid, ncvarid, ndim
882!-
883! LOCAL
884!-
885  INTEGER :: iret
886!---------------------------------------------------------------------
887  fid = ncids(fid_in)
888!-
889  w_sta(1:4) = (/      1,      1,  1,  1 /)
890  w_len(1:2) = (/    iim,    jjm /)
891  w_dim(1:2) = (/ nlonid, nlatid /)
892!-
893  IF ( (llm > 0).AND.(ttm > 0) ) THEN
894    ndim = 4
895    w_len(3:4) = (/    llm,    ttm /)
896    w_dim(3:4) = (/ zdimid, tdimid /)
897  ELSE IF (llm > 0) THEN
898    ndim = 3
899    w_dim(3) = zdimid
900    w_len(3) = llm
901  ELSE IF (ttm > 0) THEN
902    ndim = 3
903    w_dim(3) = tdimid
904    w_len(3) = ttm
905  ELSE
906    ndim = 2
907  ENDIF
908!-
909  iret = NF90_REDEF   (fid)
910  iret = NF90_DEF_VAR (fid,varname,NF90_FLOAT,w_dim(1:ndim),ncvarid)
911  iret = NF90_PUT_ATT (fid,ncvarid,'short_name',TRIM(varname))
912  iret = NF90_ENDDEF  (fid)
913!--------------------------
914END  SUBROUTINE flinput_mat
915!-
916!===
917!-
918SUBROUTINE flinput_scal &
919  (fid_in, varname, iim, nlonid, jjm, nlatid, &
920                    llm, zdimid, ttm, tdimid, var)
921!---------------------------------------------------------------------
922  IMPLICIT NONE
923!-
924  INTEGER :: fid_in
925  CHARACTER(LEN=*) :: varname
926  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
927  REAL :: var
928!-
929! LOCAL
930!-
931  INTEGER :: fid, iret
932!---------------------------------------------------------------------
933  fid = ncids(fid_in)
934!-
935  iret = NF90_REDEF   (fid)
936  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, varname, REAL(var,KIND=4))
937  iret = NF90_ENDDEF  (fid)
938!---------------------------
939END  SUBROUTINE flinput_scal
940!-
941!===
942!-
943SUBROUTINE flinget_r1d &
944  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
945!---------------------------------------------------------------------
946  IMPLICIT NONE
947!-
948  INTEGER :: fid_in
949  CHARACTER(LEN=*) :: varname
950  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
951  REAL :: var(:)
952!-
953  INTEGER :: jl, ji
954  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
955  LOGICAL :: check = .FALSE.
956!---------------------------------------------------------------------
957  IF (.NOT.ALLOCATED(buff_tmp)) THEN
958    IF (check) WRITE(*,*) &
959      "flinget_r1d : allocate buff_tmp for buff_sz = ",SIZE(var)
960    ALLOCATE (buff_tmp(SIZE(var)))
961  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
962    IF (check) WRITE(*,*) &
963      "flinget_r1d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
964    DEALLOCATE (buff_tmp)
965    ALLOCATE (buff_tmp(SIZE(var)))
966  ENDIF
967!-
968  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
969                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
970!-
971  jl=0
972  DO ji=1,SIZE(var,1)
973    jl=jl+1
974    var(ji) = buff_tmp(jl)
975  ENDDO
976!-------------------------
977END SUBROUTINE flinget_r1d
978!-
979!===
980!-
981SUBROUTINE flinget_r2d &
982  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
983!---------------------------------------------------------------------
984  IMPLICIT NONE
985!-
986  INTEGER :: fid_in
987  CHARACTER(LEN=*) :: varname
988  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
989  REAL :: var(:,:)
990!-
991  INTEGER :: jl, jj, ji
992  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
993  LOGICAL :: check = .FALSE.
994!---------------------------------------------------------------------
995  IF (.NOT.ALLOCATED(buff_tmp)) THEN
996    IF (check) WRITE(*,*) &
997      "flinget_r2d : allocate buff_tmp for buff_sz = ",SIZE(var)
998    ALLOCATE (buff_tmp(SIZE(var)))
999  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1000    IF (check) WRITE(*,*) &
1001      "flinget_r2d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1002    DEALLOCATE (buff_tmp)
1003    ALLOCATE (buff_tmp(SIZE(var)))
1004  ENDIF
1005!-
1006  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1007                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1008!-
1009  jl=0
1010  DO jj=1,SIZE(var,2)
1011    DO ji=1,SIZE(var,1)
1012      jl=jl+1
1013      var(ji,jj) = buff_tmp(jl)
1014    ENDDO
1015  ENDDO
1016!-------------------------
1017END SUBROUTINE flinget_r2d
1018!-
1019!===
1020!-
1021SUBROUTINE flinget_r2d_zoom2d &
1022  (fid_in,varname,iim,jjm,llm,ttm, &
1023   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1024!---------------------------------------------------------------------
1025  IMPLICIT NONE
1026!-
1027  INTEGER :: fid_in
1028  CHARACTER(LEN=*) :: varname
1029  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1030  REAL :: var(:,:)
1031!-
1032  INTEGER :: jl, jj, ji
1033  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1034  LOGICAL :: check = .FALSE.
1035!---------------------------------------------------------------------
1036  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1037    IF (check) WRITE(*,*) &
1038      "flinget_r2d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1039    ALLOCATE (buff_tmp(SIZE(var)))
1040  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1041    IF (check) WRITE(*,*) &
1042      "flinget_r2d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1043    DEALLOCATE (buff_tmp)
1044    ALLOCATE (buff_tmp(SIZE(var)))
1045  ENDIF
1046!-
1047  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1048                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1049!-
1050  jl=0
1051  DO jj=1,SIZE(var,2)
1052    DO ji=1,SIZE(var,1)
1053      jl=jl+1
1054      var(ji,jj) = buff_tmp(jl)
1055    ENDDO
1056  ENDDO
1057!--------------------------------
1058END SUBROUTINE flinget_r2d_zoom2d
1059!-
1060!===
1061!-
1062SUBROUTINE flinget_r3d &
1063  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
1064!---------------------------------------------------------------------
1065  IMPLICIT NONE
1066!-
1067  INTEGER :: fid_in
1068  CHARACTER(LEN=*) :: varname
1069  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1070  REAL :: var(:,:,:)
1071!-
1072  INTEGER :: jl, jk, jj, ji
1073  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1074  LOGICAL :: check = .FALSE.
1075!---------------------------------------------------------------------
1076  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1077    IF (check) WRITE(*,*) &
1078      "flinget_r3d : allocate buff_tmp for buff_sz = ",SIZE(var)
1079    ALLOCATE (buff_tmp(SIZE(var)))
1080  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1081    IF (check) WRITE(*,*) &
1082      "flinget_r3d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1083    DEALLOCATE (buff_tmp)
1084    ALLOCATE (buff_tmp(SIZE(var)))
1085  ENDIF
1086!-
1087  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1088                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1089!-
1090  jl=0
1091  DO jk=1,SIZE(var,3)
1092    DO jj=1,SIZE(var,2)
1093      DO ji=1,SIZE(var,1)
1094        jl=jl+1
1095        var(ji,jj,jk) = buff_tmp(jl)
1096      ENDDO
1097    ENDDO
1098  ENDDO
1099!-------------------------
1100END SUBROUTINE flinget_r3d
1101!-
1102!===
1103!-
1104SUBROUTINE flinget_r3d_zoom2d &
1105  (fid_in,varname,iim,jjm,llm,ttm, &
1106   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1107!---------------------------------------------------------------------
1108  IMPLICIT NONE
1109!-
1110  INTEGER :: fid_in
1111  CHARACTER(LEN=*) :: varname
1112  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1113  REAL :: var(:,:,:)
1114!-
1115  INTEGER :: jl, jk, jj, ji
1116  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1117  LOGICAL :: check = .FALSE.
1118!---------------------------------------------------------------------
1119  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1120    IF (check) WRITE(*,*) &
1121      "flinget_r3d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1122    ALLOCATE (buff_tmp(SIZE(var)))
1123  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1124    IF (check) WRITE(*,*) &
1125      "flinget_r3d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1126    DEALLOCATE (buff_tmp)
1127    ALLOCATE (buff_tmp(SIZE(var)))
1128  ENDIF
1129!-
1130  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1131                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1132!-
1133  jl=0
1134  DO jk=1,SIZE(var,3)
1135    DO jj=1,SIZE(var,2)
1136      DO ji=1,SIZE(var,1)
1137        jl=jl+1
1138        var(ji,jj,jk) = buff_tmp(jl)
1139      ENDDO
1140    ENDDO
1141  ENDDO
1142!--------------------------------
1143END SUBROUTINE flinget_r3d_zoom2d
1144!-
1145!===
1146!-
1147SUBROUTINE flinget_r4d &
1148  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
1149!---------------------------------------------------------------------
1150  IMPLICIT NONE
1151!-
1152  INTEGER :: fid_in
1153  CHARACTER(LEN=*) :: varname
1154  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1155  REAL :: var(:,:,:,:)
1156!-
1157  INTEGER :: jl, jk, jj, ji, jm
1158  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1159  LOGICAL :: check = .FALSE.
1160!---------------------------------------------------------------------
1161  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1162    IF (check) WRITE(*,*) &
1163      "flinget_r4d : allocate buff_tmp for buff_sz = ",SIZE(var)
1164    ALLOCATE (buff_tmp(SIZE(var)))
1165  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1166    IF (check) WRITE(*,*) &
1167      "flinget_r4d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1168    DEALLOCATE (buff_tmp)
1169    ALLOCATE (buff_tmp(SIZE(var)))
1170  ENDIF
1171!-
1172  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1173                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1174!-
1175  jl=0
1176  DO jm=1,SIZE(var,4)
1177    DO jk=1,SIZE(var,3)
1178      DO jj=1,SIZE(var,2)
1179        DO ji=1,SIZE(var,1)
1180          jl=jl+1
1181          var(ji,jj,jk,jm) = buff_tmp(jl)
1182        ENDDO
1183      ENDDO
1184    ENDDO
1185  ENDDO
1186!-------------------------
1187END SUBROUTINE flinget_r4d
1188!-
1189!===
1190!-
1191SUBROUTINE flinget_r4d_zoom2d &
1192  (fid_in,varname,iim,jjm,llm,ttm, &
1193   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1194!---------------------------------------------------------------------
1195  IMPLICIT NONE
1196!-
1197  INTEGER :: fid_in
1198  CHARACTER(LEN=*) :: varname
1199  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1200  REAL :: var(:,:,:,:)
1201!-
1202  INTEGER :: jl, jk, jj, ji, jm
1203  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1204  LOGICAL :: check = .FALSE.
1205!---------------------------------------------------------------------
1206  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1207    IF (check) WRITE(*,*) &
1208      "flinget_r4d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1209    ALLOCATE (buff_tmp(SIZE(var)))
1210  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1211    IF (check) WRITE(*,*) &
1212      "flinget_r4d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1213    DEALLOCATE (buff_tmp)
1214    ALLOCATE (buff_tmp(SIZE(var)))
1215  ENDIF
1216!-
1217  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1218                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1219!-
1220  jl=0
1221  DO jm=1,SIZE(var,4)
1222    DO jk=1,SIZE(var,3)
1223      DO jj=1,SIZE(var,2)
1224        DO ji=1,SIZE(var,1)
1225          jl=jl+1
1226          var(ji,jj,jk,jm) = buff_tmp(jl)
1227        ENDDO
1228      ENDDO
1229    ENDDO
1230  ENDDO
1231!--------------------------------
1232END SUBROUTINE flinget_r4d_zoom2d
1233!-
1234!===
1235!-
1236SUBROUTINE flinget_mat &
1237  (fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
1238   itau_fin, iideb, iilen, jjdeb, jjlen, var)
1239!---------------------------------------------------------------------
1240!- This subroutine will read the variable named varname from
1241!- the file previously opened by flinopen and identified by fid
1242!-
1243!- It is checked that the dimensions of the variable to be read
1244!- correspond to what the user requested when he specified
1245!- iim, jjm and llm. The only exception which is allowed is
1246!- for compressed data where the horizontal grid is not expected
1247!- to be iim x jjm.
1248!-
1249!- If variable is of size zero a global attribute is read.
1250!- This global attribute will be of type real
1251!-
1252!- INPUT
1253!-
1254!- fid      : File ID returned by flinopen
1255!- varname  : Name of the variable to be read from the file
1256!- iim      : | These three variables give the size of the variables
1257!- jjm      : | to be read. It will be verified that the variables
1258!- llm      : | fits in there.
1259!- ttm      : |
1260!- itau_dep : Time step at which we will start to read
1261!- itau_fin : Time step until which we are going to read
1262!-            For the moment this is done on indexes
1263!-            but it should be in the physical space.
1264!-            If there is no time-axis in the file then use a
1265!-            itau_fin < itau_dep, this will tell flinget not to
1266!-            expect a time-axis in the file.
1267!- iideb    : index i for zoom
1268!- iilen    : length of zoom
1269!- jjdeb    : index j for zoom
1270!- jjlen    : length of zoom
1271!-
1272!- OUTPUT
1273!-
1274!- var      : array that will contain the data
1275!---------------------------------------------------------------------
1276  IMPLICIT NONE
1277!-
1278! ARGUMENTS
1279!-
1280  INTEGER :: fid_in
1281  CHARACTER(LEN=*) :: varname
1282  INTEGER :: iim, jjm, llm, ttm
1283  INTEGER :: itau_dep, itau_fin, iideb, iilen, jjdeb, jjlen
1284  REAL :: var(:)
1285!-
1286! LOCAL
1287!-
1288  INTEGER :: iret, fid
1289  INTEGER :: vid, cvid, clen
1290  CHARACTER(LEN=70) :: str1
1291  CHARACTER(LEN=250) :: att_n, tmp_n
1292  CHARACTER(LEN=5) :: axs_l
1293  INTEGER :: tmp_i
1294  REAL,SAVE :: mis_v=0.
1295  REAL :: tmp_r
1296  INTEGER :: ndims, x_typ, nb_atts
1297  INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: dimids
1298  INTEGER :: i, iv, nvars, i2d, cnd
1299  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: var_tmp
1300  LOGICAL :: uncompress = .FALSE.
1301  LOGICAL :: check = .FALSE.
1302!---------------------------------------------------------------------
1303  fid = ncids(fid_in)
1304!-
1305  IF (check) THEN
1306    WRITE(*,*) &
1307    'flinget_mat : fid_in, fid, varname :', fid_in, fid, TRIM(varname)
1308    WRITE(*,*) &
1309    'flinget_mat : iim, jjm, llm, ttm, itau_dep, itau_fin :', &
1310    iim, jjm, llm, ttm, itau_dep, itau_fin
1311    WRITE(*,*) &
1312    'flinget_mat : iideb, iilen, jjdeb, jjlen :', &
1313    iideb, iilen, jjdeb, jjlen
1314  ENDIF
1315!-
1316  uncompress = .FALSE.
1317!-
1318! 1.0 We get first all the details on this variable from the file
1319!-
1320  nvars = ncnbva(fid_in)
1321!-
1322  vid = -1
1323  iret = NF90_INQ_VARID (fid, varname, vid)
1324!-
1325  IF (vid < 0 .OR. iret /= NF90_NOERR) THEN
1326    CALL histerr (3,'flinget', &
1327      'Variable '//TRIM(varname)//' not found in file',' ',' ')
1328  ENDIF
1329!-
1330  iret = NF90_INQUIRE_VARIABLE (fid, vid, &
1331           ndims=ndims, dimids=dimids, nAtts=nb_atts)
1332  IF (check) THEN
1333    WRITE(*,*) &
1334    'flinget_mat : fid, vid :', fid, vid
1335    WRITE(*,*) &
1336    'flinget_mat : ndims, dimids(1:ndims), nb_atts :', &
1337    ndims, dimids(1:ndims), nb_atts
1338  ENDIF
1339!-
1340  w_dim(:) = 0
1341  DO i=1,ndims
1342    iret  = NF90_INQUIRE_DIMENSION (fid, dimids(i), len=w_dim(i))
1343  ENDDO
1344  IF (check) WRITE(*,*) &
1345    'flinget_mat : w_dim :', w_dim(1:ndims)
1346!-
1347  mis_v = 0.0; axs_l = ' ';
1348!-
1349  IF (nb_atts > 0) THEN
1350    IF (check) THEN
1351      WRITE(*,*) 'flinget_mat : attributes for variable :'
1352    ENDIF
1353  ENDIF
1354  DO i=1,nb_atts
1355    iret = NF90_INQ_ATTNAME (fid, vid, i, att_n)
1356    iret = NF90_INQUIRE_ATTRIBUTE (fid, vid, att_n, xtype=x_typ)
1357    CALL strlowercase (att_n)
1358    IF      (    (x_typ == NF90_INT).OR.(x_typ == NF90_SHORT) &
1359             .OR.(x_typ == NF90_BYTE) ) THEN
1360      iret = NF90_GET_ATT (fid, vid, att_n, tmp_i)
1361      IF (check) THEN
1362        WRITE(*,*) '   ',TRIM(att_n),' : ',tmp_i
1363      ENDIF
1364    ELSE IF ( (x_typ == NF90_FLOAT).OR.(x_typ == NF90_DOUBLE) ) THEN
1365      iret = NF90_GET_ATT (fid, vid, att_n, tmp_r)
1366      IF (check) THEN
1367        WRITE(*,*) '   ',TRIM(att_n),' : ',tmp_r
1368      ENDIF
1369      IF (index(att_n,'missing_value') > 0) THEN
1370        mis_v = tmp_r
1371      ENDIF
1372    ELSE
1373      tmp_n = ''
1374      iret = NF90_GET_ATT (fid, vid, att_n, tmp_n)
1375      IF (check) THEN
1376        WRITE(*,*) '   ',TRIM(att_n),' : ',TRIM(tmp_n)
1377      ENDIF
1378      IF (index(att_n,'axis') > 0) THEN
1379        axs_l = tmp_n
1380      ENDIF
1381    ENDIF
1382  ENDDO
1383!?
1384!!!!!!!!!! We will need a verification on the type of the variable
1385!?
1386!-
1387! 2.0 The dimensions are analysed to determine what is to be read
1388!-
1389! 2.1 the longitudes
1390!-
1391  IF ( w_dim(1) /= iim .OR. w_dim(2) /= jjm) THEN
1392!---
1393!-- There is a possibility that we have to deal with a compressed axis !
1394!---
1395    iret = NF90_INQUIRE_DIMENSION (fid, dimids(1), &
1396             name=tmp_n, len=clen)
1397    iret = NF90_INQ_VARID (fid, tmp_n, cvid)
1398!---
1399    IF (check) WRITE(*,*) &
1400      'Dimname, iret , NF90_NOERR : ',TRIM(tmp_n),iret,NF90_NOERR
1401!---
1402!-- If we have an axis which has the same name
1403!-- as the dimension we can see if it is compressed
1404!---
1405!-- TODO TODO for zoom2d
1406!---
1407    IF (iret == NF90_NOERR) THEN
1408      iret = NF90_GET_ATT (fid, cvid, 'compress', str1)
1409!-----
1410      IF (iret == NF90_NOERR) THEN
1411        iret = NF90_INQUIRE_VARIABLE (fid,cvid,xtype=x_typ,ndims=cnd)
1412!-------
1413        IF ( cnd /= 1 .AND. x_typ /= NF90_INT) THEN
1414          CALL histerr (3,'flinget', &
1415            'Variable '//TRIM(tmp_n)//' can not be a compressed axis', &
1416            'Either it has too many dimensions'// &
1417            ' or it is not of type integer', ' ')
1418        ELSE
1419!---------
1420!-------- Let us see if we already have that index table
1421!---------
1422          IF (    (cind_len /= clen).OR.(cind_vid /= cvid) &
1423              .OR.(cind_fid /= fid) ) THEN
1424            IF (ALLOCATED(cindex))   DEALLOCATE(cindex)
1425            ALLOCATE(cindex(clen))
1426            cind_len = clen
1427            cind_vid = cvid
1428            cind_fid = fid
1429            iret = NF90_GET_VAR (fid, cvid, cindex)
1430          ENDIF
1431!---------
1432!-------- In any case we need to set the slab of data to be read
1433!---------
1434          uncompress = .TRUE.
1435          w_sta(1) = 1
1436          w_len(1) = clen
1437          i2d = 1
1438        ENDIF
1439      ELSE
1440        str1 = 'The horizontal dimensions of '//varname
1441        CALL histerr (3,'flinget',str1, &
1442          'is not compressed and does not'// &
1443          ' correspond to the requested size',' ')
1444      ENDIF
1445    ELSE
1446      IF (w_dim(1) /= iim) THEN
1447        str1 = 'The longitude dimension of '//varname
1448        CALL histerr (3,'flinget',str1, &
1449          'in the file is not equal to the dimension', &
1450          'that should be read')
1451      ENDIF
1452      IF (w_dim(2) /= jjm) THEN
1453        str1 = 'The latitude dimension of '//varname
1454        CALL histerr (3,'flinget',str1, &
1455          'in the file is not equal to the dimension', &
1456          'that should be read')
1457      ENDIF
1458    ENDIF
1459  ELSE
1460    w_sta(1:2) = (/ iideb, jjdeb /)
1461    w_len(1:2) = (/ iilen, jjlen /)
1462    i2d = 2
1463  ENDIF
1464!-
1465! 2.3 Now the difficult part, the 3rd dimension which can be
1466! time or levels.
1467!-
1468! Priority is given to the time axis if only three axes are present.
1469!-
1470  IF (ndims > i2d) THEN
1471!---
1472!-- 2.3.1 We have a vertical axis
1473!---
1474    IF (llm == 1 .AND. ndims == i2d+2 .OR. llm == w_dim(i2d+1)) THEN
1475!-----
1476      IF (w_dim(i2d+1) /= llm) THEN
1477        CALL histerr (3,'flinget', &
1478          'The vertical dimension of '//varname, &
1479          'in the file is not equal to the dimension', &
1480          'that should be read')
1481      ELSE
1482        w_sta(i2d+1) = 1
1483        IF (llm > 0) THEN
1484          w_len(i2d+1) = llm
1485        ELSE
1486          w_len(i2d+1) = w_sta(i2d+1)
1487        ENDIF
1488      ENDIF
1489!-----
1490      IF ((itau_fin-itau_dep) >= 0) THEN
1491        IF      (ndims /= i2d+2) THEN
1492          CALL histerr (3,'flinget', &
1493            'You attempt to read a time slab', &
1494            'but there is no time axis on this variable', varname)
1495        ELSE IF ((itau_fin - itau_dep) <= w_dim(i2d+2)) THEN
1496          w_sta(i2d+2) = itau_dep
1497          w_len(i2d+2) = itau_fin-itau_dep+1
1498        ELSE
1499          CALL histerr (3,'flinget', &
1500            'The time step you try to read is not', &
1501            'in the file (1)', varname)
1502        ENDIF
1503      ELSE IF (ndims == i2d+2 .AND. w_dim(i2d+2) > 1) THEN
1504        CALL histerr (3,'flinget', &
1505          'There is a time axis in the file but no', &
1506          'time step give in the call', varname)
1507      ELSE
1508        w_sta(i2d+2) = 1
1509        w_len(i2d+2) = 1
1510      ENDIF
1511    ELSE
1512!-----
1513!---- 2.3.2 We do not have any vertical axis
1514!-----
1515      IF (ndims == i2d+2) THEN
1516        CALL histerr (3,'flinget', &
1517          'The file contains 4 dimensions', &
1518          'but only 3 are requestes for variable ', varname)
1519      ENDIF
1520      IF ((itau_fin-itau_dep) >= 0) THEN
1521        IF (ndims == i2d+1) THEN
1522          IF ((itau_fin-itau_dep) < w_dim(i2d+1) ) THEN
1523            w_sta(i2d+1) = itau_dep
1524            w_len(i2d+1) = itau_fin-itau_dep+1
1525          ELSE
1526            CALL histerr (3,'flinget', &
1527              'The time step you try to read is not', &
1528              'in the file (2)', varname)
1529          ENDIF
1530        ELSE
1531          CALL histerr (3,'flinget', &
1532            'From your input you sould have 3 dimensions', &
1533            'in the file but there are 4', varname)
1534        ENDIF
1535      ELSE
1536        IF (ndims == i2d+1 .AND. w_dim(i2d+1) > 1) THEN
1537          CALL histerr (3,'flinget', &
1538            'There is a time axis in the file but no', &
1539            'time step given in the call', varname)
1540        ELSE
1541          w_sta(i2d+1) = 1
1542          w_len(i2d+1) = 1
1543        ENDIF
1544      ENDIF
1545    ENDIF
1546  ELSE
1547!---
1548!-- 2.3.3 We do not have any vertical axis
1549!---
1550    w_sta(i2d+1:i2d+2) = (/ 0, 0 /)
1551    w_len(i2d+1:i2d+2) = (/ 0, 0 /)
1552  ENDIF
1553!-
1554! 3.0 Reading the data
1555!-
1556  IF (check) WRITE(*,*) &
1557    'flinget_mat 3.0 : ', uncompress, w_sta, w_len
1558!---
1559  IF (uncompress) THEN
1560!---
1561    IF (ALLOCATED(var_tmp)) THEN
1562      IF (SIZE(var_tmp) < clen) THEN
1563        DEALLOCATE(var_tmp)
1564        ALLOCATE(var_tmp(clen))
1565      ENDIF
1566    ELSE
1567      ALLOCATE(var_tmp(clen))
1568    ENDIF
1569!---
1570    iret = NF90_GET_VAR (fid, vid, var_tmp, &
1571             start=w_sta(:), count=w_len(:))
1572!---
1573    var(:) = mis_v
1574    var(cindex(:)) = var_tmp(:)
1575!---
1576  ELSE
1577    iret = NF90_GET_VAR (fid, vid, var, &
1578             start=w_sta(:), count=w_len(:))
1579  ENDIF
1580!-
1581  IF (check) WRITE(*,*) 'flinget_mat 3.1 : ',NF90_STRERROR (iret)
1582!--------------------------
1583END  SUBROUTINE flinget_mat
1584!-
1585!===
1586!-
1587SUBROUTINE flinget_scal &
1588  (fid_in, varname, iim, jjm, llm, ttm, itau_dep, itau_fin, var)
1589!---------------------------------------------------------------------
1590!- This subroutine will read the variable named varname from
1591!- the file previously opened by flinopen and identified by fid
1592!-
1593!- If variable is of size zero a global attribute is read. This
1594!- global attribute will be of type real
1595!-
1596!- INPUT
1597!-
1598!- fid      : File ID returned by flinopen
1599!- varname  : Name of the variable to be read from the file
1600!- iim      : | These three variables give the size of the variables
1601!- jjm      : | to be read. It will be verified that the variables
1602!- llm      : | fits in there.
1603!- ttm      : |
1604!- itau_dep : Time step at which we will start to read
1605!- itau_fin : Time step until which we are going to read
1606!-           For the moment this is done on indeces but it should be
1607!-           in the physical space
1608!-           If there is no time-axis in the file then use a
1609!-           itau_fin < itau_dep, this will tell flinget not to
1610!-           expect a time-axis in the file.
1611!-
1612!- OUTPUT
1613!-
1614!- var      : scalar that will contain the data
1615!---------------------------------------------------------------------
1616  IMPLICIT NONE
1617!-
1618! ARGUMENTS
1619!-
1620  INTEGER :: fid_in
1621  CHARACTER(LEN=*) :: varname
1622  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1623  REAL :: var
1624!-
1625! LOCAL
1626!-
1627  INTEGER :: iret, fid
1628!-
1629  LOGICAL :: check = .FALSE.
1630!---------------------------------------------------------------------
1631  fid = ncids(fid_in)
1632!-
1633! 1.0 Reading a global attribute
1634!-
1635  iret = NF90_GET_ATT (fid, NF90_GLOBAL, varname, var)
1636!---------------------------
1637END  SUBROUTINE flinget_scal
1638!-
1639!===
1640!-
1641SUBROUTINE flinfindcood (fid_in, axtype, vid, ndim)
1642!---------------------------------------------------------------------
1643!- This subroutine explores the file in order to find
1644!- the coordinate according to a number of rules
1645!---------------------------------------------------------------------
1646  IMPLICIT NONE
1647!-
1648! ARGUMENTS
1649!-
1650  INTEGER :: fid_in, vid, ndim
1651  CHARACTER(LEN=3) :: axtype
1652!-
1653! LOCAL
1654!-
1655  INTEGER :: iv, iret, dimnb
1656  CHARACTER(LEN=40) :: dimname, dimuni1, dimuni2, dimuni3
1657  CHARACTER(LEN=30) :: str1
1658  LOGICAL :: found_rule = .FALSE.
1659!---------------------------------------------------------------------
1660  vid = -1
1661!-
1662! Make sure all strings are invalid
1663!-
1664  dimname = '?-?'
1665  dimuni1 = '?-?'
1666  dimuni2 = '?-?'
1667  dimuni3 = '?-?'
1668!-
1669! First rule : we look for the correct units
1670! lon : east
1671! lat : north
1672! We make an exact check as it would be too easy to mistake
1673! some units by just comparing the substrings.
1674!-
1675  SELECTCASE(axtype)
1676  CASE ('lon')
1677    dimuni1 = 'degree_e'
1678    dimuni2 = 'degrees_e'
1679    found_rule = .TRUE.
1680  CASE('lat')
1681    dimuni1 = 'degree_n'
1682    dimuni2 = 'degrees_n'
1683    found_rule = .TRUE.
1684  CASE('lev')
1685    dimuni1 = 'm'
1686    dimuni2 = 'km'
1687    dimuni3 = 'hpa'
1688    found_rule = .TRUE.
1689  CASE DEFAULT
1690    found_rule = .FALSE.
1691  END SELECT
1692!-
1693  IF (found_rule) THEN
1694    iv = 0
1695    DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1696      iv = iv+1
1697      str1 = ''
1698      iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1)
1699      IF (iret == NF90_NOERR) THEN
1700        CALL strlowercase (str1)
1701        IF (    (INDEX(str1, TRIM(dimuni1)) == 1) &
1702            .OR.(INDEX(str1, TRIM(dimuni2)) == 1) &
1703            .OR.(INDEX(str1, TRIM(dimuni3)) == 1) ) THEN
1704          vid = iv
1705          iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim)
1706        ENDIF
1707      ENDIF
1708    ENDDO
1709  ENDIF
1710!-
1711! Second rule : we find specific names :
1712! lon : nav_lon
1713! lat : nav_lat
1714! Here we can check if we find the substring as the
1715! names are more specific.
1716!-
1717  SELECTCASE(axtype)
1718  CASE ('lon')
1719    dimname = 'nav_lon lon longitude'
1720    found_rule = .TRUE.
1721  CASE('lat')
1722    dimname = 'nav_lat lat latitude'
1723    found_rule = .TRUE.
1724  CASE('lev')
1725    dimname = 'plev level depth deptht'
1726    found_rule = .TRUE.
1727  CASE DEFAULT
1728    found_rule = .FALSE.
1729  END SELECT
1730!-
1731  IF (found_rule) THEN
1732    iv = 0
1733    DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1734      iv = iv+1
1735      str1=''
1736      iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
1737                                    name=str1, ndims=ndim)
1738      IF (INDEX(dimname,TRIM(str1)) >= 1) THEN
1739        vid = iv
1740      ENDIF
1741    ENDDO
1742  ENDIF
1743!-
1744! Third rule : we find a variable with the same name as the dimension
1745! lon = 1
1746! lat = 2
1747! lev = 3
1748!-
1749  IF (vid < 0) THEN
1750    SELECTCASE(axtype)
1751    CASE ('lon')
1752      dimnb = 1
1753      found_rule = .TRUE.
1754    CASE('lat')
1755      dimnb = 2
1756      found_rule = .TRUE.
1757    CASE('lev')
1758      dimnb = 3
1759      found_rule = .TRUE.
1760    CASE DEFAULT
1761      found_rule = .FALSE.
1762    END SELECT
1763!---
1764    IF (found_rule) THEN
1765      iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname)
1766      iv = 0
1767      DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1768        iv = iv+1
1769        str1=''
1770        iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
1771                                      name=str1, ndims=ndim)
1772        IF (INDEX(dimname,TRIM(str1)) == 1) THEN
1773          vid = iv
1774        ENDIF
1775      ENDDO
1776    ENDIF
1777  ENDIF
1778!-
1779! Stop the program if no coordinate was found
1780!-
1781  IF (vid < 0) THEN
1782    CALL histerr (3,'flinfindcood', &
1783           'No coordinate axis was found in the file', &
1784           'The data in this file can not be used', axtype)
1785  ENDIF
1786!--------------------------
1787END SUBROUTINE flinfindcood
1788!-
1789!===
1790!-
1791SUBROUTINE flinclo (fid_in)
1792!---------------------------------------------------------------------
1793  IMPLICIT NONE
1794!-
1795  INTEGER :: fid_in
1796!-
1797  INTEGER :: iret
1798!---------------------------------------------------------------------
1799  iret = NF90_CLOSE (ncids(fid_in))
1800  ncfileopen(fid_in) = .FALSE.
1801!---------------------
1802END SUBROUTINE flinclo
1803!-
1804!===
1805!-
1806SUBROUTINE flinquery_var(fid_in, varname, exists)
1807!---------------------------------------------------------------------
1808!- Queries the existance of a variable in the file.
1809!---------------------------------------------------------------------
1810  IMPLICIT NONE
1811!-
1812  INTEGER :: fid_in
1813  CHARACTER(LEN=*) varname
1814  LOGICAL :: exists
1815!-
1816  INTEGER :: iret, fid, vid
1817!---------------------------------------------------------------------
1818  fid = ncids(fid_in)
1819  vid = -1
1820  iret = NF90_INQ_VARID (fid, varname, vid)
1821!-
1822  exists = ( (vid >= 0).AND.(iret == NF90_NOERR) )
1823!---------------------------
1824END SUBROUTINE flinquery_var
1825!-
1826!===
1827!-
1828SUBROUTINE flininspect (fid_in)
1829!---------------------------------------------------------------------
1830  IMPLICIT NONE
1831!-
1832! fid : File id to inspect
1833!-
1834  INTEGER :: fid_in
1835!-
1836!- LOCAL
1837!-
1838  INTEGER :: iim, jjm, llm, ttm, fid_out
1839  INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim
1840  INTEGER :: iv, in, lll
1841  INTEGER :: xid, yid, zid, tid
1842  INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: idimid
1843  CHARACTER(LEN=80) :: name
1844  CHARACTER(LEN=30) :: axname
1845!---------------------------------------------------------------------
1846  fid = ncids(fid_in)
1847!-
1848  iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, &
1849                       nAttributes=nb_atts, unlimitedDimId=id_unlim)
1850!-
1851  WRITE (*,*) 'IOIPSL ID                   : ',fid_in
1852  WRITE (*,*) 'NetCDF ID                   : ',fid
1853  WRITE (*,*) 'Number of dimensions        : ',ndims
1854  WRITE (*,*) 'Number of variables         : ',nvars
1855  WRITE (*,*) 'Number of global attributes : ',nb_atts
1856  WRITE (*,*) 'ID unlimited                : ',id_unlim
1857!-
1858  xid = -1; iim = 0;
1859  yid = -1; jjm = 0;
1860  zid = -1; llm = 0;
1861  tid = -1; ttm = 0;
1862!-
1863  DO iv=1,ndims
1864!---
1865    iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll)
1866    CALL strlowercase (axname)
1867    axname = ADJUSTL(axname)
1868!---
1869    WRITE (*,*) 'Dimension number : ',iv
1870    WRITE (*,*) 'Dimension name   : ',TRIM(axname)
1871!---
1872    IF      (    (INDEX(axname,'x') == 1) &
1873             .OR.(INDEX(axname,'lon') == 1)) THEN
1874      xid = iv; iim = lll;
1875      WRITE (*,*) 'Dimension X size   : ',iim
1876    ELSE IF (    (INDEX(axname,'y') == 1) &
1877             .OR.(INDEX(axname,'lat') == 1)) THEN
1878      yid = iv; jjm = lll;
1879      WRITE (*,*) 'Dimension Y size   : ',jjm
1880    ELSE IF (    (INDEX(axname,'lev') == 1) &
1881             .OR.(INDEX(axname,'plev') == 1) &
1882             .OR.(INDEX(axname,'z') == 1) &
1883             .OR.(INDEX(axname,'depth') == 1)) THEN
1884      zid = iv; llm = lll;
1885      WRITE (*,*) 'Dimension Z size   : ',llm
1886    ELSE IF (    (INDEX(axname,'tstep') == 1) &
1887             .OR.(INDEX(axname,'time_counter') == 1)) THEN
1888!---- For the time we certainly need to allow for other names
1889      tid = iv; ttm = lll;
1890    ELSE IF (ndims == 1) THEN
1891!---- Nothing was found and ndims=1 then we have a vector of data
1892      xid = 1; iim = lll;
1893    ENDIF
1894!---
1895  ENDDO
1896!-
1897! Keep all this information
1898!-
1899  nbfiles = nbfiles+1
1900!-
1901  IF (nbfiles > nbfile_max) THEN
1902    CALL histerr(3,'flininspect', &
1903      'Too many files. Please increase nbfil_max', &
1904      'in program flincom.F90.',' ')
1905  ENDIF
1906!-
1907  ncids(nbfiles) = fid
1908  ncnbd(nbfiles) = ndims
1909!-
1910  ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /)
1911!-
1912  ncfunli(nbfiles) = id_unlim
1913  ncnba(nbfiles)   = nb_atts
1914  ncnbva(nbfiles)  = nvars
1915  ncfileopen(nbfiles) = .TRUE.
1916!-
1917  fid_out = nbfiles
1918!-
1919  DO in=1,nvars
1920    iret = NF90_INQUIRE_VARIABLE (fid, in, &
1921             name=name, ndims=ndims, dimids=idimid, nAtts=nb_atts)
1922    WRITE (*,*) 'Variable number  ------------ > ', in
1923    WRITE (*,*) 'Variable name        : ', TRIM(name)
1924    WRITE (*,*) 'Number of dimensions : ', ndims
1925    WRITE (*,*) 'Dimensions ID''s     : ', idimid(1:ndims)
1926    WRITE (*,*) 'Number of attributes : ', nb_atts
1927  ENDDO
1928!-------------------------
1929END SUBROUTINE flininspect
1930!-
1931!===
1932!-
1933END MODULE flincom
Note: See TracBrowser for help on using the repository browser.