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

Last change on this file since 368 was 368, checked in by bellier, 16 years ago

Suppressed an unused variable.

  • Property svn:keywords set to Id
File size: 57.0 KB
Line 
1!$Id$
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 :: 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
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,start=(/ 1 /),count=(/ ttm /))
602    itaus(1:ttm) = NINT(vec_tmp(1:ttm))
603    DEALLOCATE(vec_tmp)
604!---
605    IF (check) WRITE(*,*) 'flinopen 5.1 Times ',itaus
606!---
607!-- Getting all the details for the time axis
608!---
609!-- Find the calendar
610    iret = NF90_GET_ATT (fid,gdtmaf_id,'calendar',calendar)
611    IF (iret == NF90_NOERR) THEN
612      CALL ioconf_calendar(calendar)
613    ENDIF
614!--
615    units = ''
616    iret = NF90_GET_ATT (fid,vid,'units',units)
617    IF (gdtt_id > 0) THEN
618      units = units(INDEX(units,'since')+6:LEN_TRIM(units))
619      READ (units,'(I4.4,5(a,I2.2))') &
620        year0, strc, month0, strc, day0, &
621               strc, hours0, strc, minutes0, strc, seci
622      sec0 = hours0*3600. + minutes0*60. + seci
623      CALL ymds2ju (year0, month0, day0, sec0, date0)
624      IF (check) &
625        WRITE(*,*) 'flinopen 5.1 gdtt_id year0 ... date0 ', &
626                   year0, month0, day0, sec0, date0
627!-----
628      iret = NF90_GET_ATT (fid, gdtt_id, 'tstep_sec', dt)
629    ELSE IF (gdtmaf_id > 0) THEN
630      units = units(INDEX(units,'since')+6:LEN_TRIM(units))
631      READ (units,'(I4.4,5(a,I2.2))') &
632        year0, strc, month0, strc, day0, &
633               strc, hours0, strc, minutes0, strc, seci
634      sec0 = hours0*3600. + minutes0*60. + seci
635      CALL ymds2ju (year0, month0, day0, sec0, date0)
636!-----
637      IF (check) &
638        WRITE(*,*) 'flinopen 5.1 gdtmaf_id year0 ... date0 ', &
639                   year0, month0, day0, sec0, date0
640    ELSE IF (old_id > 0) THEN
641      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'delta_tstep_sec', dt)
642      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'day0', r_day)
643      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'sec0', sec)
644      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'year0', r_year)
645      iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'month0', r_month)
646!-----
647      day = INT(r_day)
648      month = INT(r_month)
649      year = INT(r_year)
650!-----
651      CALL ymds2ju (year, month, day, sec, date0)
652    ENDIF
653  ENDIF
654!-
655  IF (check) WRITE(*,*) 'flinopen 6.0 File opened', date0, dt
656!---------------------------
657END SUBROUTINE flinopen_work
658!-
659!===
660!-
661SUBROUTINE flininfo (filename, iim, jjm, llm, ttm, fid_out)
662!---------------------------------------------------------------------
663!- This subroutine allows to get some information.
664!- It is usualy done within flinopen but the user may want to call
665!- it before in order to allocate the space needed to extract the
666!- data from the file.
667!---------------------------------------------------------------------
668  IMPLICIT NONE
669!-
670! ARGUMENTS
671!-
672  CHARACTER(LEN=*) :: filename
673  INTEGER :: iim, jjm, llm, ttm, fid_out
674!-
675! LOCAL
676!-
677  INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim
678  INTEGER :: iv, lll
679  INTEGER :: xid, yid, zid, tid
680  CHARACTER(LEN=80) :: name
681  CHARACTER(LEN=30) :: axname
682!-
683  LOGICAL :: check = .FALSE.
684!---------------------------------------------------------------------
685  lll = LEN_TRIM(filename)
686  IF (filename(lll-2:lll) /= '.nc') THEN
687    name = filename(1:lll)//'.nc'
688  ELSE
689    name = filename(1:lll)
690  ENDIF
691!-
692  iret = NF90_OPEN (name, NF90_NOWRITE, fid)
693  IF (iret /= NF90_NOERR) THEN
694    CALL histerr(3, 'flininfo','Could not open file :',TRIM(name),' ')
695  ENDIF
696!-
697  iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, &
698                      nAttributes=nb_atts, unlimitedDimId=id_unlim)
699!-
700  xid = -1; iim = 0;
701  yid = -1; jjm = 0;
702  zid = -1; llm = 0;
703  tid = -1; ttm = 0;
704!-
705  DO iv=1,ndims
706!---
707    iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll)
708    CALL strlowercase (axname)
709    axname = ADJUSTL(axname)
710!---
711    IF (check) WRITE(*,*) &
712      'flininfo - getting axname',iv,axname,lll
713!---
714    IF      (    (INDEX(axname,'x') == 1) &
715             .OR.(INDEX(axname,'lon') == 1) ) THEN
716      xid = iv; iim = lll;
717    ELSE IF (    (INDEX(axname,'y') == 1) &
718             .OR.(INDEX(axname,'lat') == 1) ) THEN
719      yid = iv; jjm = lll;
720    ELSE IF (    (INDEX(axname,'lev') == 1) &
721             .OR.(INDEX(axname,'plev') == 1) &
722             .OR.(INDEX(axname,'z') == 1) &
723             .OR.(INDEX(axname,'depth') == 1) ) THEN
724      zid = iv; llm = lll;
725    ELSE IF (    (INDEX(axname,'tstep') == 1) &
726             .OR.(INDEX(axname,'time_counter') == 1) ) THEN
727!---- For the time we certainly need to allow for other names
728      tid = iv; ttm = lll;
729    ELSE IF (ndims == 1) THEN
730!---- Nothing was found and ndims=1 then we have a vector of data
731      xid = 1; iim = lll;
732    ENDIF
733!---
734  ENDDO
735!-
736! Keep all this information
737!-
738  nbfiles = nbfiles+1
739!-
740  IF (nbfiles > nbfile_max) THEN
741    CALL histerr (3,'flininfo', &
742      'Too many files. Please increase nbfil_max', &
743      'in program flincom.F90.',' ')
744  ENDIF
745!-
746  ncids(nbfiles) = fid
747  ncnbd(nbfiles) = ndims
748!-
749  ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /)
750!-
751  ncfunli(nbfiles) = id_unlim
752  ncnba(nbfiles)   = nb_atts
753  ncnbva(nbfiles)  = nvars
754  ncfileopen(nbfiles) = .TRUE.
755!-
756  fid_out = nbfiles
757!----------------------
758END SUBROUTINE flininfo
759!-
760!===
761!-
762SUBROUTINE flinput_r1d &
763  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
764!---------------------------------------------------------------------
765  IMPLICIT NONE
766!-
767  INTEGER :: fid_in
768  CHARACTER(LEN=*) :: varname
769  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
770  REAL :: var(:)
771!-
772  INTEGER :: fid, ncvarid, ndim, iret
773  LOGICAL :: check = .FALSE.
774!---------------------------------------------------------------------
775  IF (check) WRITE(*,*) &
776     "flinput_r1d : SIZE(var) = ",SIZE(var)
777!-
778  CALL flinput_mat &
779    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
780     fid,ncvarid,ndim)
781!-
782  iret = NF90_PUT_VAR (fid, ncvarid, var, &
783           start=w_sta(1:ndim), count=w_len(1:ndim))
784!-------------------------
785END SUBROUTINE flinput_r1d
786!-
787!===
788!-
789SUBROUTINE flinput_r2d &
790  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
791!---------------------------------------------------------------------
792  IMPLICIT NONE
793!-
794  INTEGER :: fid_in
795  CHARACTER(LEN=*) :: varname
796  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
797  REAL :: var(:,:)
798!-
799  INTEGER :: fid, ncvarid, ndim, iret
800  LOGICAL :: check = .FALSE.
801!---------------------------------------------------------------------
802  IF (check) WRITE(*,*) &
803     "flinput_r2d : SIZE(var) = ",SIZE(var)
804!-
805  CALL flinput_mat &
806    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
807     fid,ncvarid,ndim)
808!-
809  iret = NF90_PUT_VAR (fid, ncvarid, var, &
810           start=w_sta(1:ndim), count=w_len(1:ndim))
811!-------------------------
812END SUBROUTINE flinput_r2d
813!-
814!===
815!-
816SUBROUTINE flinput_r3d &
817  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
818!---------------------------------------------------------------------
819  IMPLICIT NONE
820!-
821  INTEGER :: fid_in
822  CHARACTER(LEN=*) :: varname
823  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
824  REAL :: var(:,:,:)
825!-
826  INTEGER :: fid, ncvarid, ndim, iret
827  LOGICAL :: check = .FALSE.
828!---------------------------------------------------------------------
829  IF (check) WRITE(*,*) &
830     "flinput_r3d : SIZE(var) = ",SIZE(var)
831!-
832  CALL flinput_mat &
833    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
834     fid,ncvarid,ndim)
835!-
836  iret = NF90_PUT_VAR (fid, ncvarid, var, &
837           start=w_sta(1:ndim), count=w_len(1:ndim))
838!-------------------------
839END SUBROUTINE flinput_r3d
840!-
841!===
842!-
843SUBROUTINE flinput_r4d &
844  (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid,var)
845!---------------------------------------------------------------------
846  IMPLICIT NONE
847!-
848  INTEGER :: fid_in
849  CHARACTER(LEN=*) :: varname
850  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
851  REAL :: var(:,:,:,:)
852!-
853  INTEGER :: fid, ncvarid, ndim, iret
854  LOGICAL :: check = .FALSE.
855!---------------------------------------------------------------------
856  IF (check) WRITE(*,*) &
857     "flinput_r4d : SIZE(var) = ",SIZE(var)
858!-
859  CALL flinput_mat &
860    (fid_in,varname,iim,nlonid,jjm,nlatid,llm,zdimid,ttm,tdimid, &
861     fid,ncvarid,ndim)
862!-
863  iret = NF90_PUT_VAR (fid, ncvarid, var, &
864           start=w_sta(1:ndim), count=w_len(1:ndim))
865!-------------------------
866END SUBROUTINE flinput_r4d
867!-
868!===
869!-
870SUBROUTINE flinput_mat &
871  (fid_in,varname,iim,nlonid,jjm,nlatid, &
872                  llm,zdimid,ttm,tdimid,fid,ncvarid,ndim)
873!---------------------------------------------------------------------
874  IMPLICIT NONE
875!-
876  INTEGER :: fid_in
877  CHARACTER(LEN=*) :: varname
878  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
879  INTEGER :: fid, ncvarid, ndim
880!-
881! LOCAL
882!-
883  INTEGER :: iret
884!---------------------------------------------------------------------
885  fid = ncids(fid_in)
886!-
887  w_sta(1:4) = (/      1,      1,  1,  1 /)
888  w_len(1:2) = (/    iim,    jjm /)
889  w_dim(1:2) = (/ nlonid, nlatid /)
890!-
891  IF ( (llm > 0).AND.(ttm > 0) ) THEN
892    ndim = 4
893    w_len(3:4) = (/    llm,    ttm /)
894    w_dim(3:4) = (/ zdimid, tdimid /)
895  ELSE IF (llm > 0) THEN
896    ndim = 3
897    w_dim(3) = zdimid
898    w_len(3) = llm
899  ELSE IF (ttm > 0) THEN
900    ndim = 3
901    w_dim(3) = tdimid
902    w_len(3) = ttm
903  ELSE
904    ndim = 2
905  ENDIF
906!-
907  iret = NF90_REDEF   (fid)
908  iret = NF90_DEF_VAR (fid,varname,NF90_FLOAT,w_dim(1:ndim),ncvarid)
909  iret = NF90_PUT_ATT (fid,ncvarid,'short_name',TRIM(varname))
910  iret = NF90_ENDDEF  (fid)
911!--------------------------
912END  SUBROUTINE flinput_mat
913!-
914!===
915!-
916SUBROUTINE flinput_scal &
917  (fid_in, varname, iim, nlonid, jjm, nlatid, &
918                    llm, zdimid, ttm, tdimid, var)
919!---------------------------------------------------------------------
920  IMPLICIT NONE
921!-
922  INTEGER :: fid_in
923  CHARACTER(LEN=*) :: varname
924  INTEGER :: iim, nlonid, jjm, nlatid, llm, zdimid, ttm, tdimid
925  REAL :: var
926!-
927! LOCAL
928!-
929  INTEGER :: fid, iret
930!---------------------------------------------------------------------
931  fid = ncids(fid_in)
932!-
933  iret = NF90_REDEF   (fid)
934  iret = NF90_PUT_ATT (fid, NF90_GLOBAL, varname, REAL(var,KIND=4))
935  iret = NF90_ENDDEF  (fid)
936!---------------------------
937END  SUBROUTINE flinput_scal
938!-
939!===
940!-
941SUBROUTINE flinget_r1d &
942  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
943!---------------------------------------------------------------------
944  IMPLICIT NONE
945!-
946  INTEGER :: fid_in
947  CHARACTER(LEN=*) :: varname
948  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
949  REAL :: var(:)
950!-
951  INTEGER :: jl, ji
952  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
953  LOGICAL :: check = .FALSE.
954!---------------------------------------------------------------------
955  IF (.NOT.ALLOCATED(buff_tmp)) THEN
956    IF (check) WRITE(*,*) &
957      "flinget_r1d : allocate buff_tmp for buff_sz = ",SIZE(var)
958    ALLOCATE (buff_tmp(SIZE(var)))
959  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
960    IF (check) WRITE(*,*) &
961      "flinget_r1d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
962    DEALLOCATE (buff_tmp)
963    ALLOCATE (buff_tmp(SIZE(var)))
964  ENDIF
965!-
966  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
967                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
968!-
969  jl=0
970  DO ji=1,SIZE(var,1)
971    jl=jl+1
972    var(ji) = buff_tmp(jl)
973  ENDDO
974!-------------------------
975END SUBROUTINE flinget_r1d
976!-
977!===
978!-
979SUBROUTINE flinget_r2d &
980  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
981!---------------------------------------------------------------------
982  IMPLICIT NONE
983!-
984  INTEGER :: fid_in
985  CHARACTER(LEN=*) :: varname
986  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
987  REAL :: var(:,:)
988!-
989  INTEGER :: jl, jj, ji
990  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
991  LOGICAL :: check = .FALSE.
992!---------------------------------------------------------------------
993  IF (.NOT.ALLOCATED(buff_tmp)) THEN
994    IF (check) WRITE(*,*) &
995      "flinget_r2d : allocate buff_tmp for buff_sz = ",SIZE(var)
996    ALLOCATE (buff_tmp(SIZE(var)))
997  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
998    IF (check) WRITE(*,*) &
999      "flinget_r2d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1000    DEALLOCATE (buff_tmp)
1001    ALLOCATE (buff_tmp(SIZE(var)))
1002  ENDIF
1003!-
1004  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1005                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1006!-
1007  jl=0
1008  DO jj=1,SIZE(var,2)
1009    DO ji=1,SIZE(var,1)
1010      jl=jl+1
1011      var(ji,jj) = buff_tmp(jl)
1012    ENDDO
1013  ENDDO
1014!-------------------------
1015END SUBROUTINE flinget_r2d
1016!-
1017!===
1018!-
1019SUBROUTINE flinget_r2d_zoom2d &
1020  (fid_in,varname,iim,jjm,llm,ttm, &
1021   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1022!---------------------------------------------------------------------
1023  IMPLICIT NONE
1024!-
1025  INTEGER :: fid_in
1026  CHARACTER(LEN=*) :: varname
1027  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1028  REAL :: var(:,:)
1029!-
1030  INTEGER :: jl, jj, ji
1031  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1032  LOGICAL :: check = .FALSE.
1033!---------------------------------------------------------------------
1034  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1035    IF (check) WRITE(*,*) &
1036      "flinget_r2d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1037    ALLOCATE (buff_tmp(SIZE(var)))
1038  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1039    IF (check) WRITE(*,*) &
1040      "flinget_r2d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1041    DEALLOCATE (buff_tmp)
1042    ALLOCATE (buff_tmp(SIZE(var)))
1043  ENDIF
1044!-
1045  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1046                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1047!-
1048  jl=0
1049  DO jj=1,SIZE(var,2)
1050    DO ji=1,SIZE(var,1)
1051      jl=jl+1
1052      var(ji,jj) = buff_tmp(jl)
1053    ENDDO
1054  ENDDO
1055!--------------------------------
1056END SUBROUTINE flinget_r2d_zoom2d
1057!-
1058!===
1059!-
1060SUBROUTINE flinget_r3d &
1061  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
1062!---------------------------------------------------------------------
1063  IMPLICIT NONE
1064!-
1065  INTEGER :: fid_in
1066  CHARACTER(LEN=*) :: varname
1067  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1068  REAL :: var(:,:,:)
1069!-
1070  INTEGER :: jl, jk, jj, ji
1071  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1072  LOGICAL :: check = .FALSE.
1073!---------------------------------------------------------------------
1074  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1075    IF (check) WRITE(*,*) &
1076      "flinget_r3d : allocate buff_tmp for buff_sz = ",SIZE(var)
1077    ALLOCATE (buff_tmp(SIZE(var)))
1078  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1079    IF (check) WRITE(*,*) &
1080      "flinget_r3d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1081    DEALLOCATE (buff_tmp)
1082    ALLOCATE (buff_tmp(SIZE(var)))
1083  ENDIF
1084!-
1085  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1086                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1087!-
1088  jl=0
1089  DO jk=1,SIZE(var,3)
1090    DO jj=1,SIZE(var,2)
1091      DO ji=1,SIZE(var,1)
1092        jl=jl+1
1093        var(ji,jj,jk) = buff_tmp(jl)
1094      ENDDO
1095    ENDDO
1096  ENDDO
1097!-------------------------
1098END SUBROUTINE flinget_r3d
1099!-
1100!===
1101!-
1102SUBROUTINE flinget_r3d_zoom2d &
1103  (fid_in,varname,iim,jjm,llm,ttm, &
1104   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1105!---------------------------------------------------------------------
1106  IMPLICIT NONE
1107!-
1108  INTEGER :: fid_in
1109  CHARACTER(LEN=*) :: varname
1110  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1111  REAL :: var(:,:,:)
1112!-
1113  INTEGER :: jl, jk, jj, ji
1114  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1115  LOGICAL :: check = .FALSE.
1116!---------------------------------------------------------------------
1117  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1118    IF (check) WRITE(*,*) &
1119      "flinget_r3d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1120    ALLOCATE (buff_tmp(SIZE(var)))
1121  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1122    IF (check) WRITE(*,*) &
1123      "flinget_r3d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1124    DEALLOCATE (buff_tmp)
1125    ALLOCATE (buff_tmp(SIZE(var)))
1126  ENDIF
1127!-
1128  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1129                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1130!-
1131  jl=0
1132  DO jk=1,SIZE(var,3)
1133    DO jj=1,SIZE(var,2)
1134      DO ji=1,SIZE(var,1)
1135        jl=jl+1
1136        var(ji,jj,jk) = buff_tmp(jl)
1137      ENDDO
1138    ENDDO
1139  ENDDO
1140!--------------------------------
1141END SUBROUTINE flinget_r3d_zoom2d
1142!-
1143!===
1144!-
1145SUBROUTINE flinget_r4d &
1146  (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var)
1147!---------------------------------------------------------------------
1148  IMPLICIT NONE
1149!-
1150  INTEGER :: fid_in
1151  CHARACTER(LEN=*) :: varname
1152  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1153  REAL :: var(:,:,:,:)
1154!-
1155  INTEGER :: jl, jk, jj, ji, jm
1156  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1157  LOGICAL :: check = .FALSE.
1158!---------------------------------------------------------------------
1159  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1160    IF (check) WRITE(*,*) &
1161      "flinget_r4d : allocate buff_tmp for buff_sz = ",SIZE(var)
1162    ALLOCATE (buff_tmp(SIZE(var)))
1163  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1164    IF (check) WRITE(*,*) &
1165      "flinget_r4d : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1166    DEALLOCATE (buff_tmp)
1167    ALLOCATE (buff_tmp(SIZE(var)))
1168  ENDIF
1169!-
1170  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1171                    itau_dep,itau_fin,1,iim,1,jjm,buff_tmp)
1172!-
1173  jl=0
1174  DO jm=1,SIZE(var,4)
1175    DO jk=1,SIZE(var,3)
1176      DO jj=1,SIZE(var,2)
1177        DO ji=1,SIZE(var,1)
1178          jl=jl+1
1179          var(ji,jj,jk,jm) = buff_tmp(jl)
1180        ENDDO
1181      ENDDO
1182    ENDDO
1183  ENDDO
1184!-------------------------
1185END SUBROUTINE flinget_r4d
1186!-
1187!===
1188!-
1189SUBROUTINE flinget_r4d_zoom2d &
1190  (fid_in,varname,iim,jjm,llm,ttm, &
1191   itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var)
1192!---------------------------------------------------------------------
1193  IMPLICIT NONE
1194!-
1195  INTEGER :: fid_in
1196  CHARACTER(LEN=*) :: varname
1197  INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen
1198  REAL :: var(:,:,:,:)
1199!-
1200  INTEGER :: jl, jk, jj, ji, jm
1201  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp
1202  LOGICAL :: check = .FALSE.
1203!---------------------------------------------------------------------
1204  IF (.NOT.ALLOCATED(buff_tmp)) THEN
1205    IF (check) WRITE(*,*) &
1206      "flinget_r4d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var)
1207    ALLOCATE (buff_tmp(SIZE(var)))
1208  ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN
1209    IF (check) WRITE(*,*) &
1210      "flinget_r4d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var)
1211    DEALLOCATE (buff_tmp)
1212    ALLOCATE (buff_tmp(SIZE(var)))
1213  ENDIF
1214!-
1215  CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, &
1216                    itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp)
1217!-
1218  jl=0
1219  DO jm=1,SIZE(var,4)
1220    DO jk=1,SIZE(var,3)
1221      DO jj=1,SIZE(var,2)
1222        DO ji=1,SIZE(var,1)
1223          jl=jl+1
1224          var(ji,jj,jk,jm) = buff_tmp(jl)
1225        ENDDO
1226      ENDDO
1227    ENDDO
1228  ENDDO
1229!--------------------------------
1230END SUBROUTINE flinget_r4d_zoom2d
1231!-
1232!===
1233!-
1234SUBROUTINE flinget_mat &
1235  (fid_in, varname, iim, jjm, llm, ttm, itau_dep, &
1236   itau_fin, iideb, iilen, jjdeb, jjlen, var)
1237!---------------------------------------------------------------------
1238!- This subroutine will read the variable named varname from
1239!- the file previously opened by flinopen and identified by fid
1240!-
1241!- It is checked that the dimensions of the variable to be read
1242!- correspond to what the user requested when he specified
1243!- iim, jjm and llm. The only exception which is allowed is
1244!- for compressed data where the horizontal grid is not expected
1245!- to be iim x jjm.
1246!-
1247!- If variable is of size zero a global attribute is read.
1248!- This global attribute will be of type real
1249!-
1250!- INPUT
1251!-
1252!- fid      : File ID returned by flinopen
1253!- varname  : Name of the variable to be read from the file
1254!- iim      : | These three variables give the size of the variables
1255!- jjm      : | to be read. It will be verified that the variables
1256!- llm      : | fits in there.
1257!- ttm      : |
1258!- itau_dep : Time step at which we will start to read
1259!- itau_fin : Time step until which we are going to read
1260!-            For the moment this is done on indexes
1261!-            but it should be in the physical space.
1262!-            If there is no time-axis in the file then use a
1263!-            itau_fin < itau_dep, this will tell flinget not to
1264!-            expect a time-axis in the file.
1265!- iideb    : index i for zoom
1266!- iilen    : length of zoom
1267!- jjdeb    : index j for zoom
1268!- jjlen    : length of zoom
1269!-
1270!- OUTPUT
1271!-
1272!- var      : array that will contain the data
1273!---------------------------------------------------------------------
1274  IMPLICIT NONE
1275!-
1276! ARGUMENTS
1277!-
1278  INTEGER :: fid_in
1279  CHARACTER(LEN=*) :: varname
1280  INTEGER :: iim, jjm, llm, ttm
1281  INTEGER :: itau_dep, itau_fin, iideb, iilen, jjdeb, jjlen
1282  REAL :: var(:)
1283!-
1284! LOCAL
1285!-
1286  INTEGER :: iret, fid
1287  INTEGER :: vid, cvid, clen
1288  CHARACTER(LEN=70) :: str1
1289  CHARACTER(LEN=250) :: att_n, tmp_n
1290  CHARACTER(LEN=5) :: axs_l
1291  INTEGER :: tmp_i
1292  REAL,SAVE :: mis_v=0.
1293  REAL :: tmp_r
1294  INTEGER :: ndims, x_typ, nb_atts
1295  INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: dimids
1296  INTEGER :: i, nvars, i2d, cnd
1297  REAL,DIMENSION(:),ALLOCATABLE,SAVE :: var_tmp
1298  LOGICAL :: uncompress = .FALSE.
1299  LOGICAL :: check = .FALSE.
1300!---------------------------------------------------------------------
1301  fid = ncids(fid_in)
1302!-
1303  IF (check) THEN
1304    WRITE(*,*) &
1305    'flinget_mat : fid_in, fid, varname :', fid_in, fid, TRIM(varname)
1306    WRITE(*,*) &
1307    'flinget_mat : iim, jjm, llm, ttm, itau_dep, itau_fin :', &
1308    iim, jjm, llm, ttm, itau_dep, itau_fin
1309    WRITE(*,*) &
1310    'flinget_mat : iideb, iilen, jjdeb, jjlen :', &
1311    iideb, iilen, jjdeb, jjlen
1312  ENDIF
1313!-
1314  uncompress = .FALSE.
1315!-
1316! 1.0 We get first all the details on this variable from the file
1317!-
1318  nvars = ncnbva(fid_in)
1319!-
1320  vid = -1
1321  iret = NF90_INQ_VARID (fid, varname, vid)
1322!-
1323  IF (vid < 0 .OR. iret /= NF90_NOERR) THEN
1324    CALL histerr (3,'flinget', &
1325      'Variable '//TRIM(varname)//' not found in file',' ',' ')
1326  ENDIF
1327!-
1328  iret = NF90_INQUIRE_VARIABLE (fid, vid, &
1329           ndims=ndims, dimids=dimids, nAtts=nb_atts)
1330  IF (check) THEN
1331    WRITE(*,*) &
1332    'flinget_mat : fid, vid :', fid, vid
1333    WRITE(*,*) &
1334    'flinget_mat : ndims, dimids(1:ndims), nb_atts :', &
1335    ndims, dimids(1:ndims), nb_atts
1336  ENDIF
1337!-
1338  w_dim(:) = 0
1339  DO i=1,ndims
1340    iret  = NF90_INQUIRE_DIMENSION (fid, dimids(i), len=w_dim(i))
1341  ENDDO
1342  IF (check) WRITE(*,*) &
1343    'flinget_mat : w_dim :', w_dim(1:ndims)
1344!-
1345  mis_v = 0.0; axs_l = ' ';
1346!-
1347  IF (nb_atts > 0) THEN
1348    IF (check) THEN
1349      WRITE(*,*) 'flinget_mat : attributes for variable :'
1350    ENDIF
1351  ENDIF
1352  DO i=1,nb_atts
1353    iret = NF90_INQ_ATTNAME (fid, vid, i, att_n)
1354    iret = NF90_INQUIRE_ATTRIBUTE (fid, vid, att_n, xtype=x_typ)
1355    CALL strlowercase (att_n)
1356    IF      (    (x_typ == NF90_INT).OR.(x_typ == NF90_SHORT) &
1357             .OR.(x_typ == NF90_BYTE) ) THEN
1358      iret = NF90_GET_ATT (fid, vid, att_n, tmp_i)
1359      IF (check) THEN
1360        WRITE(*,*) '   ',TRIM(att_n),' : ',tmp_i
1361      ENDIF
1362    ELSE IF ( (x_typ == NF90_FLOAT).OR.(x_typ == NF90_DOUBLE) ) THEN
1363      iret = NF90_GET_ATT (fid, vid, att_n, tmp_r)
1364      IF (check) THEN
1365        WRITE(*,*) '   ',TRIM(att_n),' : ',tmp_r
1366      ENDIF
1367      IF (index(att_n,'missing_value') > 0) THEN
1368        mis_v = tmp_r
1369      ENDIF
1370    ELSE
1371      tmp_n = ''
1372      iret = NF90_GET_ATT (fid, vid, att_n, tmp_n)
1373      IF (check) THEN
1374        WRITE(*,*) '   ',TRIM(att_n),' : ',TRIM(tmp_n)
1375      ENDIF
1376      IF (index(att_n,'axis') > 0) THEN
1377        axs_l = tmp_n
1378      ENDIF
1379    ENDIF
1380  ENDDO
1381!?
1382!!!!!!!!!! We will need a verification on the type of the variable
1383!?
1384!-
1385! 2.0 The dimensions are analysed to determine what is to be read
1386!-
1387! 2.1 the longitudes
1388!-
1389  IF ( w_dim(1) /= iim .OR. w_dim(2) /= jjm) THEN
1390!---
1391!-- There is a possibility that we have to deal with a compressed axis !
1392!---
1393    iret = NF90_INQUIRE_DIMENSION (fid, dimids(1), &
1394             name=tmp_n, len=clen)
1395    iret = NF90_INQ_VARID (fid, tmp_n, cvid)
1396!---
1397    IF (check) WRITE(*,*) &
1398      'Dimname, iret , NF90_NOERR : ',TRIM(tmp_n),iret,NF90_NOERR
1399!---
1400!-- If we have an axis which has the same name
1401!-- as the dimension we can see if it is compressed
1402!---
1403!-- TODO TODO for zoom2d
1404!---
1405    IF (iret == NF90_NOERR) THEN
1406      iret = NF90_GET_ATT (fid, cvid, 'compress', str1)
1407!-----
1408      IF (iret == NF90_NOERR) THEN
1409        iret = NF90_INQUIRE_VARIABLE (fid,cvid,xtype=x_typ,ndims=cnd)
1410!-------
1411        IF ( cnd /= 1 .AND. x_typ /= NF90_INT) THEN
1412          CALL histerr (3,'flinget', &
1413            'Variable '//TRIM(tmp_n)//' can not be a compressed axis', &
1414            'Either it has too many dimensions'// &
1415            ' or it is not of type integer', ' ')
1416        ELSE
1417!---------
1418!-------- Let us see if we already have that index table
1419!---------
1420          IF (    (cind_len /= clen).OR.(cind_vid /= cvid) &
1421              .OR.(cind_fid /= fid) ) THEN
1422            IF (ALLOCATED(cindex))   DEALLOCATE(cindex)
1423            ALLOCATE(cindex(clen))
1424            cind_len = clen
1425            cind_vid = cvid
1426            cind_fid = fid
1427            iret = NF90_GET_VAR (fid, cvid, cindex)
1428          ENDIF
1429!---------
1430!-------- In any case we need to set the slab of data to be read
1431!---------
1432          uncompress = .TRUE.
1433          w_sta(1) = 1
1434          w_len(1) = clen
1435          i2d = 1
1436        ENDIF
1437      ELSE
1438        str1 = 'The horizontal dimensions of '//varname
1439        CALL histerr (3,'flinget',str1, &
1440          'is not compressed and does not'// &
1441          ' correspond to the requested size',' ')
1442      ENDIF
1443    ELSE
1444      IF (w_dim(1) /= iim) THEN
1445        str1 = 'The longitude dimension of '//varname
1446        CALL histerr (3,'flinget',str1, &
1447          'in the file is not equal to the dimension', &
1448          'that should be read')
1449      ENDIF
1450      IF (w_dim(2) /= jjm) THEN
1451        str1 = 'The latitude dimension of '//varname
1452        CALL histerr (3,'flinget',str1, &
1453          'in the file is not equal to the dimension', &
1454          'that should be read')
1455      ENDIF
1456    ENDIF
1457  ELSE
1458    w_sta(1:2) = (/ iideb, jjdeb /)
1459    w_len(1:2) = (/ iilen, jjlen /)
1460    i2d = 2
1461  ENDIF
1462!-
1463! 2.3 Now the difficult part, the 3rd dimension which can be
1464! time or levels.
1465!-
1466! Priority is given to the time axis if only three axes are present.
1467!-
1468  IF (ndims > i2d) THEN
1469!---
1470!-- 2.3.1 We have a vertical axis
1471!---
1472    IF (llm == 1 .AND. ndims == i2d+2 .OR. llm == w_dim(i2d+1)) THEN
1473!-----
1474      IF (w_dim(i2d+1) /= llm) THEN
1475        CALL histerr (3,'flinget', &
1476          'The vertical dimension of '//varname, &
1477          'in the file is not equal to the dimension', &
1478          'that should be read')
1479      ELSE
1480        w_sta(i2d+1) = 1
1481        IF (llm > 0) THEN
1482          w_len(i2d+1) = llm
1483        ELSE
1484          w_len(i2d+1) = w_sta(i2d+1)
1485        ENDIF
1486      ENDIF
1487!-----
1488      IF ((itau_fin-itau_dep) >= 0) THEN
1489        IF      (ndims /= i2d+2) THEN
1490          CALL histerr (3,'flinget', &
1491            'You attempt to read a time slab', &
1492            'but there is no time axis on this variable', varname)
1493        ELSE IF ((itau_fin - itau_dep) <= w_dim(i2d+2)) THEN
1494          w_sta(i2d+2) = itau_dep
1495          w_len(i2d+2) = itau_fin-itau_dep+1
1496        ELSE
1497          CALL histerr (3,'flinget', &
1498            'The time step you try to read is not', &
1499            'in the file (1)', varname)
1500        ENDIF
1501      ELSE IF (ndims == i2d+2 .AND. w_dim(i2d+2) > 1) THEN
1502        CALL histerr (3,'flinget', &
1503          'There is a time axis in the file but no', &
1504          'time step give in the call', varname)
1505      ELSE
1506        w_sta(i2d+2) = 1
1507        w_len(i2d+2) = 1
1508      ENDIF
1509    ELSE
1510!-----
1511!---- 2.3.2 We do not have any vertical axis
1512!-----
1513      IF (ndims == i2d+2) THEN
1514        CALL histerr (3,'flinget', &
1515          'The file contains 4 dimensions', &
1516          'but only 3 are requestes for variable ', varname)
1517      ENDIF
1518      IF ((itau_fin-itau_dep) >= 0) THEN
1519        IF (ndims == i2d+1) THEN
1520          IF ((itau_fin-itau_dep) < w_dim(i2d+1) ) THEN
1521            w_sta(i2d+1) = itau_dep
1522            w_len(i2d+1) = itau_fin-itau_dep+1
1523          ELSE
1524            CALL histerr (3,'flinget', &
1525              'The time step you try to read is not', &
1526              'in the file (2)', varname)
1527          ENDIF
1528        ELSE
1529          CALL histerr (3,'flinget', &
1530            'From your input you sould have 3 dimensions', &
1531            'in the file but there are 4', varname)
1532        ENDIF
1533      ELSE
1534        IF (ndims == i2d+1 .AND. w_dim(i2d+1) > 1) THEN
1535          CALL histerr (3,'flinget', &
1536            'There is a time axis in the file but no', &
1537            'time step given in the call', varname)
1538        ELSE
1539          w_sta(i2d+1) = 1
1540          w_len(i2d+1) = 1
1541        ENDIF
1542      ENDIF
1543    ENDIF
1544  ELSE
1545!---
1546!-- 2.3.3 We do not have any vertical axis
1547!---
1548    w_sta(i2d+1:i2d+2) = (/ 0, 0 /)
1549    w_len(i2d+1:i2d+2) = (/ 0, 0 /)
1550  ENDIF
1551!-
1552! 3.0 Reading the data
1553!-
1554  IF (check) WRITE(*,*) &
1555    'flinget_mat 3.0 : ', uncompress, w_sta, w_len
1556!---
1557  IF (uncompress) THEN
1558!---
1559    IF (ALLOCATED(var_tmp)) THEN
1560      IF (SIZE(var_tmp) < clen) THEN
1561        DEALLOCATE(var_tmp)
1562        ALLOCATE(var_tmp(clen))
1563      ENDIF
1564    ELSE
1565      ALLOCATE(var_tmp(clen))
1566    ENDIF
1567!---
1568    iret = NF90_GET_VAR (fid, vid, var_tmp, &
1569             start=w_sta(:), count=w_len(:))
1570!---
1571    var(:) = mis_v
1572    var(cindex(:)) = var_tmp(:)
1573!---
1574  ELSE
1575    iret = NF90_GET_VAR (fid, vid, var, &
1576             start=w_sta(:), count=w_len(:))
1577  ENDIF
1578!-
1579  IF (check) WRITE(*,*) 'flinget_mat 3.1 : ',NF90_STRERROR (iret)
1580!--------------------------
1581END  SUBROUTINE flinget_mat
1582!-
1583!===
1584!-
1585SUBROUTINE flinget_scal &
1586  (fid_in, varname, iim, jjm, llm, ttm, itau_dep, itau_fin, var)
1587!---------------------------------------------------------------------
1588!- This subroutine will read the variable named varname from
1589!- the file previously opened by flinopen and identified by fid
1590!-
1591!- If variable is of size zero a global attribute is read. This
1592!- global attribute will be of type real
1593!-
1594!- INPUT
1595!-
1596!- fid      : File ID returned by flinopen
1597!- varname  : Name of the variable to be read from the file
1598!- iim      : | These three variables give the size of the variables
1599!- jjm      : | to be read. It will be verified that the variables
1600!- llm      : | fits in there.
1601!- ttm      : |
1602!- itau_dep : Time step at which we will start to read
1603!- itau_fin : Time step until which we are going to read
1604!-           For the moment this is done on indeces but it should be
1605!-           in the physical space
1606!-           If there is no time-axis in the file then use a
1607!-           itau_fin < itau_dep, this will tell flinget not to
1608!-           expect a time-axis in the file.
1609!-
1610!- OUTPUT
1611!-
1612!- var      : scalar that will contain the data
1613!---------------------------------------------------------------------
1614  IMPLICIT NONE
1615!-
1616! ARGUMENTS
1617!-
1618  INTEGER :: fid_in
1619  CHARACTER(LEN=*) :: varname
1620  INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin
1621  REAL :: var
1622!-
1623! LOCAL
1624!-
1625  INTEGER :: iret, fid
1626!-
1627  LOGICAL :: check = .FALSE.
1628!---------------------------------------------------------------------
1629  IF (check) THEN
1630    WRITE (*,*) 'flinget_scal in file with id ',fid_in
1631  ENDIF
1632!-
1633  fid = ncids(fid_in)
1634!-
1635! 1.0 Reading a global attribute
1636!-
1637  iret = NF90_GET_ATT (fid, NF90_GLOBAL, varname, var)
1638!---------------------------
1639END  SUBROUTINE flinget_scal
1640!-
1641!===
1642!-
1643SUBROUTINE flinfindcood (fid_in, axtype, vid, ndim)
1644!---------------------------------------------------------------------
1645!- This subroutine explores the file in order to find
1646!- the coordinate according to a number of rules
1647!---------------------------------------------------------------------
1648  IMPLICIT NONE
1649!-
1650! ARGUMENTS
1651!-
1652  INTEGER :: fid_in, vid, ndim
1653  CHARACTER(LEN=3) :: axtype
1654!-
1655! LOCAL
1656!-
1657  INTEGER :: iv, iret, dimnb
1658  CHARACTER(LEN=40) :: dimname, dimuni1, dimuni2, dimuni3
1659  CHARACTER(LEN=30) :: str1
1660  LOGICAL :: found_rule = .FALSE.
1661!---------------------------------------------------------------------
1662  vid = -1
1663!-
1664! Make sure all strings are invalid
1665!-
1666  dimname = '?-?'
1667  dimuni1 = '?-?'
1668  dimuni2 = '?-?'
1669  dimuni3 = '?-?'
1670!-
1671! First rule : we look for the correct units
1672! lon : east
1673! lat : north
1674! We make an exact check as it would be too easy to mistake
1675! some units by just comparing the substrings.
1676!-
1677  SELECTCASE(axtype)
1678  CASE ('lon')
1679    dimuni1 = 'degree_e'
1680    dimuni2 = 'degrees_e'
1681    found_rule = .TRUE.
1682  CASE('lat')
1683    dimuni1 = 'degree_n'
1684    dimuni2 = 'degrees_n'
1685    found_rule = .TRUE.
1686  CASE('lev')
1687    dimuni1 = 'm'
1688    dimuni2 = 'km'
1689    dimuni3 = 'hpa'
1690    found_rule = .TRUE.
1691  CASE DEFAULT
1692    found_rule = .FALSE.
1693  END SELECT
1694!-
1695  IF (found_rule) THEN
1696    iv = 0
1697    DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1698      iv = iv+1
1699      str1 = ''
1700      iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1)
1701      IF (iret == NF90_NOERR) THEN
1702        CALL strlowercase (str1)
1703        IF (    (INDEX(str1, TRIM(dimuni1)) == 1) &
1704            .OR.(INDEX(str1, TRIM(dimuni2)) == 1) &
1705            .OR.(INDEX(str1, TRIM(dimuni3)) == 1) ) THEN
1706          vid = iv
1707          iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim)
1708        ENDIF
1709      ENDIF
1710    ENDDO
1711  ENDIF
1712!-
1713! Second rule : we find specific names :
1714! lon : nav_lon
1715! lat : nav_lat
1716! Here we can check if we find the substring as the
1717! names are more specific.
1718!-
1719  SELECTCASE(axtype)
1720  CASE ('lon')
1721    dimname = 'nav_lon lon longitude'
1722    found_rule = .TRUE.
1723  CASE('lat')
1724    dimname = 'nav_lat lat latitude'
1725    found_rule = .TRUE.
1726  CASE('lev')
1727    dimname = 'plev level depth deptht'
1728    found_rule = .TRUE.
1729  CASE DEFAULT
1730    found_rule = .FALSE.
1731  END SELECT
1732!-
1733  IF (found_rule) THEN
1734    iv = 0
1735    DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1736      iv = iv+1
1737      str1=''
1738      iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
1739                                    name=str1, ndims=ndim)
1740      IF (INDEX(dimname,TRIM(str1)) >= 1) THEN
1741        vid = iv
1742      ENDIF
1743    ENDDO
1744  ENDIF
1745!-
1746! Third rule : we find a variable with the same name as the dimension
1747! lon = 1
1748! lat = 2
1749! lev = 3
1750!-
1751  IF (vid < 0) THEN
1752    SELECTCASE(axtype)
1753    CASE ('lon')
1754      dimnb = 1
1755      found_rule = .TRUE.
1756    CASE('lat')
1757      dimnb = 2
1758      found_rule = .TRUE.
1759    CASE('lev')
1760      dimnb = 3
1761      found_rule = .TRUE.
1762    CASE DEFAULT
1763      found_rule = .FALSE.
1764    END SELECT
1765!---
1766    IF (found_rule) THEN
1767      iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname)
1768      iv = 0
1769      DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )
1770        iv = iv+1
1771        str1=''
1772        iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
1773                                      name=str1, ndims=ndim)
1774        IF (INDEX(dimname,TRIM(str1)) == 1) THEN
1775          vid = iv
1776        ENDIF
1777      ENDDO
1778    ENDIF
1779  ENDIF
1780!-
1781! Stop the program if no coordinate was found
1782!-
1783  IF (vid < 0) THEN
1784    CALL histerr (3,'flinfindcood', &
1785           'No coordinate axis was found in the file', &
1786           'The data in this file can not be used', axtype)
1787  ENDIF
1788!--------------------------
1789END SUBROUTINE flinfindcood
1790!-
1791!===
1792!-
1793SUBROUTINE flinclo (fid_in)
1794!---------------------------------------------------------------------
1795  IMPLICIT NONE
1796!-
1797  INTEGER :: fid_in
1798!-
1799  INTEGER :: iret
1800!---------------------------------------------------------------------
1801  iret = NF90_CLOSE (ncids(fid_in))
1802  ncfileopen(fid_in) = .FALSE.
1803!---------------------
1804END SUBROUTINE flinclo
1805!-
1806!===
1807!-
1808SUBROUTINE flinquery_var(fid_in, varname, exists)
1809!---------------------------------------------------------------------
1810!- Queries the existance of a variable in the file.
1811!---------------------------------------------------------------------
1812  IMPLICIT NONE
1813!-
1814  INTEGER :: fid_in
1815  CHARACTER(LEN=*) varname
1816  LOGICAL :: exists
1817!-
1818  INTEGER :: iret, fid, vid
1819!---------------------------------------------------------------------
1820  fid = ncids(fid_in)
1821  vid = -1
1822  iret = NF90_INQ_VARID (fid, varname, vid)
1823!-
1824  exists = ( (vid >= 0).AND.(iret == NF90_NOERR) )
1825!---------------------------
1826END SUBROUTINE flinquery_var
1827!-
1828!===
1829!-
1830SUBROUTINE flininspect (fid_in)
1831!---------------------------------------------------------------------
1832  IMPLICIT NONE
1833!-
1834! fid : File id to inspect
1835!-
1836  INTEGER :: fid_in
1837!-
1838!- LOCAL
1839!-
1840  INTEGER :: iim, jjm, llm, ttm, fid_out
1841  INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim
1842  INTEGER :: iv, in, lll
1843  INTEGER :: xid, yid, zid, tid
1844  INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: idimid
1845  CHARACTER(LEN=80) :: name
1846  CHARACTER(LEN=30) :: axname
1847!---------------------------------------------------------------------
1848  fid = ncids(fid_in)
1849!-
1850  iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, &
1851                       nAttributes=nb_atts, unlimitedDimId=id_unlim)
1852!-
1853  WRITE (*,*) 'IOIPSL ID                   : ',fid_in
1854  WRITE (*,*) 'NetCDF ID                   : ',fid
1855  WRITE (*,*) 'Number of dimensions        : ',ndims
1856  WRITE (*,*) 'Number of variables         : ',nvars
1857  WRITE (*,*) 'Number of global attributes : ',nb_atts
1858  WRITE (*,*) 'ID unlimited                : ',id_unlim
1859!-
1860  xid = -1; iim = 0;
1861  yid = -1; jjm = 0;
1862  zid = -1; llm = 0;
1863  tid = -1; ttm = 0;
1864!-
1865  DO iv=1,ndims
1866!---
1867    iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll)
1868    CALL strlowercase (axname)
1869    axname = ADJUSTL(axname)
1870!---
1871    WRITE (*,*) 'Dimension number : ',iv
1872    WRITE (*,*) 'Dimension name   : ',TRIM(axname)
1873!---
1874    IF      (    (INDEX(axname,'x') == 1) &
1875             .OR.(INDEX(axname,'lon') == 1)) THEN
1876      xid = iv; iim = lll;
1877      WRITE (*,*) 'Dimension X size   : ',iim
1878    ELSE IF (    (INDEX(axname,'y') == 1) &
1879             .OR.(INDEX(axname,'lat') == 1)) THEN
1880      yid = iv; jjm = lll;
1881      WRITE (*,*) 'Dimension Y size   : ',jjm
1882    ELSE IF (    (INDEX(axname,'lev') == 1) &
1883             .OR.(INDEX(axname,'plev') == 1) &
1884             .OR.(INDEX(axname,'z') == 1) &
1885             .OR.(INDEX(axname,'depth') == 1)) THEN
1886      zid = iv; llm = lll;
1887      WRITE (*,*) 'Dimension Z size   : ',llm
1888    ELSE IF (    (INDEX(axname,'tstep') == 1) &
1889             .OR.(INDEX(axname,'time_counter') == 1)) THEN
1890!---- For the time we certainly need to allow for other names
1891      tid = iv; ttm = lll;
1892    ELSE IF (ndims == 1) THEN
1893!---- Nothing was found and ndims=1 then we have a vector of data
1894      xid = 1; iim = lll;
1895    ENDIF
1896!---
1897  ENDDO
1898!-
1899! Keep all this information
1900!-
1901  nbfiles = nbfiles+1
1902!-
1903  IF (nbfiles > nbfile_max) THEN
1904    CALL histerr(3,'flininspect', &
1905      'Too many files. Please increase nbfil_max', &
1906      'in program flincom.F90.',' ')
1907  ENDIF
1908!-
1909  ncids(nbfiles) = fid
1910  ncnbd(nbfiles) = ndims
1911!-
1912  ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /)
1913!-
1914  ncfunli(nbfiles) = id_unlim
1915  ncnba(nbfiles)   = nb_atts
1916  ncnbva(nbfiles)  = nvars
1917  ncfileopen(nbfiles) = .TRUE.
1918!-
1919  fid_out = nbfiles
1920!-
1921  DO in=1,nvars
1922    iret = NF90_INQUIRE_VARIABLE (fid, in, &
1923             name=name, ndims=ndims, dimids=idimid, nAtts=nb_atts)
1924    WRITE (*,*) 'Variable number  ------------ > ', in
1925    WRITE (*,*) 'Variable name        : ', TRIM(name)
1926    WRITE (*,*) 'Number of dimensions : ', ndims
1927    WRITE (*,*) 'Dimensions ID''s     : ', idimid(1:ndims)
1928    WRITE (*,*) 'Number of attributes : ', nb_atts
1929  ENDDO
1930!-------------------------
1931END SUBROUTINE flininspect
1932!-
1933!===
1934!-
1935END MODULE flincom
Note: See TracBrowser for help on using the repository browser.