1 |
guez |
30 |
MODULE flincom |
2 |
|
|
|
3 |
|
|
! From flincom.f90, version 2.2 2006/03/07 09:21:51 |
4 |
|
|
|
5 |
|
|
IMPLICIT NONE |
6 |
|
|
|
7 |
|
|
PRIVATE |
8 |
guez |
32 |
PUBLIC flinclo, flinopen_nozoom, flininfo, ncids |
9 |
guez |
30 |
|
10 |
guez |
32 |
! This is the data we keep on each file we open: |
11 |
|
|
INTEGER, PARAMETER:: nbfile_max = 200 |
12 |
|
|
INTEGER, SAVE:: ncids(nbfile_max) |
13 |
|
|
INTEGER, SAVE:: ncnbva(nbfile_max) |
14 |
|
|
LOGICAL, SAVE:: ncfileopen(nbfile_max)=.FALSE. |
15 |
guez |
30 |
|
16 |
|
|
CONTAINS |
17 |
|
|
|
18 |
guez |
32 |
SUBROUTINE flinopen_nozoom(iim, jjm, llm, lon, lat, lev, & |
19 |
guez |
30 |
ttm, itaus, date0, dt, fid_out) |
20 |
|
|
|
21 |
|
|
!- The routine will open an input file |
22 |
|
|
!- INPUT |
23 |
|
|
!- There is no test of the content of the file against the input |
24 |
|
|
! from the model |
25 |
|
|
|
26 |
|
|
!- iim : size in the x direction in the file (longitude) |
27 |
|
|
!- jjm : size in the y direction |
28 |
|
|
!- llm : number of levels |
29 |
|
|
!- (llm = 0 means no axis to be expected) |
30 |
|
|
|
31 |
|
|
!- WARNING : |
32 |
|
|
!- It is for the user to check |
33 |
|
|
!- that the dimensions of lon lat and lev are correct when passed to |
34 |
|
|
!- flinopen. This can be done after the call when iim and jjm have |
35 |
|
|
!- been retrieved from the netCDF file. In F90 this problem will |
36 |
|
|
!- be solved with an internal assign |
37 |
|
|
!- IF iim, jjm, llm or ttm are parameters in the calling program |
38 |
|
|
!- it will create a segmentation fault |
39 |
|
|
|
40 |
|
|
!- ttm : size of time axis |
41 |
|
|
|
42 |
|
|
!- OUTPUT |
43 |
|
|
|
44 |
|
|
!- lon : array of (iim,jjm), |
45 |
|
|
!- that contains the longitude of each point |
46 |
|
|
!- lat : same for latitude |
47 |
|
|
!- lev : An array of llm for the latitude |
48 |
|
|
!- itaus : Time steps within this file |
49 |
|
|
!- date0 : Julian date at which itau = 0 |
50 |
|
|
!- dt : length of the time steps of the data |
51 |
|
|
|
52 |
|
|
!--------------------------------------------------------------------- |
53 |
|
|
|
54 |
guez |
32 |
USE calendar, ONLY : ymds2ju, ioconf_calendar |
55 |
|
|
USE errioipsl, ONLY : histerr |
56 |
|
|
USE netcdf, ONLY : nf90_get_att, nf90_get_var, nf90_global, & |
57 |
|
|
nf90_inquire_variable |
58 |
|
|
|
59 |
guez |
30 |
IMPLICIT NONE |
60 |
|
|
|
61 |
|
|
! ARGUMENTS |
62 |
|
|
|
63 |
|
|
INTEGER, intent(in) :: iim, jjm, llm, ttm |
64 |
|
|
real, intent(out):: lon(iim,jjm), lat(iim,jjm), lev(llm) |
65 |
|
|
INTEGER, intent(out):: itaus(ttm) |
66 |
|
|
REAL, intent(out):: date0, dt |
67 |
|
|
|
68 |
|
|
INTEGER, intent(in):: fid_out |
69 |
|
|
! (file ID which is later used to read the data) |
70 |
|
|
|
71 |
|
|
! LOCAL |
72 |
|
|
|
73 |
|
|
INTEGER :: iret, vid, fid, nbdim, i |
74 |
|
|
INTEGER :: gdtt_id, old_id, iv, gdtmaf_id |
75 |
|
|
CHARACTER(LEN=250) :: name |
76 |
guez |
32 |
CHARACTER(LEN=80) :: units, my_calendar |
77 |
guez |
30 |
INTEGER :: year, month, day |
78 |
|
|
REAL :: r_year, r_month, r_day |
79 |
|
|
INTEGER :: year0, month0, day0, hours0, minutes0, seci |
80 |
|
|
REAL :: sec, sec0 |
81 |
|
|
CHARACTER :: strc |
82 |
|
|
|
83 |
|
|
REAL,DIMENSION(:),ALLOCATABLE :: vec_tmp |
84 |
|
|
|
85 |
|
|
!--------------------------------------------------------------------- |
86 |
|
|
|
87 |
|
|
IF ( (fid_out < 1).OR.(fid_out > nbfile_max) ) THEN |
88 |
|
|
! Either the fid_out has not been initialized (0 or very large) |
89 |
|
|
! then we have to open anyway. Else we only need to open the file |
90 |
|
|
! if it has not been opened before. |
91 |
|
|
print *, "Call flinfo before flinopen" |
92 |
|
|
stop 1 |
93 |
|
|
end IF |
94 |
|
|
IF (.NOT.ncfileopen(fid_out)) THEN |
95 |
|
|
print *, "Call flinfo before flinopen" |
96 |
|
|
stop 1 |
97 |
|
|
end IF |
98 |
|
|
|
99 |
|
|
!-- The user has already opened the file |
100 |
|
|
!-- and we trust that he knows the dimensions |
101 |
|
|
|
102 |
|
|
fid = ncids(fid_out) |
103 |
|
|
|
104 |
|
|
! 2.0 get the sizes and names of the different coordinates |
105 |
|
|
! and do a first set of verification. |
106 |
|
|
|
107 |
|
|
! 3.0 Check if we are realy talking about the same coodinate system |
108 |
|
|
! if not then we get the lon, lat and lev variables from the file |
109 |
|
|
|
110 |
|
|
!-- 4.0 extracting the coordinates |
111 |
|
|
!--- |
112 |
|
|
CALL flinfindcood (fid_out, 'lon', vid, nbdim) |
113 |
|
|
IF (nbdim == 2) THEN |
114 |
|
|
iret = NF90_GET_VAR (fid, vid, lon, & |
115 |
|
|
start=(/ 1, 1 /), count=(/ iim, jjm /)) |
116 |
|
|
ELSE |
117 |
|
|
ALLOCATE(vec_tmp(iim)) |
118 |
|
|
iret = NF90_GET_VAR (fid, vid, vec_tmp, & |
119 |
|
|
start=(/ 1 /), count=(/ iim /)) |
120 |
|
|
DO i=1,jjm |
121 |
|
|
lon(:,i) = vec_tmp(:) |
122 |
|
|
ENDDO |
123 |
|
|
DEALLOCATE(vec_tmp) |
124 |
|
|
ENDIF |
125 |
|
|
!--- |
126 |
|
|
CALL flinfindcood (fid_out, 'lat', vid, nbdim) |
127 |
|
|
IF (nbdim == 2) THEN |
128 |
|
|
iret = NF90_GET_VAR (fid, vid, lat, & |
129 |
|
|
start=(/ 1, 1 /), count=(/ iim, jjm /)) |
130 |
|
|
ELSE |
131 |
|
|
ALLOCATE(vec_tmp(jjm)) |
132 |
|
|
iret = NF90_GET_VAR (fid, vid, vec_tmp, & |
133 |
|
|
start=(/ 1 /), count=(/ jjm /)) |
134 |
|
|
DO i=1,iim |
135 |
|
|
lat(i,:) = vec_tmp(:) |
136 |
|
|
ENDDO |
137 |
|
|
DEALLOCATE(vec_tmp) |
138 |
|
|
ENDIF |
139 |
|
|
!--- |
140 |
|
|
IF (llm > 0) THEN |
141 |
|
|
CALL flinfindcood (fid_out, 'lev', vid, nbdim) |
142 |
|
|
IF (nbdim == 1) THEN |
143 |
|
|
iret = NF90_GET_VAR (fid, vid, lev, & |
144 |
|
|
start=(/ 1 /), count=(/ llm /)) |
145 |
|
|
ELSE |
146 |
|
|
CALL histerr (3,'flinopen', & |
147 |
|
|
'Can not handle vertical coordinates that have more',& |
148 |
|
|
'than 1 dimension',' ') |
149 |
|
|
ENDIF |
150 |
|
|
ENDIF |
151 |
|
|
|
152 |
|
|
! 5.0 Get all the details for the time if possible needed |
153 |
|
|
|
154 |
|
|
IF (ttm > 0) THEN |
155 |
|
|
!--- |
156 |
|
|
!-- 5.1 Find the time axis. Prefered method is the 'timestep since' |
157 |
|
|
!--- |
158 |
|
|
gdtmaf_id = -1 |
159 |
|
|
gdtt_id = -1 |
160 |
|
|
old_id = -1 |
161 |
|
|
DO iv=1,ncnbva(fid_out) |
162 |
|
|
name='' |
163 |
|
|
iret = NF90_INQUIRE_VARIABLE (fid, iv, name=name) |
164 |
|
|
units='' |
165 |
|
|
iret = NF90_GET_ATT (fid, iv, 'units', units) |
166 |
|
|
IF (INDEX(units,'seconds since') > 0) gdtmaf_id = iv |
167 |
|
|
IF (INDEX(units,'timesteps since') > 0) gdtt_id = iv |
168 |
|
|
IF (INDEX(name, 'tstep') > 0) old_id = iv |
169 |
|
|
ENDDO |
170 |
|
|
!--- |
171 |
|
|
IF (gdtt_id > 0) THEN |
172 |
|
|
vid = gdtt_id |
173 |
|
|
ELSE IF (gdtmaf_id > 0) THEN |
174 |
|
|
vid = gdtmaf_id |
175 |
|
|
ELSE IF (old_id > 0) THEN |
176 |
|
|
vid = old_id |
177 |
|
|
ELSE |
178 |
|
|
CALL histerr (3, 'flinopen', 'No time axis found',' ',' ') |
179 |
|
|
ENDIF |
180 |
|
|
!--- |
181 |
|
|
ALLOCATE(vec_tmp(ttm)) |
182 |
|
|
iret = NF90_GET_VAR (fid, vid, vec_tmp, & |
183 |
|
|
start=(/ 1 /), count=(/ ttm /)) |
184 |
|
|
itaus(1:ttm) = NINT(vec_tmp(1:ttm)) |
185 |
|
|
DEALLOCATE(vec_tmp) |
186 |
|
|
!--- |
187 |
|
|
!-- Getting all the details for the time axis |
188 |
|
|
!--- |
189 |
|
|
!-- Find the calendar |
190 |
guez |
32 |
my_calendar='XXXX' |
191 |
|
|
iret = NF90_GET_ATT (fid, gdtmaf_id, 'calendar', my_calendar) |
192 |
|
|
IF ( INDEX(my_calendar,'XXXX') < 1 ) THEN |
193 |
|
|
CALL ioconf_calendar(my_calendar) |
194 |
guez |
30 |
ENDIF |
195 |
|
|
!-- |
196 |
|
|
units = '' |
197 |
|
|
iret = NF90_GET_ATT (fid, vid, 'units', units) |
198 |
|
|
IF (gdtt_id > 0) THEN |
199 |
|
|
units = units(INDEX(units,'since')+6:LEN_TRIM(units)) |
200 |
|
|
READ (units,'(I4.4,5(a,I2.2))') & |
201 |
|
|
year0, strc, month0, strc, day0, & |
202 |
|
|
strc, hours0, strc, minutes0, strc, seci |
203 |
|
|
sec0 = hours0*3600. + minutes0*60. + seci |
204 |
|
|
CALL ymds2ju (year0, month0, day0, sec0, date0) |
205 |
|
|
iret = NF90_GET_ATT (fid, gdtt_id, 'tstep_sec', dt) |
206 |
|
|
ELSE IF (gdtmaf_id > 0) THEN |
207 |
|
|
units = units(INDEX(units,'since')+6:LEN_TRIM(units)) |
208 |
|
|
READ (units,'(I4.4,5(a,I2.2))') & |
209 |
|
|
year0, strc, month0, strc, day0, & |
210 |
|
|
strc, hours0, strc, minutes0, strc, seci |
211 |
|
|
sec0 = hours0*3600. + minutes0*60. + seci |
212 |
|
|
CALL ymds2ju (year0, month0, day0, sec0, date0) |
213 |
|
|
!----- |
214 |
|
|
ELSE IF (old_id > 0) THEN |
215 |
|
|
iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'delta_tstep_sec', dt) |
216 |
|
|
iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'day0', r_day) |
217 |
|
|
iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'sec0', sec) |
218 |
|
|
iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'year0', r_year) |
219 |
|
|
iret = NF90_GET_ATT (fid, NF90_GLOBAL, 'month0', r_month) |
220 |
|
|
!----- |
221 |
|
|
day = INT(r_day) |
222 |
|
|
month = INT(r_month) |
223 |
|
|
year = INT(r_year) |
224 |
|
|
!----- |
225 |
|
|
CALL ymds2ju (year, month, day, sec, date0) |
226 |
|
|
ENDIF |
227 |
|
|
ENDIF |
228 |
|
|
|
229 |
|
|
END SUBROUTINE flinopen_nozoom |
230 |
|
|
|
231 |
|
|
!=== |
232 |
|
|
|
233 |
|
|
SUBROUTINE flininfo (filename, iim, jjm, llm, ttm, fid_out) |
234 |
|
|
!--------------------------------------------------------------------- |
235 |
|
|
!- This subroutine allows to get some information. |
236 |
|
|
!- It is usualy done within flinopen but the user may want to call |
237 |
|
|
!- it before in order to allocate the space needed to extract the |
238 |
|
|
!- data from the file. |
239 |
|
|
!--------------------------------------------------------------------- |
240 |
guez |
32 |
USE strlowercase_m, ONLY : strlowercase |
241 |
|
|
USE errioipsl, ONLY : histerr |
242 |
|
|
USE netcdf, ONLY : nf90_inquire, nf90_inquire_dimension, nf90_noerr, & |
243 |
|
|
nf90_nowrite, nf90_open |
244 |
|
|
|
245 |
guez |
30 |
IMPLICIT NONE |
246 |
|
|
|
247 |
|
|
! ARGUMENTS |
248 |
|
|
|
249 |
|
|
CHARACTER(LEN=*), intent(in):: filename |
250 |
|
|
INTEGER, intent(out):: iim, jjm, llm, ttm, fid_out |
251 |
|
|
|
252 |
|
|
! LOCAL |
253 |
|
|
|
254 |
guez |
32 |
INTEGER, SAVE :: nbfiles = 0 |
255 |
|
|
INTEGER, SAVE :: ncdims(nbfile_max,4) |
256 |
guez |
30 |
INTEGER :: iret, fid, ndims, nvars, nb_atts, id_unlim |
257 |
|
|
INTEGER :: iv, lll |
258 |
|
|
CHARACTER(LEN=80) :: name |
259 |
|
|
CHARACTER(LEN=30) :: axname |
260 |
|
|
|
261 |
|
|
LOGICAL :: check = .FALSE. |
262 |
|
|
!--------------------------------------------------------------------- |
263 |
|
|
lll = LEN_TRIM(filename) |
264 |
|
|
IF (filename(lll-2:lll) /= '.nc') THEN |
265 |
|
|
name = filename(1:lll)//'.nc' |
266 |
|
|
ELSE |
267 |
|
|
name = filename(1:lll) |
268 |
|
|
ENDIF |
269 |
|
|
|
270 |
|
|
iret = NF90_OPEN (name, NF90_NOWRITE, fid) |
271 |
|
|
IF (iret /= NF90_NOERR) THEN |
272 |
|
|
CALL histerr(3, 'flininfo','Could not open file :',TRIM(name),' ') |
273 |
|
|
ENDIF |
274 |
|
|
|
275 |
|
|
iret = NF90_INQUIRE (fid, nDimensions=ndims, nVariables=nvars, & |
276 |
|
|
nAttributes=nb_atts, unlimitedDimId=id_unlim) |
277 |
|
|
|
278 |
|
|
iim = 0; |
279 |
|
|
jjm = 0; |
280 |
|
|
llm = 0; |
281 |
|
|
ttm = 0; |
282 |
|
|
|
283 |
|
|
DO iv=1,ndims |
284 |
|
|
!--- |
285 |
|
|
iret = NF90_INQUIRE_DIMENSION (fid, iv, name=axname, len=lll) |
286 |
|
|
CALL strlowercase (axname) |
287 |
|
|
axname = ADJUSTL(axname) |
288 |
|
|
!--- |
289 |
|
|
IF (check) WRITE(*,*) & |
290 |
|
|
'flininfo - getting axname',iv,axname,lll |
291 |
|
|
!--- |
292 |
|
|
IF ( (INDEX(axname,'x') == 1) & |
293 |
|
|
.OR.(INDEX(axname,'lon') == 1) ) THEN |
294 |
|
|
iim = lll; |
295 |
|
|
ELSE IF ( (INDEX(axname,'y') == 1) & |
296 |
|
|
.OR.(INDEX(axname,'lat') == 1) ) THEN |
297 |
|
|
jjm = lll; |
298 |
|
|
ELSE IF ( (INDEX(axname,'lev') == 1) & |
299 |
|
|
.OR.(INDEX(axname,'plev') == 1) & |
300 |
|
|
.OR.(INDEX(axname,'z') == 1) & |
301 |
|
|
.OR.(INDEX(axname,'depth') == 1) ) THEN |
302 |
|
|
llm = lll; |
303 |
|
|
ELSE IF ( (INDEX(axname,'tstep') == 1) & |
304 |
|
|
.OR.(INDEX(axname,'time_counter') == 1) ) THEN |
305 |
|
|
!---- For the time we certainly need to allow for other names |
306 |
|
|
ttm = lll; |
307 |
|
|
ELSE IF (ndims == 1) THEN |
308 |
|
|
!---- Nothing was found and ndims=1 then we have a vector of data |
309 |
|
|
iim = lll; |
310 |
|
|
ENDIF |
311 |
|
|
!--- |
312 |
|
|
ENDDO |
313 |
|
|
|
314 |
|
|
! Keep all this information |
315 |
|
|
|
316 |
|
|
nbfiles = nbfiles+1 |
317 |
|
|
|
318 |
|
|
IF (nbfiles > nbfile_max) THEN |
319 |
|
|
CALL histerr (3,'flininfo', & |
320 |
|
|
'Too many files. Please increase nbfil_max', & |
321 |
|
|
'in program flincom.F90.',' ') |
322 |
|
|
ENDIF |
323 |
|
|
|
324 |
|
|
ncids(nbfiles) = fid |
325 |
|
|
ncdims(nbfiles,1:4) = (/ iim, jjm, llm, ttm /) |
326 |
|
|
|
327 |
|
|
ncnbva(nbfiles) = nvars |
328 |
|
|
ncfileopen(nbfiles) = .TRUE. |
329 |
|
|
|
330 |
|
|
fid_out = nbfiles |
331 |
|
|
!---------------------- |
332 |
|
|
END SUBROUTINE flininfo |
333 |
|
|
|
334 |
|
|
!=== |
335 |
|
|
|
336 |
|
|
SUBROUTINE flinfindcood (fid_in, axtype, vid, ndim) |
337 |
|
|
!--------------------------------------------------------------------- |
338 |
|
|
!- This subroutine explores the file in order to find |
339 |
|
|
!- the coordinate according to a number of rules |
340 |
|
|
!--------------------------------------------------------------------- |
341 |
guez |
32 |
USE strlowercase_m, ONLY : strlowercase |
342 |
|
|
USE errioipsl, ONLY : histerr |
343 |
|
|
USE netcdf, ONLY : nf90_get_att, nf90_inquire_dimension, & |
344 |
|
|
nf90_inquire_variable, nf90_noerr |
345 |
|
|
|
346 |
guez |
30 |
IMPLICIT NONE |
347 |
|
|
|
348 |
|
|
! ARGUMENTS |
349 |
|
|
|
350 |
|
|
INTEGER, intent(in):: fid_in |
351 |
|
|
integer vid, ndim |
352 |
|
|
CHARACTER(LEN=3) :: axtype |
353 |
|
|
|
354 |
|
|
! LOCAL |
355 |
|
|
|
356 |
|
|
INTEGER :: iv, iret, dimnb |
357 |
|
|
CHARACTER(LEN=40) :: dimname, dimuni1, dimuni2, dimuni3 |
358 |
|
|
CHARACTER(LEN=30) :: str1 |
359 |
|
|
LOGICAL :: found_rule = .FALSE. |
360 |
|
|
!--------------------------------------------------------------------- |
361 |
|
|
vid = -1 |
362 |
|
|
|
363 |
|
|
! Make sure all strings are invalid |
364 |
|
|
|
365 |
|
|
dimname = '?-?' |
366 |
|
|
dimuni1 = '?-?' |
367 |
|
|
dimuni2 = '?-?' |
368 |
|
|
dimuni3 = '?-?' |
369 |
|
|
|
370 |
|
|
! First rule : we look for the correct units |
371 |
|
|
! lon : east |
372 |
|
|
! lat : north |
373 |
|
|
! We make an exact check as it would be too easy to mistake |
374 |
|
|
! some units by just comparing the substrings. |
375 |
|
|
|
376 |
|
|
SELECTCASE(axtype) |
377 |
|
|
CASE ('lon') |
378 |
|
|
dimuni1 = 'degree_e' |
379 |
|
|
dimuni2 = 'degrees_e' |
380 |
|
|
found_rule = .TRUE. |
381 |
|
|
CASE('lat') |
382 |
|
|
dimuni1 = 'degree_n' |
383 |
|
|
dimuni2 = 'degrees_n' |
384 |
|
|
found_rule = .TRUE. |
385 |
|
|
CASE('lev') |
386 |
|
|
dimuni1 = 'm' |
387 |
|
|
dimuni2 = 'km' |
388 |
|
|
dimuni3 = 'hpa' |
389 |
|
|
found_rule = .TRUE. |
390 |
|
|
CASE DEFAULT |
391 |
|
|
found_rule = .FALSE. |
392 |
|
|
END SELECT |
393 |
|
|
|
394 |
|
|
IF (found_rule) THEN |
395 |
|
|
iv = 0 |
396 |
|
|
DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
397 |
|
|
iv = iv+1 |
398 |
|
|
str1 = '' |
399 |
|
|
iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1) |
400 |
|
|
IF (iret == NF90_NOERR) THEN |
401 |
|
|
CALL strlowercase (str1) |
402 |
|
|
IF ( (INDEX(str1, TRIM(dimuni1)) == 1) & |
403 |
|
|
.OR.(INDEX(str1, TRIM(dimuni2)) == 1) & |
404 |
|
|
.OR.(INDEX(str1, TRIM(dimuni3)) == 1) ) THEN |
405 |
|
|
vid = iv |
406 |
|
|
iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim) |
407 |
|
|
ENDIF |
408 |
|
|
ENDIF |
409 |
|
|
ENDDO |
410 |
|
|
ENDIF |
411 |
|
|
|
412 |
|
|
! Second rule : we find specific names : |
413 |
|
|
! lon : nav_lon |
414 |
|
|
! lat : nav_lat |
415 |
|
|
! Here we can check if we find the substring as the |
416 |
|
|
! names are more specific. |
417 |
|
|
|
418 |
|
|
SELECTCASE(axtype) |
419 |
|
|
CASE ('lon') |
420 |
|
|
dimname = 'nav_lon lon longitude' |
421 |
|
|
found_rule = .TRUE. |
422 |
|
|
CASE('lat') |
423 |
|
|
dimname = 'nav_lat lat latitude' |
424 |
|
|
found_rule = .TRUE. |
425 |
|
|
CASE('lev') |
426 |
|
|
dimname = 'plev level depth deptht' |
427 |
|
|
found_rule = .TRUE. |
428 |
|
|
CASE DEFAULT |
429 |
|
|
found_rule = .FALSE. |
430 |
|
|
END SELECT |
431 |
|
|
|
432 |
|
|
IF (found_rule) THEN |
433 |
|
|
iv = 0 |
434 |
|
|
DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
435 |
|
|
iv = iv+1 |
436 |
|
|
str1='' |
437 |
|
|
iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, & |
438 |
|
|
name=str1, ndims=ndim) |
439 |
|
|
IF (INDEX(dimname,TRIM(str1)) >= 1) THEN |
440 |
|
|
vid = iv |
441 |
|
|
ENDIF |
442 |
|
|
ENDDO |
443 |
|
|
ENDIF |
444 |
|
|
|
445 |
|
|
! Third rule : we find a variable with the same name as the dimension |
446 |
|
|
! lon = 1 |
447 |
|
|
! lat = 2 |
448 |
|
|
! lev = 3 |
449 |
|
|
|
450 |
|
|
IF (vid < 0) THEN |
451 |
|
|
SELECTCASE(axtype) |
452 |
|
|
CASE ('lon') |
453 |
|
|
dimnb = 1 |
454 |
|
|
found_rule = .TRUE. |
455 |
|
|
CASE('lat') |
456 |
|
|
dimnb = 2 |
457 |
|
|
found_rule = .TRUE. |
458 |
|
|
CASE('lev') |
459 |
|
|
dimnb = 3 |
460 |
|
|
found_rule = .TRUE. |
461 |
|
|
CASE DEFAULT |
462 |
|
|
found_rule = .FALSE. |
463 |
|
|
END SELECT |
464 |
|
|
!--- |
465 |
|
|
IF (found_rule) THEN |
466 |
|
|
iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname) |
467 |
|
|
iv = 0 |
468 |
|
|
DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) ) |
469 |
|
|
iv = iv+1 |
470 |
|
|
str1='' |
471 |
|
|
iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, & |
472 |
|
|
name=str1, ndims=ndim) |
473 |
|
|
IF (INDEX(dimname,TRIM(str1)) == 1) THEN |
474 |
|
|
vid = iv |
475 |
|
|
ENDIF |
476 |
|
|
ENDDO |
477 |
|
|
ENDIF |
478 |
|
|
ENDIF |
479 |
|
|
|
480 |
|
|
! Stop the program if no coordinate was found |
481 |
|
|
|
482 |
|
|
IF (vid < 0) THEN |
483 |
|
|
CALL histerr (3,'flinfindcood', & |
484 |
|
|
'No coordinate axis was found in the file', & |
485 |
|
|
'The data in this file can not be used', axtype) |
486 |
|
|
ENDIF |
487 |
|
|
!-------------------------- |
488 |
|
|
END SUBROUTINE flinfindcood |
489 |
|
|
|
490 |
|
|
!=== |
491 |
|
|
|
492 |
|
|
SUBROUTINE flinclo (fid_in) |
493 |
|
|
!--------------------------------------------------------------------- |
494 |
guez |
32 |
USE netcdf, ONLY : nf90_close |
495 |
guez |
30 |
IMPLICIT NONE |
496 |
|
|
|
497 |
|
|
INTEGER :: fid_in |
498 |
|
|
|
499 |
|
|
INTEGER :: iret |
500 |
|
|
!--------------------------------------------------------------------- |
501 |
|
|
iret = NF90_CLOSE (ncids(fid_in)) |
502 |
|
|
ncfileopen(fid_in) = .FALSE. |
503 |
|
|
!--------------------- |
504 |
|
|
END SUBROUTINE flinclo |
505 |
|
|
|
506 |
|
|
END MODULE flincom |