New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
flincom.f90 in utils/tools/DOMAINcfg/src – NEMO

source: utils/tools/DOMAINcfg/src/flincom.f90 @ 14931

Last change on this file since 14931 was 14623, checked in by ldebreu, 3 years ago

AGFdomcfg: 1) Update DOMAINcfg to be compliant with the removal of halo cells 2) Update most of the LBC ... subroutines to a recent NEMO 4 version #2638

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