1 | dnl $Id :$
|
---|
2 | dnl
|
---|
3 | dnl This file is used for the automatic generation of the file
|
---|
4 | dnl ncdf_getvar.f90.
|
---|
5 | dnl
|
---|
6 | dnl AUTHOR
|
---|
7 | dnl C. Marquardt, Darmstadt, Germany <christian@marquardt.sc>
|
---|
8 | dnl
|
---|
9 | dnl COPYRIGHT
|
---|
10 | dnl
|
---|
11 | dnl Copyright (c) 2005 Christian Marquardt <christian@marquardt.sc>
|
---|
12 | dnl
|
---|
13 | dnl All rights reserved.
|
---|
14 | dnl
|
---|
15 | dnl Permission is hereby granted, free of charge, to any person obtaining
|
---|
16 | dnl a copy of this software and associated documentation files (the
|
---|
17 | dnl "Software"), to deal in the Software without restriction, including
|
---|
18 | dnl without limitation the rights to use, copy, modify, merge, publish,
|
---|
19 | dnl distribute, sublicense, and/or sell copies of the Software, and to
|
---|
20 | dnl permit persons to whom the Software is furnished to do so, subject to
|
---|
21 | dnl the following conditions:
|
---|
22 | dnl
|
---|
23 | dnl The above copyright notice and this permission notice shall be
|
---|
24 | dnl included in all copies or substantial portions of the Software as well
|
---|
25 | dnl as in supporting documentation.
|
---|
26 | dnl
|
---|
27 | dnl THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
|
---|
28 | dnl EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
|
---|
29 | dnl MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
|
---|
30 | dnl NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE
|
---|
31 | dnl LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
|
---|
32 | dnl OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION
|
---|
33 | dnl WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
|
---|
34 | dnl
|
---|
35 | subroutine NCDF_AFUN (name, values, ncfile, ncid, rec, start, count, units, range)
|
---|
36 |
|
---|
37 | use typeSizes
|
---|
38 | use messages
|
---|
39 | use unitconvert, only: ut_convert
|
---|
40 | use ncdf, not_this => NCDF_AFUN
|
---|
41 |
|
---|
42 | implicit none
|
---|
43 |
|
---|
44 | character(len = *), intent( in) :: name
|
---|
45 | TYPE, dimension(COLONS), &
|
---|
46 | INTENT_OR_POINTER :: values
|
---|
47 | TYPE, dimension(2), optional :: range
|
---|
48 | character(len = *), optional :: ncfile
|
---|
49 | integer, optional :: ncid
|
---|
50 | integer, optional :: rec
|
---|
51 | integer, dimension(:), optional :: start
|
---|
52 | integer, dimension(:), optional :: count
|
---|
53 | character(len = *), optional :: units
|
---|
54 |
|
---|
55 | integer :: status, varid, ncid_local, groupid
|
---|
56 | integer :: ndims, dimrec, i, j
|
---|
57 | integer, dimension(NF90_MAX_VAR_DIMS) :: strt, cnts, dimids
|
---|
58 | character(len = NF90_MAX_NAME) :: ncdf_units
|
---|
59 | TYPE, dimension(COLONS), &
|
---|
60 | allocatable :: rvalues
|
---|
61 | TYPE, dimension(2) :: rrange
|
---|
62 | TYPE &
|
---|
63 | :: scale_factor, add_offset
|
---|
64 | logical :: have_scale, have_offset
|
---|
65 |
|
---|
66 | logical :: have_range
|
---|
67 | character(len = NF90_MAX_NAME) :: vname
|
---|
68 | logical :: havegroup
|
---|
69 |
|
---|
70 | ! See if this is the current netcdf
|
---|
71 | ! ---------------------------------
|
---|
72 |
|
---|
73 | if (present(ncfile)) then
|
---|
74 | if (ncfile == ncdf_ncname) then
|
---|
75 | ncid_local = ncdf_ncid
|
---|
76 | else
|
---|
77 | status = nf90_open(ncfile, nf90_share, ncid)
|
---|
78 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
---|
79 | endif
|
---|
80 |
|
---|
81 | else if (present(ncid)) then
|
---|
82 |
|
---|
83 | ncid_local = ncid
|
---|
84 |
|
---|
85 | else
|
---|
86 |
|
---|
87 | ncid_local = ncdf_ncid
|
---|
88 |
|
---|
89 | endif
|
---|
90 |
|
---|
91 | j = 0
|
---|
92 | ifelse(POINTER, `yes',,dnl
|
---|
93 | ! Allocate memory
|
---|
94 | ! ---------------
|
---|
95 |
|
---|
96 | cnts(:NUMDIMS) = shape(values)
|
---|
97 | allocate(rvalues(ALLOC_ARGS`'))
|
---|
98 | )
|
---|
99 |
|
---|
100 | ! Get group id if necessary and replace ncid_local with it
|
---|
101 | ! --------------------------------------------------------
|
---|
102 |
|
---|
103 | status = ncdf_getgroupid(ncid_local, name, vname, groupid, havegroup)
|
---|
104 | if (status /= nf90_noerr) then
|
---|
105 | WRITE ( *, FMT="(A)" ) "ERROR: Group ID not found: "// name
|
---|
106 | call ncdf_error_handler(status)
|
---|
107 | return
|
---|
108 | endif
|
---|
109 | ncid_local = groupid
|
---|
110 |
|
---|
111 | ! Get variable ID
|
---|
112 | ! ---------------
|
---|
113 |
|
---|
114 | status = nf90_inq_varid(ncid_local, vname, varid)
|
---|
115 | if (status /= nf90_noerr) then
|
---|
116 | WRITE ( *, FMT="(A)" ) "ERROR: Variable not found: "// name
|
---|
117 | call ncdf_error_handler(status)
|
---|
118 | endif
|
---|
119 |
|
---|
120 | ! Obtain some information about the variables dimensionality
|
---|
121 | ! ----------------------------------------------------------
|
---|
122 |
|
---|
123 | status = nf90_inquire_variable(ncid_local, varid, &
|
---|
124 | ndims = ndims, dimids = dimids)
|
---|
125 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
---|
126 |
|
---|
127 | ! Obtain scaling factors for the variable
|
---|
128 | ! ---------------------------------------
|
---|
129 |
|
---|
130 | status = nf90_get_att(ncid_local, varid, 'scale_factor', scale_factor)
|
---|
131 | if (status == nf90_enotatt) then
|
---|
132 | have_scale = .false.
|
---|
133 | status = nf90_noerr
|
---|
134 | else
|
---|
135 | have_scale = .true.
|
---|
136 | endif
|
---|
137 | if (status /= nf90_noerr) then
|
---|
138 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute scale_factor not found for variable: "// name
|
---|
139 | call ncdf_error_handler(status)
|
---|
140 | endif
|
---|
141 |
|
---|
142 | status = nf90_get_att(ncid_local, varid, 'add_offset', add_offset)
|
---|
143 | if (status == nf90_enotatt) then
|
---|
144 | have_offset = .false.
|
---|
145 | status = nf90_noerr
|
---|
146 | else
|
---|
147 | have_offset = .true.
|
---|
148 | endif
|
---|
149 | if (status /= nf90_noerr) then
|
---|
150 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute add_offset not found for variable: "// name
|
---|
151 | call ncdf_error_handler(status)
|
---|
152 | endif
|
---|
153 |
|
---|
154 | ! Prepare start and count arrays - these are the defaults
|
---|
155 | ! -------------------------------------------------------
|
---|
156 |
|
---|
157 | strt = 1
|
---|
158 | cnts = 0
|
---|
159 |
|
---|
160 | do i = 1, NUMDIMS
|
---|
161 | status = nf90_inquire_dimension(ncid_local, dimids(i), len = cnts(i))
|
---|
162 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
---|
163 | ifelse(POINTER, `yes',,dnl
|
---|
164 | cnts(i) = min(cnts(i), size(values, i))
|
---|
165 | )dnl
|
---|
166 | enddo
|
---|
167 |
|
---|
168 | ! Special cases: record is given...
|
---|
169 | ! ---------------------------------
|
---|
170 |
|
---|
171 | if (present(rec)) then
|
---|
172 |
|
---|
173 | ! ...see if an unlimited (record) dimension is available...
|
---|
174 |
|
---|
175 | status = nf90_inquire(ncid_local, unlimitedDimID = dimrec)
|
---|
176 | if (status /= nf90_noerr .or. dimrec == -1) &
|
---|
177 | call ncdf_error_handler(status)
|
---|
178 |
|
---|
179 | ! ...make sure this is true for the variable in question...
|
---|
180 |
|
---|
181 | i = 1
|
---|
182 | do while (i <= ndims)
|
---|
183 | if (dimids(i) == dimrec) exit
|
---|
184 | i = i + 1
|
---|
185 | enddo
|
---|
186 |
|
---|
187 | ! ...and set the start array.
|
---|
188 |
|
---|
189 | if (i <= ndims) then
|
---|
190 | strt(i) = rec
|
---|
191 | cnts(i) = 1
|
---|
192 | else
|
---|
193 | call ncdf_error_handler(NF90_ENORECVARS)
|
---|
194 | endif
|
---|
195 |
|
---|
196 | ! ...or start and count...
|
---|
197 |
|
---|
198 | else if (present(start) .and. present(count)) then
|
---|
199 |
|
---|
200 | strt(1:size(start)) = start
|
---|
201 | cnts(1:size(count)) = count
|
---|
202 |
|
---|
203 | ! ...or only start...
|
---|
204 |
|
---|
205 | else if (present(start)) then
|
---|
206 |
|
---|
207 | strt(1:size(start)) = start
|
---|
208 | cnts(1:ndims) = strt(1:ndims) - cnts(1:ndims) + 1
|
---|
209 |
|
---|
210 | ! ...or only count...
|
---|
211 |
|
---|
212 | else if (present(count)) then
|
---|
213 |
|
---|
214 | cnts(1:size(count)) = count
|
---|
215 |
|
---|
216 | endif
|
---|
217 |
|
---|
218 | ifelse(PUTORGET,`get',dnl
|
---|
219 | ! Check start and count arrays
|
---|
220 | ! ----------------------------
|
---|
221 |
|
---|
222 | j = 0
|
---|
223 | do i = `1, ndims'
|
---|
224 | if (cnts(i) > 1) then
|
---|
225 | j = j + 1
|
---|
226 | endif
|
---|
227 | enddo
|
---|
228 |
|
---|
229 | if (j > NUMDIMS) then
|
---|
230 | call message_set_routine('ncdf_getvar')
|
---|
231 | call message(msg_fatal, 'Data array does not have enough dimensions to hold requested data.')
|
---|
232 | endif
|
---|
233 | )dnl
|
---|
234 | ifelse(POINTER, `yes',dnl
|
---|
235 | ! Allocate Memory
|
---|
236 | ! ---------------
|
---|
237 |
|
---|
238 | allocate(rvalues(ALLOC_ARGS`'))
|
---|
239 | allocate(values(ALLOC_ARGS`'))
|
---|
240 | )dnl
|
---|
241 | ifelse(PUTORGET,`put',dnl
|
---|
242 | ! Convert units or copy data values
|
---|
243 | ! ---------------------------------
|
---|
244 |
|
---|
245 | if (present(units)) then
|
---|
246 | ncdf_units(:) = ' '
|
---|
247 | status = nf90_get_att(ncid_local, varid, 'units', ncdf_units)
|
---|
248 |
|
---|
249 | if (status /= nf90_noerr) then
|
---|
250 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute units not found for variable: "// name
|
---|
251 | call ncdf_error_handler(status)
|
---|
252 | endif
|
---|
253 |
|
---|
254 | call ut_convert(values, units, rvalues, ncdf_units)
|
---|
255 |
|
---|
256 | ! Only convert valid data
|
---|
257 | status = nf90_get_att(ncid_local, varid, 'valid_range', rrange)
|
---|
258 | if (status /= nf90_noerr) then
|
---|
259 | have_range=.false.
|
---|
260 | else
|
---|
261 | have_range=.true.
|
---|
262 | where(rvalues < rrange(1) .or. rvalues > rrange(2)) rvalues = values
|
---|
263 | endif
|
---|
264 |
|
---|
265 | else
|
---|
266 | rvalues = values
|
---|
267 | endif
|
---|
268 |
|
---|
269 | ! Scale variables if necessary
|
---|
270 | ! ----------------------------
|
---|
271 |
|
---|
272 | if (have_scale .and. have_offset) then
|
---|
273 | rvalues = (rvalues - add_offset) / scale_factor
|
---|
274 | else if (have_scale .and. (.not. have_offset)) then
|
---|
275 | rvalues = rvalues / scale_factor
|
---|
276 | else if ((.not. have_scale) .and. have_offset) then
|
---|
277 | rvalues = rvalues - add_offset
|
---|
278 | endif
|
---|
279 | )
|
---|
280 | ! Read/write values
|
---|
281 | ! -------------------
|
---|
282 |
|
---|
283 | status = `nf90_'PUTORGET`_var(ncid_local, varid, rvalues, start = strt, count = cnts)'
|
---|
284 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
---|
285 |
|
---|
286 | ifelse(PUTORGET,`get',dnl
|
---|
287 | ! Scale variables if necessary
|
---|
288 | ! ----------------------------
|
---|
289 |
|
---|
290 | if (have_scale .and. have_offset) then
|
---|
291 | rvalues = rvalues * scale_factor + add_offset
|
---|
292 | else if (have_scale .and. (.not. have_offset)) then
|
---|
293 | rvalues = rvalues * scale_factor
|
---|
294 | else if ((.not. have_scale) .and. have_offset) then
|
---|
295 | rvalues = rvalues + add_offset
|
---|
296 | endif
|
---|
297 |
|
---|
298 | ! Obtain range values - use defaults if not present
|
---|
299 | ! -------------------
|
---|
300 |
|
---|
301 | if (present(range))then
|
---|
302 | status = nf90_get_att(ncid_local, varid, 'valid_range', rrange)
|
---|
303 | if (status /= nf90_noerr) then
|
---|
304 | rrange = range
|
---|
305 | have_range = .false.
|
---|
306 | else
|
---|
307 | have_range = .true.
|
---|
308 | endif
|
---|
309 | endif
|
---|
310 |
|
---|
311 | ! Convert units - if present & not the same as target and data valid
|
---|
312 | ! -------------
|
---|
313 |
|
---|
314 | if (present(units)) then
|
---|
315 | ncdf_units(:) = ' '
|
---|
316 | status = nf90_get_att(ncid_local, varid, 'units', ncdf_units)
|
---|
317 | if (status /= nf90_noerr) then
|
---|
318 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute units not found for variable: "// name
|
---|
319 | call ncdf_error_handler(status)
|
---|
320 | endif
|
---|
321 |
|
---|
322 | call ut_convert(rvalues, ncdf_units, values, units)
|
---|
323 |
|
---|
324 | if(present(range) .and. have_range) then
|
---|
325 | call ut_convert(rrange, ncdf_units, range, units)
|
---|
326 |
|
---|
327 | ! Only convert valid data
|
---|
328 | where(rvalues < rrange(1) .or. rvalues > rrange(2)) values = rvalues
|
---|
329 | endif
|
---|
330 |
|
---|
331 | else
|
---|
332 | values = rvalues
|
---|
333 | if(present(range)) range = rrange
|
---|
334 | endif
|
---|
335 |
|
---|
336 | ! Add to counter of number of variables read
|
---|
337 | ! ------------------------------------------
|
---|
338 |
|
---|
339 | if (.not. (havegroup)) ncdf_read(varid) = .true.
|
---|
340 |
|
---|
341 | )dnl
|
---|
342 |
|
---|
343 | ! Dummy lines to avoid warnings - ignore
|
---|
344 | ! --------------------------------------
|
---|
345 | if (present(range)) then
|
---|
346 | continue
|
---|
347 | endif
|
---|
348 |
|
---|
349 | ! Deallocate memory
|
---|
350 | ! -----------------
|
---|
351 |
|
---|
352 | deallocate(rvalues)
|
---|
353 |
|
---|
354 | end subroutine NCDF_AFUN
|
---|