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

Last change on this file since 1378 was 1378, checked in by mmaipsl, 13 years ago

Enhancement : use ipslout number from errioipsl to redirect all prints of IOIPSL
in the local process when use with parallelization.
This variable ipslout can be modified with ipslnlf function of errioipsl module.

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