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

Last change on this file since 5793 was 4863, checked in by jgipsl, 5 years ago

Following changes have been done by A.Jornet/LSCE. No change is results and no change in usage have been seen. Some more error checking might stop the model for example if dimensions are not correct in call to histcom module.

Restcom:

  • Define a new var size length (20 to 100 )→ pbs found without no errors
  • Raise an error when var name is too long
  • Deallocate any buffer at the end of all restput/restcget calls → buffers only increase size. After loading/saving nothing is done with this memory

Histcom:

  • Raise an error if given history declared variables do not match with given dimensions from histwrite

getincom and stringop:

  • Enable any length character for the run.def → useful for long filepaths

flincom

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