1 | MODULE 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 | !- |
---|
120 | CONTAINS |
---|
121 | !- |
---|
122 | !=== |
---|
123 | !- |
---|
124 | SUBROUTINE 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 | !--------------------- |
---|
295 | END SUBROUTINE flincre |
---|
296 | !- |
---|
297 | !=== |
---|
298 | !- |
---|
299 | SUBROUTINE 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 | !----------------------------- |
---|
328 | END SUBROUTINE flinopen_zoom2d |
---|
329 | !- |
---|
330 | !=== |
---|
331 | !- |
---|
332 | SUBROUTINE 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 | !------------------------- |
---|
355 | END SUBROUTINE flinopen_nozoom |
---|
356 | !- |
---|
357 | !=== |
---|
358 | !- |
---|
359 | SUBROUTINE 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 | !--------------------------- |
---|
670 | END SUBROUTINE flinopen_work |
---|
671 | !- |
---|
672 | !=== |
---|
673 | !- |
---|
674 | SUBROUTINE 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=80) :: 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 | !---------------------- |
---|
774 | END SUBROUTINE flininfo |
---|
775 | !- |
---|
776 | !=== |
---|
777 | !- |
---|
778 | SUBROUTINE 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 | !------------------------- |
---|
803 | END SUBROUTINE flinput_r1d |
---|
804 | !- |
---|
805 | !=== |
---|
806 | !- |
---|
807 | SUBROUTINE 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 | !------------------------- |
---|
832 | END SUBROUTINE flinput_r2d |
---|
833 | !- |
---|
834 | !=== |
---|
835 | !- |
---|
836 | SUBROUTINE 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 | !------------------------- |
---|
861 | END SUBROUTINE flinput_r3d |
---|
862 | !- |
---|
863 | !=== |
---|
864 | !- |
---|
865 | SUBROUTINE 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 | !------------------------- |
---|
890 | END SUBROUTINE flinput_r4d |
---|
891 | !- |
---|
892 | !=== |
---|
893 | !- |
---|
894 | SUBROUTINE 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 | !-------------------------- |
---|
936 | END SUBROUTINE flinput_mat |
---|
937 | !- |
---|
938 | !=== |
---|
939 | !- |
---|
940 | SUBROUTINE 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 | !--------------------------- |
---|
961 | END SUBROUTINE flinput_scal |
---|
962 | !- |
---|
963 | !=== |
---|
964 | !- |
---|
965 | SUBROUTINE 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 | END SUBROUTINE flinget_r1d |
---|
1002 | !- |
---|
1003 | !=== |
---|
1004 | !- |
---|
1005 | SUBROUTINE flinget_r2d & |
---|
1006 | (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var) |
---|
1007 | !--------------------------------------------------------------------- |
---|
1008 | IMPLICIT NONE |
---|
1009 | !- |
---|
1010 | INTEGER :: fid_in |
---|
1011 | CHARACTER(LEN=*) :: varname |
---|
1012 | INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin |
---|
1013 | REAL :: var(:,:) |
---|
1014 | !- |
---|
1015 | INTEGER :: jl, jj, ji |
---|
1016 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1017 | LOGICAL :: l_dbg |
---|
1018 | !--------------------------------------------------------------------- |
---|
1019 | CALL ipsldbg (old_status=l_dbg) |
---|
1020 | |
---|
1021 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1022 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1023 | "flinget_r2d : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1024 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1025 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1026 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1027 | "flinget_r2d : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1028 | DEALLOCATE (buff_tmp) |
---|
1029 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1030 | ENDIF |
---|
1031 | !- |
---|
1032 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1033 | itau_dep,itau_fin,1,iim,1,jjm,buff_tmp) |
---|
1034 | !- |
---|
1035 | jl=0 |
---|
1036 | DO jj=1,SIZE(var,2) |
---|
1037 | DO ji=1,SIZE(var,1) |
---|
1038 | jl=jl+1 |
---|
1039 | var(ji,jj) = buff_tmp(jl) |
---|
1040 | ENDDO |
---|
1041 | ENDDO |
---|
1042 | !------------------------- |
---|
1043 | END SUBROUTINE flinget_r2d |
---|
1044 | !- |
---|
1045 | !=== |
---|
1046 | !- |
---|
1047 | SUBROUTINE flinget_r2d_zoom2d & |
---|
1048 | (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1049 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var) |
---|
1050 | !--------------------------------------------------------------------- |
---|
1051 | IMPLICIT NONE |
---|
1052 | !- |
---|
1053 | INTEGER :: fid_in |
---|
1054 | CHARACTER(LEN=*) :: varname |
---|
1055 | INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen |
---|
1056 | REAL :: var(:,:) |
---|
1057 | !- |
---|
1058 | INTEGER :: jl, jj, ji |
---|
1059 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1060 | LOGICAL :: l_dbg |
---|
1061 | !--------------------------------------------------------------------- |
---|
1062 | CALL ipsldbg (old_status=l_dbg) |
---|
1063 | |
---|
1064 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1065 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1066 | "flinget_r2d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1067 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1068 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1069 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1070 | "flinget_r2d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1071 | DEALLOCATE (buff_tmp) |
---|
1072 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1073 | ENDIF |
---|
1074 | !- |
---|
1075 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1076 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp) |
---|
1077 | !- |
---|
1078 | jl=0 |
---|
1079 | DO jj=1,SIZE(var,2) |
---|
1080 | DO ji=1,SIZE(var,1) |
---|
1081 | jl=jl+1 |
---|
1082 | var(ji,jj) = buff_tmp(jl) |
---|
1083 | ENDDO |
---|
1084 | ENDDO |
---|
1085 | !-------------------------------- |
---|
1086 | END SUBROUTINE flinget_r2d_zoom2d |
---|
1087 | !- |
---|
1088 | !=== |
---|
1089 | !- |
---|
1090 | SUBROUTINE flinget_r3d & |
---|
1091 | (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var) |
---|
1092 | !--------------------------------------------------------------------- |
---|
1093 | IMPLICIT NONE |
---|
1094 | !- |
---|
1095 | INTEGER :: fid_in |
---|
1096 | CHARACTER(LEN=*) :: varname |
---|
1097 | INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin |
---|
1098 | REAL :: var(:,:,:) |
---|
1099 | !- |
---|
1100 | INTEGER :: jl, jk, jj, ji |
---|
1101 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1102 | LOGICAL :: l_dbg |
---|
1103 | !--------------------------------------------------------------------- |
---|
1104 | CALL ipsldbg (old_status=l_dbg) |
---|
1105 | |
---|
1106 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1107 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1108 | "flinget_r3d : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1109 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1110 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1111 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1112 | "flinget_r3d : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1113 | DEALLOCATE (buff_tmp) |
---|
1114 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1115 | ENDIF |
---|
1116 | !- |
---|
1117 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1118 | itau_dep,itau_fin,1,iim,1,jjm,buff_tmp) |
---|
1119 | !- |
---|
1120 | jl=0 |
---|
1121 | DO jk=1,SIZE(var,3) |
---|
1122 | DO jj=1,SIZE(var,2) |
---|
1123 | DO ji=1,SIZE(var,1) |
---|
1124 | jl=jl+1 |
---|
1125 | var(ji,jj,jk) = buff_tmp(jl) |
---|
1126 | ENDDO |
---|
1127 | ENDDO |
---|
1128 | ENDDO |
---|
1129 | !------------------------- |
---|
1130 | END SUBROUTINE flinget_r3d |
---|
1131 | !- |
---|
1132 | !=== |
---|
1133 | !- |
---|
1134 | SUBROUTINE flinget_r3d_zoom2d & |
---|
1135 | (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1136 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var) |
---|
1137 | !--------------------------------------------------------------------- |
---|
1138 | IMPLICIT NONE |
---|
1139 | !- |
---|
1140 | INTEGER :: fid_in |
---|
1141 | CHARACTER(LEN=*) :: varname |
---|
1142 | INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen |
---|
1143 | REAL :: var(:,:,:) |
---|
1144 | !- |
---|
1145 | INTEGER :: jl, jk, jj, ji |
---|
1146 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1147 | LOGICAL :: l_dbg |
---|
1148 | !--------------------------------------------------------------------- |
---|
1149 | CALL ipsldbg (old_status=l_dbg) |
---|
1150 | |
---|
1151 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1152 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1153 | "flinget_r3d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1154 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1155 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1156 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1157 | "flinget_r3d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1158 | DEALLOCATE (buff_tmp) |
---|
1159 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1160 | ENDIF |
---|
1161 | !- |
---|
1162 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1163 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp) |
---|
1164 | !- |
---|
1165 | jl=0 |
---|
1166 | DO jk=1,SIZE(var,3) |
---|
1167 | DO jj=1,SIZE(var,2) |
---|
1168 | DO ji=1,SIZE(var,1) |
---|
1169 | jl=jl+1 |
---|
1170 | var(ji,jj,jk) = buff_tmp(jl) |
---|
1171 | ENDDO |
---|
1172 | ENDDO |
---|
1173 | ENDDO |
---|
1174 | !-------------------------------- |
---|
1175 | END SUBROUTINE flinget_r3d_zoom2d |
---|
1176 | !- |
---|
1177 | !=== |
---|
1178 | !- |
---|
1179 | SUBROUTINE flinget_r4d & |
---|
1180 | (fid_in,varname,iim,jjm,llm,ttm,itau_dep,itau_fin,var) |
---|
1181 | !--------------------------------------------------------------------- |
---|
1182 | IMPLICIT NONE |
---|
1183 | !- |
---|
1184 | INTEGER :: fid_in |
---|
1185 | CHARACTER(LEN=*) :: varname |
---|
1186 | INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin |
---|
1187 | REAL :: var(:,:,:,:) |
---|
1188 | !- |
---|
1189 | INTEGER :: jl, jk, jj, ji, jm |
---|
1190 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1191 | LOGICAL :: l_dbg |
---|
1192 | !--------------------------------------------------------------------- |
---|
1193 | CALL ipsldbg (old_status=l_dbg) |
---|
1194 | |
---|
1195 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1196 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1197 | "flinget_r4d : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1198 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1199 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1200 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1201 | "flinget_r4d : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1202 | DEALLOCATE (buff_tmp) |
---|
1203 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1204 | ENDIF |
---|
1205 | !- |
---|
1206 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1207 | itau_dep,itau_fin,1,iim,1,jjm,buff_tmp) |
---|
1208 | !- |
---|
1209 | jl=0 |
---|
1210 | DO jm=1,SIZE(var,4) |
---|
1211 | DO jk=1,SIZE(var,3) |
---|
1212 | DO jj=1,SIZE(var,2) |
---|
1213 | DO ji=1,SIZE(var,1) |
---|
1214 | jl=jl+1 |
---|
1215 | var(ji,jj,jk,jm) = buff_tmp(jl) |
---|
1216 | ENDDO |
---|
1217 | ENDDO |
---|
1218 | ENDDO |
---|
1219 | ENDDO |
---|
1220 | !------------------------- |
---|
1221 | END SUBROUTINE flinget_r4d |
---|
1222 | !- |
---|
1223 | !=== |
---|
1224 | !- |
---|
1225 | SUBROUTINE flinget_r4d_zoom2d & |
---|
1226 | (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1227 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,var) |
---|
1228 | !--------------------------------------------------------------------- |
---|
1229 | IMPLICIT NONE |
---|
1230 | !- |
---|
1231 | INTEGER :: fid_in |
---|
1232 | CHARACTER(LEN=*) :: varname |
---|
1233 | INTEGER :: iim,jjm,llm,ttm,itau_dep,itau_fin,iideb,jjdeb,iilen,jjlen |
---|
1234 | REAL :: var(:,:,:,:) |
---|
1235 | !- |
---|
1236 | INTEGER :: jl, jk, jj, ji, jm |
---|
1237 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: buff_tmp |
---|
1238 | LOGICAL :: l_dbg |
---|
1239 | !--------------------------------------------------------------------- |
---|
1240 | CALL ipsldbg (old_status=l_dbg) |
---|
1241 | |
---|
1242 | IF (.NOT.ALLOCATED(buff_tmp)) THEN |
---|
1243 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1244 | "flinget_r4d_zoom : allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1245 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1246 | ELSE IF (SIZE(var) > SIZE(buff_tmp)) THEN |
---|
1247 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1248 | "flinget_r4d_zoom : re-allocate buff_tmp for buff_sz = ",SIZE(var) |
---|
1249 | DEALLOCATE (buff_tmp) |
---|
1250 | ALLOCATE (buff_tmp(SIZE(var))) |
---|
1251 | ENDIF |
---|
1252 | !- |
---|
1253 | CALL flinget_mat (fid_in,varname,iim,jjm,llm,ttm, & |
---|
1254 | itau_dep,itau_fin,iideb,iilen,jjdeb,jjlen,buff_tmp) |
---|
1255 | !- |
---|
1256 | jl=0 |
---|
1257 | DO jm=1,SIZE(var,4) |
---|
1258 | DO jk=1,SIZE(var,3) |
---|
1259 | DO jj=1,SIZE(var,2) |
---|
1260 | DO ji=1,SIZE(var,1) |
---|
1261 | jl=jl+1 |
---|
1262 | var(ji,jj,jk,jm) = buff_tmp(jl) |
---|
1263 | ENDDO |
---|
1264 | ENDDO |
---|
1265 | ENDDO |
---|
1266 | ENDDO |
---|
1267 | !-------------------------------- |
---|
1268 | END SUBROUTINE flinget_r4d_zoom2d |
---|
1269 | !- |
---|
1270 | !=== |
---|
1271 | !- |
---|
1272 | SUBROUTINE flinget_mat & |
---|
1273 | (fid_in, varname, iim, jjm, llm, ttm, itau_dep, & |
---|
1274 | itau_fin, iideb, iilen, jjdeb, jjlen, var) |
---|
1275 | !--------------------------------------------------------------------- |
---|
1276 | !- This subroutine will read the variable named varname from |
---|
1277 | !- the file previously opened by flinopen and identified by fid |
---|
1278 | !- |
---|
1279 | !- It is checked that the dimensions of the variable to be read |
---|
1280 | !- correspond to what the user requested when he specified |
---|
1281 | !- iim, jjm and llm. The only exception which is allowed is |
---|
1282 | !- for compressed data where the horizontal grid is not expected |
---|
1283 | !- to be iim x jjm. |
---|
1284 | !- |
---|
1285 | !- If variable is of size zero a global attribute is read. |
---|
1286 | !- This global attribute will be of type real |
---|
1287 | !- |
---|
1288 | !- INPUT |
---|
1289 | !- |
---|
1290 | !- fid : File ID returned by flinopen |
---|
1291 | !- varname : Name of the variable to be read from the file |
---|
1292 | !- iim : | These three variables give the size of the variables |
---|
1293 | !- jjm : | to be read. It will be verified that the variables |
---|
1294 | !- llm : | fits in there. |
---|
1295 | !- ttm : | |
---|
1296 | !- itau_dep : Time step at which we will start to read |
---|
1297 | !- itau_fin : Time step until which we are going to read |
---|
1298 | !- For the moment this is done on indexes |
---|
1299 | !- but it should be in the physical space. |
---|
1300 | !- If there is no time-axis in the file then use a |
---|
1301 | !- itau_fin < itau_dep, this will tell flinget not to |
---|
1302 | !- expect a time-axis in the file. |
---|
1303 | !- iideb : index i for zoom |
---|
1304 | !- iilen : length of zoom |
---|
1305 | !- jjdeb : index j for zoom |
---|
1306 | !- jjlen : length of zoom |
---|
1307 | !- |
---|
1308 | !- OUTPUT |
---|
1309 | !- |
---|
1310 | !- var : array that will contain the data |
---|
1311 | !--------------------------------------------------------------------- |
---|
1312 | IMPLICIT NONE |
---|
1313 | !- |
---|
1314 | ! ARGUMENTS |
---|
1315 | !- |
---|
1316 | INTEGER, INTENT(IN) :: fid_in |
---|
1317 | CHARACTER(LEN=*), INTENT(IN) :: varname |
---|
1318 | INTEGER, INTENT(IN) :: iim, jjm, llm, ttm |
---|
1319 | INTEGER, INTENT(IN) :: itau_dep, itau_fin, iideb, iilen, jjdeb, jjlen |
---|
1320 | REAL, INTENT(OUT) :: var(:) |
---|
1321 | !- |
---|
1322 | ! LOCAL |
---|
1323 | !- |
---|
1324 | INTEGER :: iret, fid |
---|
1325 | INTEGER :: vid, cvid, clen |
---|
1326 | CHARACTER(LEN=70) :: str1 |
---|
1327 | CHARACTER(LEN=250) :: att_n, tmp_n |
---|
1328 | CHARACTER(LEN=5) :: axs_l |
---|
1329 | INTEGER :: tmp_i |
---|
1330 | REAL,SAVE :: mis_v=0. |
---|
1331 | REAL :: tmp_r |
---|
1332 | INTEGER :: ndims, x_typ, nb_atts |
---|
1333 | INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: dimids |
---|
1334 | INTEGER :: i, nvars, i2d, cnd |
---|
1335 | REAL,DIMENSION(:),ALLOCATABLE,SAVE :: var_tmp |
---|
1336 | INTEGER :: itau_len |
---|
1337 | LOGICAL :: uncompress = .FALSE. |
---|
1338 | INTEGER :: il, ip, i2p, it |
---|
1339 | !- |
---|
1340 | LOGICAL :: l_dbg |
---|
1341 | !--------------------------------------------------------------------- |
---|
1342 | CALL ipsldbg (old_status=l_dbg) |
---|
1343 | !- |
---|
1344 | fid = ncids(fid_in) |
---|
1345 | !- |
---|
1346 | IF (l_dbg) THEN |
---|
1347 | WRITE(ipslout,*) & |
---|
1348 | 'flinget_mat : fid_in, fid, varname :', fid_in, fid, TRIM(varname) |
---|
1349 | WRITE(ipslout,*) & |
---|
1350 | 'flinget_mat : iim, jjm, llm, ttm, itau_dep, itau_fin :', & |
---|
1351 | iim, jjm, llm, ttm, itau_dep, itau_fin |
---|
1352 | WRITE(ipslout,*) & |
---|
1353 | 'flinget_mat : iideb, iilen, jjdeb, jjlen :', & |
---|
1354 | iideb, iilen, jjdeb, jjlen |
---|
1355 | ENDIF |
---|
1356 | !- |
---|
1357 | uncompress = .FALSE. |
---|
1358 | !- |
---|
1359 | ! 1.0 We get first all the details on this variable from the file |
---|
1360 | !- |
---|
1361 | nvars = ncnbva(fid_in) |
---|
1362 | !- |
---|
1363 | vid = -1 |
---|
1364 | iret = NF90_INQ_VARID (fid, varname, vid) |
---|
1365 | !- |
---|
1366 | IF (vid < 0 .OR. iret /= NF90_NOERR) THEN |
---|
1367 | CALL histerr (3,'flinget', & |
---|
1368 | 'Variable '//TRIM(varname)//' not found in file',' ',' ') |
---|
1369 | ENDIF |
---|
1370 | !- |
---|
1371 | iret = NF90_INQUIRE_VARIABLE (fid, vid, & |
---|
1372 | ndims=ndims, dimids=dimids, nAtts=nb_atts) |
---|
1373 | IF (l_dbg) THEN |
---|
1374 | WRITE(ipslout,*) & |
---|
1375 | 'flinget_mat : fid, vid :', fid, vid |
---|
1376 | WRITE(ipslout,*) & |
---|
1377 | 'flinget_mat : ndims, dimids(1:ndims), nb_atts :', & |
---|
1378 | ndims, dimids(1:ndims), nb_atts |
---|
1379 | ENDIF |
---|
1380 | !- |
---|
1381 | w_dim(:) = 0 |
---|
1382 | DO i=1,ndims |
---|
1383 | iret = NF90_INQUIRE_DIMENSION (fid, dimids(i), len=w_dim(i)) |
---|
1384 | ENDDO |
---|
1385 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1386 | 'flinget_mat : w_dim :', w_dim(1:ndims) |
---|
1387 | !- |
---|
1388 | mis_v = 0.0; axs_l = ' '; |
---|
1389 | !- |
---|
1390 | IF (nb_atts > 0) THEN |
---|
1391 | IF (l_dbg) THEN |
---|
1392 | WRITE(ipslout,*) 'flinget_mat : attributes for variable :' |
---|
1393 | ENDIF |
---|
1394 | ENDIF |
---|
1395 | DO i=1,nb_atts |
---|
1396 | iret = NF90_INQ_ATTNAME (fid, vid, i, att_n) |
---|
1397 | iret = NF90_INQUIRE_ATTRIBUTE (fid, vid, att_n, xtype=x_typ) |
---|
1398 | CALL strlowercase (att_n) |
---|
1399 | IF ( (x_typ == NF90_INT).OR.(x_typ == NF90_SHORT) & |
---|
1400 | .OR.(x_typ == NF90_BYTE) ) THEN |
---|
1401 | iret = NF90_GET_ATT (fid, vid, att_n, tmp_i) |
---|
1402 | IF (l_dbg) THEN |
---|
1403 | WRITE(ipslout,*) ' ',TRIM(att_n),' : ',tmp_i |
---|
1404 | ENDIF |
---|
1405 | ELSE IF ( (x_typ == NF90_FLOAT).OR.(x_typ == NF90_DOUBLE) ) THEN |
---|
1406 | iret = NF90_GET_ATT (fid, vid, att_n, tmp_r) |
---|
1407 | IF (l_dbg) THEN |
---|
1408 | WRITE(ipslout,*) ' ',TRIM(att_n),' : ',tmp_r |
---|
1409 | ENDIF |
---|
1410 | IF (index(att_n,'missing_value') > 0) THEN |
---|
1411 | mis_v = tmp_r |
---|
1412 | ENDIF |
---|
1413 | ELSE |
---|
1414 | tmp_n = '' |
---|
1415 | iret = NF90_GET_ATT (fid, vid, att_n, tmp_n) |
---|
1416 | IF (l_dbg) THEN |
---|
1417 | WRITE(ipslout,*) ' ',TRIM(att_n),' : ',TRIM(tmp_n) |
---|
1418 | ENDIF |
---|
1419 | IF (index(att_n,'axis') > 0) THEN |
---|
1420 | axs_l = tmp_n |
---|
1421 | ENDIF |
---|
1422 | ENDIF |
---|
1423 | ENDDO |
---|
1424 | !? |
---|
1425 | !!!!!!!!!! We will need a verification on the type of the variable |
---|
1426 | !? |
---|
1427 | !- |
---|
1428 | ! 2.0 The dimensions are analysed to determine what is to be read |
---|
1429 | !- |
---|
1430 | ! 2.1 the longitudes |
---|
1431 | !- |
---|
1432 | IF ( w_dim(1) /= iim .OR. w_dim(2) /= jjm) THEN |
---|
1433 | !--- |
---|
1434 | !-- There is a possibility that we have to deal with a compressed axis ! |
---|
1435 | !--- |
---|
1436 | iret = NF90_INQUIRE_DIMENSION (fid, dimids(1), & |
---|
1437 | name=tmp_n, len=clen) |
---|
1438 | iret = NF90_INQ_VARID (fid, tmp_n, cvid) |
---|
1439 | !--- |
---|
1440 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1441 | 'Dimname, iret , NF90_NOERR : ',TRIM(tmp_n),iret,NF90_NOERR |
---|
1442 | !--- |
---|
1443 | !-- If we have an axis which has the same name |
---|
1444 | !-- as the dimension we can see if it is compressed |
---|
1445 | !--- |
---|
1446 | !-- TODO TODO for zoom2d |
---|
1447 | !--- |
---|
1448 | IF (iret == NF90_NOERR) THEN |
---|
1449 | iret = NF90_GET_ATT (fid, cvid, 'compress', str1) |
---|
1450 | !----- |
---|
1451 | IF (iret == NF90_NOERR) THEN |
---|
1452 | iret = NF90_INQUIRE_VARIABLE (fid,cvid,xtype=x_typ,ndims=cnd) |
---|
1453 | !------- |
---|
1454 | IF ( cnd /= 1 .AND. x_typ /= NF90_INT) THEN |
---|
1455 | CALL histerr (3,'flinget', & |
---|
1456 | 'Variable '//TRIM(tmp_n)//' can not be a compressed axis', & |
---|
1457 | 'Either it has too many dimensions'// & |
---|
1458 | ' or it is not of type integer', ' ') |
---|
1459 | ELSE |
---|
1460 | !--------- |
---|
1461 | !-------- Let us see if we already have that index table |
---|
1462 | !--------- |
---|
1463 | IF ( (cind_len /= clen).OR.(cind_vid /= cvid) & |
---|
1464 | .OR.(cind_fid /= fid) ) THEN |
---|
1465 | IF (ALLOCATED(cindex)) DEALLOCATE(cindex) |
---|
1466 | ALLOCATE(cindex(clen)) |
---|
1467 | cind_len = clen |
---|
1468 | cind_vid = cvid |
---|
1469 | cind_fid = fid |
---|
1470 | iret = NF90_GET_VAR (fid, cvid, cindex) |
---|
1471 | ENDIF |
---|
1472 | !--------- |
---|
1473 | !-------- In any case we need to set the slab of data to be read |
---|
1474 | !--------- |
---|
1475 | uncompress = .TRUE. |
---|
1476 | w_sta(1) = 1 |
---|
1477 | w_len(1) = clen |
---|
1478 | i2d = 1 |
---|
1479 | ENDIF |
---|
1480 | ELSE |
---|
1481 | str1 = 'The horizontal dimensions of '//varname |
---|
1482 | CALL histerr (3,'flinget',str1, & |
---|
1483 | 'is not compressed and does not'// & |
---|
1484 | ' correspond to the requested size',' ') |
---|
1485 | ENDIF |
---|
1486 | ELSE |
---|
1487 | IF (w_dim(1) /= iim) THEN |
---|
1488 | str1 = 'The longitude dimension of '//varname |
---|
1489 | CALL histerr (3,'flinget',str1, & |
---|
1490 | 'in the file is not equal to the dimension', & |
---|
1491 | 'that should be read') |
---|
1492 | ENDIF |
---|
1493 | IF (w_dim(2) /= jjm) THEN |
---|
1494 | str1 = 'The latitude dimension of '//varname |
---|
1495 | CALL histerr (3,'flinget',str1, & |
---|
1496 | 'in the file is not equal to the dimension', & |
---|
1497 | 'that should be read') |
---|
1498 | ENDIF |
---|
1499 | ENDIF |
---|
1500 | ELSE |
---|
1501 | w_sta(1:2) = (/ iideb, jjdeb /) |
---|
1502 | w_len(1:2) = (/ iilen, jjlen /) |
---|
1503 | i2d = 2 |
---|
1504 | ENDIF |
---|
1505 | !- |
---|
1506 | ! 2.3 Now the difficult part, the 3rd dimension which can be |
---|
1507 | ! time or levels. |
---|
1508 | !- |
---|
1509 | ! Priority is given to the time axis if only three axes are present. |
---|
1510 | !- |
---|
1511 | IF (ndims > i2d) THEN |
---|
1512 | !--- |
---|
1513 | !-- 2.3.1 We have a vertical axis |
---|
1514 | !--- |
---|
1515 | IF (llm == 1 .AND. ndims == i2d+2 .OR. llm == w_dim(i2d+1)) THEN |
---|
1516 | !----- |
---|
1517 | IF (w_dim(i2d+1) /= llm) THEN |
---|
1518 | CALL histerr (3,'flinget', & |
---|
1519 | 'The vertical dimension of '//varname, & |
---|
1520 | 'in the file is not equal to the dimension', & |
---|
1521 | 'that should be read') |
---|
1522 | ELSE |
---|
1523 | w_sta(i2d+1) = 1 |
---|
1524 | IF (llm > 0) THEN |
---|
1525 | w_len(i2d+1) = llm |
---|
1526 | ELSE |
---|
1527 | w_len(i2d+1) = w_sta(i2d+1) |
---|
1528 | ENDIF |
---|
1529 | ENDIF |
---|
1530 | !----- |
---|
1531 | IF ((itau_fin-itau_dep) >= 0) THEN |
---|
1532 | IF (ndims /= i2d+2) THEN |
---|
1533 | CALL histerr (3,'flinget', & |
---|
1534 | 'You attempt to read a time slab', & |
---|
1535 | 'but there is no time axis on this variable', varname) |
---|
1536 | ELSE IF ((itau_fin - itau_dep) <= w_dim(i2d+2)) THEN |
---|
1537 | w_sta(i2d+2) = itau_dep |
---|
1538 | w_len(i2d+2) = itau_fin-itau_dep+1 |
---|
1539 | ELSE |
---|
1540 | CALL histerr (3,'flinget', & |
---|
1541 | 'The time step you try to read is not', & |
---|
1542 | 'in the file (1)', varname) |
---|
1543 | ENDIF |
---|
1544 | ELSE IF (ndims == i2d+2 .AND. w_dim(i2d+2) > 1) THEN |
---|
1545 | CALL histerr (3,'flinget', & |
---|
1546 | 'There is a time axis in the file but no', & |
---|
1547 | 'time step give in the call', varname) |
---|
1548 | ELSE |
---|
1549 | w_sta(i2d+2) = 1 |
---|
1550 | w_len(i2d+2) = 1 |
---|
1551 | ENDIF |
---|
1552 | ELSE |
---|
1553 | !----- |
---|
1554 | !---- 2.3.2 We do not have any vertical axis |
---|
1555 | !----- |
---|
1556 | IF (ndims == i2d+2) THEN |
---|
1557 | CALL histerr (3,'flinget', & |
---|
1558 | 'The file contains 4 dimensions', & |
---|
1559 | 'but only 3 are requestes for variable ', varname) |
---|
1560 | ENDIF |
---|
1561 | IF ((itau_fin-itau_dep) >= 0) THEN |
---|
1562 | IF (ndims == i2d+1) THEN |
---|
1563 | IF ((itau_fin-itau_dep) < w_dim(i2d+1) ) THEN |
---|
1564 | w_sta(i2d+1) = itau_dep |
---|
1565 | w_len(i2d+1) = itau_fin-itau_dep+1 |
---|
1566 | ELSE |
---|
1567 | CALL histerr (3,'flinget', & |
---|
1568 | 'The time step you try to read is not', & |
---|
1569 | 'in the file (2)', varname) |
---|
1570 | ENDIF |
---|
1571 | ELSE |
---|
1572 | CALL histerr (3,'flinget', & |
---|
1573 | 'From your input you sould have 3 dimensions', & |
---|
1574 | 'in the file but there are 4', varname) |
---|
1575 | ENDIF |
---|
1576 | ELSE |
---|
1577 | IF (ndims == i2d+1 .AND. w_dim(i2d+1) > 1) THEN |
---|
1578 | CALL histerr (3,'flinget', & |
---|
1579 | 'There is a time axis in the file but no', & |
---|
1580 | 'time step given in the call', varname) |
---|
1581 | ELSE |
---|
1582 | w_sta(i2d+1) = 1 |
---|
1583 | w_len(i2d+1) = 1 |
---|
1584 | ENDIF |
---|
1585 | ENDIF |
---|
1586 | ENDIF |
---|
1587 | ELSE |
---|
1588 | !--- |
---|
1589 | !-- 2.3.3 We do not have any vertical axis |
---|
1590 | !--- |
---|
1591 | w_sta(i2d+1:i2d+2) = (/ 0, 0 /) |
---|
1592 | w_len(i2d+1:i2d+2) = (/ 0, 0 /) |
---|
1593 | ENDIF |
---|
1594 | !- |
---|
1595 | ! 3.0 Reading the data |
---|
1596 | !- |
---|
1597 | IF (l_dbg) WRITE(ipslout,*) & |
---|
1598 | 'flinget_mat 3.0 : ', uncompress, w_sta, w_len |
---|
1599 | !--- |
---|
1600 | var(:) = mis_v |
---|
1601 | IF (uncompress) THEN |
---|
1602 | !--- |
---|
1603 | IF (ALLOCATED(var_tmp)) THEN |
---|
1604 | IF (SIZE(var_tmp) < PRODUCT(w_len(:),mask=(w_len>1))) THEN |
---|
1605 | DEALLOCATE(var_tmp) |
---|
1606 | ALLOCATE(var_tmp(PRODUCT(w_len(:),mask=(w_len>1)))) |
---|
1607 | ENDIF |
---|
1608 | ELSE |
---|
1609 | ALLOCATE(var_tmp(PRODUCT(w_len(:),mask=(w_len>1)))) |
---|
1610 | ENDIF |
---|
1611 | !--- |
---|
1612 | iret = NF90_GET_VAR (fid, vid, var_tmp, & |
---|
1613 | start=w_sta(:), count=w_len(:)) |
---|
1614 | !--- |
---|
1615 | itau_len=itau_fin-itau_dep+1 |
---|
1616 | IF (l_dbg) WRITE(ipslout,*) 'flinget_mat 3.0 : clen, itau_len ',clen,itau_len |
---|
1617 | var(:) = mis_v |
---|
1618 | IF (itau_len > 0) THEN |
---|
1619 | DO it=1,itau_len |
---|
1620 | DO il=1,clen |
---|
1621 | ip = il + (it-1)*clen |
---|
1622 | i2p = cindex(il)+(it-1)*iim*jjm |
---|
1623 | var(i2p) = var_tmp(ip) |
---|
1624 | ENDDO |
---|
1625 | ENDDO |
---|
1626 | ELSE |
---|
1627 | var(cindex(:)) = var_tmp(:) |
---|
1628 | ENDIF |
---|
1629 | !--- |
---|
1630 | ELSE |
---|
1631 | iret = NF90_GET_VAR (fid, vid, var, & |
---|
1632 | start=w_sta(:), count=w_len(:)) |
---|
1633 | ENDIF |
---|
1634 | !- |
---|
1635 | IF (l_dbg) WRITE(ipslout,*) 'flinget_mat 3.1 : ',NF90_STRERROR (iret) |
---|
1636 | !-------------------------- |
---|
1637 | END SUBROUTINE flinget_mat |
---|
1638 | !- |
---|
1639 | !=== |
---|
1640 | !- |
---|
1641 | SUBROUTINE flinget_scal & |
---|
1642 | (fid_in, varname, iim, jjm, llm, ttm, itau_dep, itau_fin, var) |
---|
1643 | !--------------------------------------------------------------------- |
---|
1644 | !- This subroutine will read the variable named varname from |
---|
1645 | !- the file previously opened by flinopen and identified by fid |
---|
1646 | !- |
---|
1647 | !- If variable is of size zero a global attribute is read. This |
---|
1648 | !- global attribute will be of type real |
---|
1649 | !- |
---|
1650 | !- INPUT |
---|
1651 | !- |
---|
1652 | !- fid : File ID returned by flinopen |
---|
1653 | !- varname : Name of the variable to be read from the file |
---|
1654 | !- iim : | These three variables give the size of the variables |
---|
1655 | !- jjm : | to be read. It will be verified that the variables |
---|
1656 | !- llm : | fits in there. |
---|
1657 | !- ttm : | |
---|
1658 | !- itau_dep : Time step at which we will start to read |
---|
1659 | !- itau_fin : Time step until which we are going to read |
---|
1660 | !- For the moment this is done on indeces but it should be |
---|
1661 | !- in the physical space |
---|
1662 | !- If there is no time-axis in the file then use a |
---|
1663 | !- itau_fin < itau_dep, this will tell flinget not to |
---|
1664 | !- expect a time-axis in the file. |
---|
1665 | !- |
---|
1666 | !- OUTPUT |
---|
1667 | !- |
---|
1668 | !- var : scalar that will contain the data |
---|
1669 | !--------------------------------------------------------------------- |
---|
1670 | IMPLICIT NONE |
---|
1671 | !- |
---|
1672 | ! ARGUMENTS |
---|
1673 | !- |
---|
1674 | INTEGER :: fid_in |
---|
1675 | CHARACTER(LEN=*) :: varname |
---|
1676 | INTEGER :: iim, jjm, llm, ttm, itau_dep, itau_fin |
---|
1677 | REAL :: var |
---|
1678 | !- |
---|
1679 | ! LOCAL |
---|
1680 | !- |
---|
1681 | INTEGER :: iret, fid, vid |
---|
1682 | INTEGER :: attlen, attnum |
---|
1683 | INTEGER :: ndims, nb_atts |
---|
1684 | INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: dimids |
---|
1685 | LOGICAL :: var_exists |
---|
1686 | !- |
---|
1687 | LOGICAL :: l_dbg |
---|
1688 | !--------------------------------------------------------------------- |
---|
1689 | CALL ipsldbg (old_status=l_dbg) |
---|
1690 | |
---|
1691 | IF (l_dbg) THEN |
---|
1692 | WRITE (ipslout,*) 'flinget_scal in file with id ',fid_in |
---|
1693 | ENDIF |
---|
1694 | !- |
---|
1695 | fid = ncids(fid_in) |
---|
1696 | iret = NF90_INQUIRE_ATTRIBUTE(fid, NF90_GLOBAL, varname, len=attlen, attnum=attnum) |
---|
1697 | !- |
---|
1698 | ! 1.0 Reading a global attribute |
---|
1699 | !- |
---|
1700 | IF ( iret == nf90_noerr ) THEN |
---|
1701 | ! |
---|
1702 | ! This seems to be a Global attribute |
---|
1703 | ! |
---|
1704 | iret = NF90_GET_ATT (fid, NF90_GLOBAL, varname, var) |
---|
1705 | ELSE |
---|
1706 | ! |
---|
1707 | ! If there was an error on the test for a global attribute it |
---|
1708 | ! is perhaps a scalar variable. |
---|
1709 | ! |
---|
1710 | vid = -1 |
---|
1711 | iret = NF90_INQ_VARID (fid, varname, vid) |
---|
1712 | ! |
---|
1713 | IF ( (vid >= 0).AND.(iret == NF90_NOERR) ) THEN |
---|
1714 | iret = NF90_INQUIRE_VARIABLE (fid, vid, & |
---|
1715 | ndims=ndims, dimids=dimids, nAtts=nb_atts) |
---|
1716 | IF ( (ndims == 0) .AND. (nb_atts >= 0) ) THEN |
---|
1717 | iret = NF90_GET_VAR(fid, vid, var) |
---|
1718 | ELSE |
---|
1719 | CALL histerr (3,'flinget_scal', & |
---|
1720 | 'The variable has coordinates and thus is probably not a scalar.', & |
---|
1721 | 'Check your netCDF file.', " ") |
---|
1722 | ENDIF |
---|
1723 | ENDIF |
---|
1724 | IF (l_dbg) THEN |
---|
1725 | WRITE(ipslout,*) "Reading a Scalar value for varibale ", varname," It has value ", var |
---|
1726 | ENDIF |
---|
1727 | ENDIF |
---|
1728 | !- |
---|
1729 | !--------------------------- |
---|
1730 | END SUBROUTINE flinget_scal |
---|
1731 | !- |
---|
1732 | !=== |
---|
1733 | !- |
---|
1734 | SUBROUTINE flinfindcood (fid_in, axtype, vid, ndim) |
---|
1735 | !--------------------------------------------------------------------- |
---|
1736 | !- This subroutine explores the file in order to find |
---|
1737 | !- the coordinate according to a number of rules |
---|
1738 | !--------------------------------------------------------------------- |
---|
1739 | IMPLICIT NONE |
---|
1740 | !- |
---|
1741 | ! ARGUMENTS |
---|
1742 | !- |
---|
1743 | INTEGER :: fid_in, vid, ndim |
---|
1744 | CHARACTER(LEN=3) :: axtype |
---|
1745 | !- |
---|
1746 | ! LOCAL |
---|
1747 | !- |
---|
1748 | INTEGER :: iv, iret, dimnb |
---|
1749 | CHARACTER(LEN=40) :: dimname, dimuni1, dimuni2, dimuni3 |
---|
1750 | CHARACTER(LEN=80) :: str1 |
---|
1751 | LOGICAL :: found_rule = .FALSE. |
---|
1752 | !--------------------------------------------------------------------- |
---|
1753 | vid = -1 |
---|
1754 | !- |
---|
1755 | ! Make sure all strings are invalid |
---|
1756 | !- |
---|
1757 | dimname = '?-?' |
---|
1758 | dimuni1 = '?-?' |
---|
1759 | dimuni2 = '?-?' |
---|
1760 | dimuni3 = '?-?' |
---|
1761 | !- |
---|
1762 | ! First rule : we look for the correct units |
---|
1763 | ! lon : east |
---|
1764 | ! lat : north |
---|
1765 | ! We make an exact check as it would be too easy to mistake |
---|
1766 | ! some units by just comparing the substrings. |
---|
1767 | !- |
---|
1768 | SELECTCASE(axtype) |
---|
1769 | CASE ('lon') |
---|
1770 | dimuni1 = 'degree_e' |
---|
1771 | dimuni2 = 'degrees_e' |
---|
1772 | found_rule = .TRUE. |
---|
1773 | CASE('lat') |
---|
1774 | dimuni1 = 'degree_n' |
---|
1775 | dimuni2 = 'degrees_n' |
---|
1776 | found_rule = .TRUE. |
---|
1777 | CASE('lev') |
---|
1778 | dimuni1 = 'm' |
---|
1779 | dimuni2 = 'km' |
---|
1780 | dimuni3 = 'hpa' |
---|
1781 | found_rule = .TRUE. |
---|
1782 | CASE DEFAULT |
---|
1783 | found_rule = .FALSE. |
---|
1784 | END SELECT |
---|
1785 | !- |
---|
1786 | IF (found_rule) THEN |
---|
1787 | iv = 0 |
---|
1788 | DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
---|
1789 | iv = iv+1 |
---|
1790 | str1 = '' |
---|
1791 | iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1) |
---|
1792 | IF (iret == NF90_NOERR) THEN |
---|
1793 | CALL strlowercase (str1) |
---|
1794 | IF ( (INDEX(str1, TRIM(dimuni1)) == 1) & |
---|
1795 | .OR.(INDEX(str1, TRIM(dimuni2)) == 1) & |
---|
1796 | .OR.(INDEX(str1, TRIM(dimuni3)) == 1) ) THEN |
---|
1797 | vid = iv |
---|
1798 | iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim) |
---|
1799 | ENDIF |
---|
1800 | ENDIF |
---|
1801 | ENDDO |
---|
1802 | ENDIF |
---|
1803 | !- |
---|
1804 | ! Second rule : we find specific names : |
---|
1805 | ! lon : nav_lon |
---|
1806 | ! lat : nav_lat |
---|
1807 | ! Here we can check if we find the substring as the |
---|
1808 | ! names are more specific. |
---|
1809 | !- |
---|
1810 | SELECTCASE(axtype) |
---|
1811 | CASE ('lon') |
---|
1812 | dimname = 'nav_lon lon longitude' |
---|
1813 | found_rule = .TRUE. |
---|
1814 | CASE('lat') |
---|
1815 | dimname = 'nav_lat lat latitude' |
---|
1816 | found_rule = .TRUE. |
---|
1817 | CASE('lev') |
---|
1818 | dimname = 'plev level depth deptht' |
---|
1819 | found_rule = .TRUE. |
---|
1820 | CASE DEFAULT |
---|
1821 | found_rule = .FALSE. |
---|
1822 | END SELECT |
---|
1823 | !- |
---|
1824 | IF (found_rule) THEN |
---|
1825 | iv = 0 |
---|
1826 | DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
---|
1827 | iv = iv+1 |
---|
1828 | str1='' |
---|
1829 | iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, & |
---|
1830 | name=str1, ndims=ndim) |
---|
1831 | IF (INDEX(dimname,TRIM(str1)) >= 1) THEN |
---|
1832 | vid = iv |
---|
1833 | ENDIF |
---|
1834 | ENDDO |
---|
1835 | ENDIF |
---|
1836 | !- |
---|
1837 | ! Third rule : we find a variable with the same name as the dimension |
---|
1838 | ! lon = 1 |
---|
1839 | ! lat = 2 |
---|
1840 | ! lev = 3 |
---|
1841 | !- |
---|
1842 | IF (vid < 0) THEN |
---|
1843 | SELECTCASE(axtype) |
---|
1844 | CASE ('lon') |
---|
1845 | dimnb = 1 |
---|
1846 | found_rule = .TRUE. |
---|
1847 | CASE('lat') |
---|
1848 | dimnb = 2 |
---|
1849 | found_rule = .TRUE. |
---|
1850 | CASE('lev') |
---|
1851 | dimnb = 3 |
---|
1852 | found_rule = .TRUE. |
---|
1853 | CASE DEFAULT |
---|
1854 | found_rule = .FALSE. |
---|
1855 | END SELECT |
---|
1856 | !--- |
---|
1857 | IF (found_rule) THEN |
---|
1858 | iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname) |
---|
1859 | iv = 0 |
---|
1860 | DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
---|
1861 | iv = iv+1 |
---|
1862 | str1='' |
---|
1863 | iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, & |
---|
1864 | name=str1, ndims=ndim) |
---|
1865 | IF (INDEX(dimname,TRIM(str1)) == 1) THEN |
---|
1866 | vid = iv |
---|
1867 | ENDIF |
---|
1868 | ENDDO |
---|
1869 | ENDIF |
---|
1870 | ENDIF |
---|
1871 | !- |
---|
1872 | ! Stop the program if no coordinate was found |
---|
1873 | !- |
---|
1874 | IF (vid < 0) THEN |
---|
1875 | CALL histerr (3,'flinfindcood', & |
---|
1876 | 'No coordinate axis was found in the file', & |
---|
1877 | 'The data in this file can not be used', axtype) |
---|
1878 | ENDIF |
---|
1879 | !-------------------------- |
---|
1880 | END SUBROUTINE flinfindcood |
---|
1881 | !- |
---|
1882 | !=== |
---|
1883 | !- |
---|
1884 | SUBROUTINE flinclo (fid_in) |
---|
1885 | !--------------------------------------------------------------------- |
---|
1886 | IMPLICIT NONE |
---|
1887 | !- |
---|
1888 | INTEGER :: fid_in |
---|
1889 | !- |
---|
1890 | INTEGER :: iret |
---|
1891 | !--------------------------------------------------------------------- |
---|
1892 | iret = NF90_CLOSE (ncids(fid_in)) |
---|
1893 | ncfileopen(fid_in) = .FALSE. |
---|
1894 | !--------------------- |
---|
1895 | END SUBROUTINE flinclo |
---|
1896 | !- |
---|
1897 | !=== |
---|
1898 | !- |
---|
1899 | SUBROUTINE flinquery_var(fid_in, varname, exists) |
---|
1900 | !--------------------------------------------------------------------- |
---|
1901 | !- Queries the existance of a variable in the file. |
---|
1902 | !--------------------------------------------------------------------- |
---|
1903 | IMPLICIT NONE |
---|
1904 | !- |
---|
1905 | INTEGER :: fid_in |
---|
1906 | CHARACTER(LEN=*) varname |
---|
1907 | LOGICAL :: exists |
---|
1908 | !- |
---|
1909 | INTEGER :: iret, fid, vid |
---|
1910 | !--------------------------------------------------------------------- |
---|
1911 | fid = ncids(fid_in) |
---|
1912 | vid = -1 |
---|
1913 | iret = NF90_INQ_VARID (fid, varname, vid) |
---|
1914 | !- |
---|
1915 | exists = ( (vid >= 0).AND.(iret == NF90_NOERR) ) |
---|
1916 | !--------------------------- |
---|
1917 | END SUBROUTINE flinquery_var |
---|
1918 | !- |
---|
1919 | !=== |
---|
1920 | !- |
---|
1921 | SUBROUTINE flininspect (fid_in) |
---|
1922 | !--------------------------------------------------------------------- |
---|
1923 | IMPLICIT NONE |
---|
1924 | !- |
---|
1925 | ! fid : File id to inspect |
---|
1926 | !- |
---|
1927 | INTEGER :: fid_in |
---|
1928 | !- |
---|
1929 | !- LOCAL |
---|
1930 | !- |
---|
1931 | INTEGER :: iim, jjm, llm, ttm, fid_out |
---|
1932 | INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim |
---|
1933 | INTEGER :: iv, in, lll |
---|
1934 | INTEGER :: xid, yid, zid, tid |
---|
1935 | INTEGER,DIMENSION(NF90_MAX_VAR_DIMS) :: idimid |
---|
1936 | CHARACTER(LEN=80) :: name |
---|
1937 | CHARACTER(LEN=30) :: axname |
---|
1938 | !--------------------------------------------------------------------- |
---|
1939 | fid = ncids(fid_in) |
---|
1940 | !- |
---|
1941 | iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, & |
---|
1942 | nAttributes=nb_atts, unlimitedDimId=id_unlim) |
---|
1943 | !- |
---|
1944 | WRITE (*,*) 'IOIPSL ID : ',fid_in |
---|
1945 | WRITE (*,*) 'NetCDF ID : ',fid |
---|
1946 | WRITE (*,*) 'Number of dimensions : ',ndims |
---|
1947 | WRITE (*,*) 'Number of variables : ',nvars |
---|
1948 | WRITE (*,*) 'Number of global attributes : ',nb_atts |
---|
1949 | WRITE (*,*) 'ID unlimited : ',id_unlim |
---|
1950 | !- |
---|
1951 | xid = -1; iim = 0; |
---|
1952 | yid = -1; jjm = 0; |
---|
1953 | zid = -1; llm = 0; |
---|
1954 | tid = -1; ttm = 0; |
---|
1955 | !- |
---|
1956 | DO iv=1,ndims |
---|
1957 | !--- |
---|
1958 | iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll) |
---|
1959 | CALL strlowercase (axname) |
---|
1960 | axname = ADJUSTL(axname) |
---|
1961 | !--- |
---|
1962 | WRITE (*,*) 'Dimension number : ',iv |
---|
1963 | WRITE (*,*) 'Dimension name : ',TRIM(axname) |
---|
1964 | !--- |
---|
1965 | IF ( (INDEX(axname,'x') == 1) & |
---|
1966 | .OR.(INDEX(axname,'lon') == 1)) THEN |
---|
1967 | xid = iv; iim = lll; |
---|
1968 | WRITE (*,*) 'Dimension X size : ',iim |
---|
1969 | ELSE IF ( (INDEX(axname,'y') == 1) & |
---|
1970 | .OR.(INDEX(axname,'lat') == 1)) THEN |
---|
1971 | yid = iv; jjm = lll; |
---|
1972 | WRITE (*,*) 'Dimension Y size : ',jjm |
---|
1973 | ELSE IF ( (INDEX(axname,'lev') == 1) & |
---|
1974 | .OR.(INDEX(axname,'plev') == 1) & |
---|
1975 | .OR.(INDEX(axname,'z') == 1) & |
---|
1976 | .OR.(INDEX(axname,'depth') == 1)) THEN |
---|
1977 | zid = iv; llm = lll; |
---|
1978 | WRITE (*,*) 'Dimension Z size : ',llm |
---|
1979 | ELSE IF ( (INDEX(axname,'tstep') == 1) & |
---|
1980 | .OR.(INDEX(axname,'time') == 1) & |
---|
1981 | .OR.(INDEX(axname,'time_counter') == 1)) THEN |
---|
1982 | !---- For the time we certainly need to allow for other names |
---|
1983 | tid = iv; ttm = lll; |
---|
1984 | ELSE IF (ndims == 1) THEN |
---|
1985 | !---- Nothing was found and ndims=1 then we have a vector of data |
---|
1986 | xid = 1; iim = lll; |
---|
1987 | ENDIF |
---|
1988 | !--- |
---|
1989 | ENDDO |
---|
1990 | !- |
---|
1991 | ! Keep all this information |
---|
1992 | !- |
---|
1993 | nbfiles = nbfiles+1 |
---|
1994 | !- |
---|
1995 | IF (nbfiles > nbfile_max) THEN |
---|
1996 | CALL histerr(3,'flininspect', & |
---|
1997 | 'Too many files. Please increase nbfil_max', & |
---|
1998 | 'in program flincom.F90.',' ') |
---|
1999 | ENDIF |
---|
2000 | !- |
---|
2001 | ncids(nbfiles) = fid |
---|
2002 | ncnbd(nbfiles) = ndims |
---|
2003 | !- |
---|
2004 | ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /) |
---|
2005 | !- |
---|
2006 | ncfunli(nbfiles) = id_unlim |
---|
2007 | ncnba(nbfiles) = nb_atts |
---|
2008 | ncnbva(nbfiles) = nvars |
---|
2009 | ncfileopen(nbfiles) = .TRUE. |
---|
2010 | !- |
---|
2011 | fid_out = nbfiles |
---|
2012 | !- |
---|
2013 | DO in=1,nvars |
---|
2014 | iret = NF90_INQUIRE_VARIABLE (fid, in, & |
---|
2015 | name=name, ndims=ndims, dimids=idimid, nAtts=nb_atts) |
---|
2016 | WRITE (*,*) 'Variable number ------------ > ', in |
---|
2017 | WRITE (*,*) 'Variable name : ', TRIM(name) |
---|
2018 | WRITE (*,*) 'Number of dimensions : ', ndims |
---|
2019 | WRITE (*,*) 'Dimensions ID''s : ', idimid(1:ndims) |
---|
2020 | WRITE (*,*) 'Number of attributes : ', nb_atts |
---|
2021 | ENDDO |
---|
2022 | !------------------------- |
---|
2023 | END SUBROUTINE flininspect |
---|
2024 | !- |
---|
2025 | !=== |
---|
2026 | !- |
---|
2027 | END MODULE flincom |
---|