1 | module nf95_get_var_m |
---|
2 | |
---|
3 | use netcdf, only: nf90_get_var, NF90_NOERR |
---|
4 | |
---|
5 | use nf95_abort_m, only: nf95_abort |
---|
6 | use check_start_count_m, only: check_start_count |
---|
7 | use nf95_get_missing_m, only: nf95_get_missing |
---|
8 | use nf95_constants, only: nf95_noerr |
---|
9 | |
---|
10 | implicit none |
---|
11 | |
---|
12 | interface nf95_get_var |
---|
13 | module procedure nf95_get_var_FourByteReal, nf95_get_var_FourByteInt, & |
---|
14 | nf95_get_var_text, nf95_get_var_1D_FourByteReal, & |
---|
15 | nf95_get_var_1D_FourByteInt, nf95_get_var_1D_EightByteReal, & |
---|
16 | nf95_get_var_2D_FourByteReal, nf95_get_var_2D_EightByteReal, & |
---|
17 | nf95_get_var_2D_twoByteInt, nf95_get_var_2D_FourByteInt, & |
---|
18 | nf95_get_var_3D_FourByteInt, nf95_get_var_3D_FourByteReal, & |
---|
19 | nf95_get_var_3D_EightByteReal, nf95_get_var_4D_FourByteReal, & |
---|
20 | nf95_get_var_4D_EightByteReal, nf95_get_var_5D_FourByteReal, & |
---|
21 | nf95_get_var_5D_EightByteReal |
---|
22 | end interface |
---|
23 | |
---|
24 | private |
---|
25 | public nf95_get_var |
---|
26 | |
---|
27 | integer ncerr_not_opt |
---|
28 | |
---|
29 | contains |
---|
30 | |
---|
31 | subroutine nf95_get_var_FourByteReal(ncid, varid, values, start, & |
---|
32 | new_missing, ncerr) |
---|
33 | |
---|
34 | use typesizes, only: FourByteReal |
---|
35 | |
---|
36 | integer, intent(in):: ncid, varid |
---|
37 | real(kind = FourByteReal), intent(out):: values |
---|
38 | integer, dimension(:), optional, intent(in):: start |
---|
39 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
40 | integer, intent(out), optional:: ncerr |
---|
41 | |
---|
42 | ! Local: |
---|
43 | character(len=*), parameter:: procedure_name = "nf95_get_var_FourByteReal" |
---|
44 | real(kind = FourByteReal) missing |
---|
45 | |
---|
46 | !------------------- |
---|
47 | |
---|
48 | include "nf95_get_var_scalar.h" |
---|
49 | |
---|
50 | end subroutine nf95_get_var_FourByteReal |
---|
51 | |
---|
52 | !*********************** |
---|
53 | |
---|
54 | subroutine nf95_get_var_FourByteInt(ncid, varid, values, start, new_missing, & |
---|
55 | ncerr) |
---|
56 | |
---|
57 | use typesizes, only: FourByteInt |
---|
58 | |
---|
59 | integer, intent(in):: ncid, varid |
---|
60 | integer(kind = FourByteInt), intent(out):: values |
---|
61 | integer, dimension(:), optional, intent(in):: start |
---|
62 | integer(kind = FourByteInt), optional, intent(in):: new_missing |
---|
63 | integer, intent(out), optional:: ncerr |
---|
64 | |
---|
65 | ! Local: |
---|
66 | character(len=*), parameter:: procedure_name = "nf95_get_var_FourByteInt" |
---|
67 | integer(kind = FourByteInt) missing |
---|
68 | |
---|
69 | !------------------- |
---|
70 | |
---|
71 | include "nf95_get_var_scalar.h" |
---|
72 | |
---|
73 | end subroutine nf95_get_var_FourByteInt |
---|
74 | |
---|
75 | !*********************** |
---|
76 | |
---|
77 | subroutine nf95_get_var_text(ncid, varid, values, start, new_missing, ncerr) |
---|
78 | |
---|
79 | integer, intent(in):: ncid, varid |
---|
80 | character, intent(out):: values |
---|
81 | integer, dimension(:), optional, intent(in):: start |
---|
82 | character, optional, intent(in):: new_missing |
---|
83 | integer, intent(out), optional:: ncerr |
---|
84 | |
---|
85 | ! Local: |
---|
86 | character(len=*), parameter:: procedure_name = "nf95_get_var_text" |
---|
87 | character missing |
---|
88 | |
---|
89 | !------------------- |
---|
90 | |
---|
91 | include "nf95_get_var_scalar.h" |
---|
92 | |
---|
93 | end subroutine nf95_get_var_text |
---|
94 | |
---|
95 | !*********************** |
---|
96 | |
---|
97 | subroutine nf95_get_var_1D_FourByteReal(ncid, varid, values, start, & |
---|
98 | count_nc, stride, map, new_missing, ncerr) |
---|
99 | |
---|
100 | use typesizes, only: FourByteReal |
---|
101 | |
---|
102 | integer, intent(in):: ncid, varid |
---|
103 | real(kind = FourByteReal), intent(out):: values(:) |
---|
104 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
105 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
106 | integer, intent(out), optional:: ncerr |
---|
107 | |
---|
108 | ! Local: |
---|
109 | character(len=*), parameter:: procedure_name = & |
---|
110 | "nf95_get_var_1D_FourByteReal" |
---|
111 | integer, parameter:: rank_values = 1 |
---|
112 | real(kind = FourByteReal) missing |
---|
113 | |
---|
114 | !------------------- |
---|
115 | |
---|
116 | include "nf95_get_var_array.h" |
---|
117 | |
---|
118 | end subroutine nf95_get_var_1D_FourByteReal |
---|
119 | |
---|
120 | !*********************** |
---|
121 | |
---|
122 | subroutine nf95_get_var_1D_FourByteInt(ncid, varid, values, start, & |
---|
123 | count_nc, stride, map, new_missing, ncerr) |
---|
124 | |
---|
125 | use typesizes, only: FourByteInt |
---|
126 | |
---|
127 | integer, intent(in):: ncid, varid |
---|
128 | integer(kind = FourByteInt), intent(out):: values(:) |
---|
129 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
130 | integer(kind = FourByteInt), optional, intent(in):: new_missing |
---|
131 | integer, intent(out), optional:: ncerr |
---|
132 | |
---|
133 | ! Local: |
---|
134 | character(len=*), parameter:: procedure_name = & |
---|
135 | "nf95_get_var_1D_FourByteInt" |
---|
136 | integer, parameter:: rank_values = 1 |
---|
137 | integer(kind = FourByteInt) missing |
---|
138 | |
---|
139 | !------------------- |
---|
140 | |
---|
141 | include "nf95_get_var_array.h" |
---|
142 | |
---|
143 | end subroutine nf95_get_var_1D_FourByteInt |
---|
144 | |
---|
145 | !*********************** |
---|
146 | |
---|
147 | subroutine nf95_get_var_1D_EightByteReal(ncid, varid, values, start, & |
---|
148 | count_nc, stride, map, new_missing, ncerr) |
---|
149 | |
---|
150 | use typesizes, only: eightByteReal |
---|
151 | |
---|
152 | integer, intent(in):: ncid, varid |
---|
153 | real(kind = EightByteReal), intent(out):: values(:) |
---|
154 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
155 | real(kind = EightByteReal), optional, intent(in):: new_missing |
---|
156 | integer, intent(out), optional:: ncerr |
---|
157 | |
---|
158 | ! Local: |
---|
159 | character(len=*), parameter:: procedure_name = & |
---|
160 | "nf95_get_var_1D_EightByteReal" |
---|
161 | integer, parameter:: rank_values = 1 |
---|
162 | real(kind = EightByteReal) missing |
---|
163 | |
---|
164 | !------------------- |
---|
165 | |
---|
166 | include "nf95_get_var_array.h" |
---|
167 | |
---|
168 | end subroutine nf95_get_var_1D_EightByteReal |
---|
169 | |
---|
170 | !*********************** |
---|
171 | |
---|
172 | subroutine nf95_get_var_2D_FourByteReal(ncid, varid, values, start, & |
---|
173 | count_nc, stride, map, new_missing, ncerr) |
---|
174 | |
---|
175 | use typesizes, only: FourByteReal |
---|
176 | |
---|
177 | integer, intent(in):: ncid, varid |
---|
178 | real(kind = FourByteReal), intent(out):: values(:, :) |
---|
179 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
180 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
181 | integer, intent(out), optional:: ncerr |
---|
182 | |
---|
183 | ! Local: |
---|
184 | character(len=*), parameter:: procedure_name = & |
---|
185 | "nf95_get_var_2D_FourByteReal" |
---|
186 | integer, parameter:: rank_values = 2 |
---|
187 | real(kind = FourByteReal) missing |
---|
188 | |
---|
189 | !------------------- |
---|
190 | |
---|
191 | include "nf95_get_var_array.h" |
---|
192 | |
---|
193 | end subroutine nf95_get_var_2D_FourByteReal |
---|
194 | |
---|
195 | !*********************** |
---|
196 | |
---|
197 | subroutine nf95_get_var_2D_EightByteReal(ncid, varid, values, start, & |
---|
198 | count_nc, stride, map, new_missing, ncerr) |
---|
199 | |
---|
200 | use typesizes, only: EightByteReal |
---|
201 | |
---|
202 | integer, intent(in):: ncid, varid |
---|
203 | real(kind = EightByteReal), intent(out):: values(:, :) |
---|
204 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
205 | real(kind = EightByteReal), optional, intent(in):: new_missing |
---|
206 | integer, intent(out), optional:: ncerr |
---|
207 | |
---|
208 | ! Local: |
---|
209 | character(len=*), parameter:: procedure_name = & |
---|
210 | "nf95_get_var_2D_EightByteReal" |
---|
211 | integer, parameter:: rank_values = 2 |
---|
212 | real(kind = EightByteReal) missing |
---|
213 | |
---|
214 | !------------------- |
---|
215 | |
---|
216 | include "nf95_get_var_array.h" |
---|
217 | |
---|
218 | end subroutine nf95_get_var_2D_EightByteReal |
---|
219 | |
---|
220 | !*********************** |
---|
221 | |
---|
222 | subroutine nf95_get_var_2D_TwoByteInt(ncid, varid, values, start, count_nc, & |
---|
223 | stride, map, new_missing, ncerr) |
---|
224 | |
---|
225 | use typesizes, only: TwoByteInt |
---|
226 | |
---|
227 | integer, intent(in):: ncid, varid |
---|
228 | integer(kind = TwoByteInt), intent(out):: values(:, :) |
---|
229 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
230 | integer(kind = TwoByteInt), optional, intent(in):: new_missing |
---|
231 | integer, intent(out), optional:: ncerr |
---|
232 | |
---|
233 | ! Local: |
---|
234 | character(len=*), parameter:: procedure_name = "nf95_get_var_2D_TwoByteInt" |
---|
235 | integer, parameter:: rank_values = 2 |
---|
236 | integer(kind = TwoByteInt) missing |
---|
237 | |
---|
238 | !------------------- |
---|
239 | |
---|
240 | include "nf95_get_var_array.h" |
---|
241 | |
---|
242 | end subroutine nf95_get_var_2D_TwoByteInt |
---|
243 | |
---|
244 | !*********************** |
---|
245 | |
---|
246 | subroutine nf95_get_var_2D_FourByteInt(ncid, varid, values, start, & |
---|
247 | count_nc, stride, map, new_missing, ncerr) |
---|
248 | |
---|
249 | use typesizes, only: FourByteInt |
---|
250 | |
---|
251 | integer, intent(in):: ncid, varid |
---|
252 | integer(kind = FourByteInt), intent(out):: values(:, :) |
---|
253 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
254 | integer(kind = FourByteInt), optional, intent(in):: new_missing |
---|
255 | integer, intent(out), optional:: ncerr |
---|
256 | |
---|
257 | ! Local: |
---|
258 | character(len=*), parameter:: procedure_name = & |
---|
259 | "nf95_get_var_2D_FourByteInt" |
---|
260 | integer, parameter:: rank_values = 2 |
---|
261 | integer(kind = FourByteInt) missing |
---|
262 | |
---|
263 | !------------------- |
---|
264 | |
---|
265 | include "nf95_get_var_array.h" |
---|
266 | |
---|
267 | end subroutine nf95_get_var_2D_FourByteInt |
---|
268 | |
---|
269 | !*********************** |
---|
270 | |
---|
271 | subroutine nf95_get_var_3D_FourByteInt(ncid, varid, values, start, & |
---|
272 | count_nc, stride, map, new_missing, ncerr) |
---|
273 | |
---|
274 | use typesizes, only: FourByteInt |
---|
275 | |
---|
276 | integer, intent(in):: ncid, varid |
---|
277 | integer(kind = FourByteInt), intent(out):: values(:, :, :) |
---|
278 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
279 | integer(kind = FourByteInt), optional, intent(in):: new_missing |
---|
280 | integer, intent(out), optional:: ncerr |
---|
281 | |
---|
282 | ! Local: |
---|
283 | character(len=*), parameter:: procedure_name = & |
---|
284 | "nf95_get_var_3D_FourByteInt" |
---|
285 | integer, parameter:: rank_values = 3 |
---|
286 | integer(kind = FourByteInt) missing |
---|
287 | |
---|
288 | !------------------- |
---|
289 | |
---|
290 | include "nf95_get_var_array.h" |
---|
291 | |
---|
292 | end subroutine nf95_get_var_3D_FourByteInt |
---|
293 | |
---|
294 | !*********************** |
---|
295 | |
---|
296 | subroutine nf95_get_var_3D_FourByteReal(ncid, varid, values, start, & |
---|
297 | count_nc, stride, map, new_missing, ncerr) |
---|
298 | |
---|
299 | use typesizes, only: FourByteReal |
---|
300 | |
---|
301 | integer, intent(in):: ncid, varid |
---|
302 | real(kind = FourByteReal), intent(out):: values(:, :, :) |
---|
303 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
304 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
305 | integer, intent(out), optional:: ncerr |
---|
306 | |
---|
307 | ! Local: |
---|
308 | character(len=*), parameter:: procedure_name = & |
---|
309 | "nf95_get_var_3D_FourByteReal" |
---|
310 | integer, parameter:: rank_values = 3 |
---|
311 | real(kind = FourByteReal) missing |
---|
312 | |
---|
313 | !------------------- |
---|
314 | |
---|
315 | include "nf95_get_var_array.h" |
---|
316 | |
---|
317 | end subroutine nf95_get_var_3D_FourByteReal |
---|
318 | |
---|
319 | !*********************** |
---|
320 | |
---|
321 | subroutine nf95_get_var_3D_EightByteReal(ncid, varid, values, start, & |
---|
322 | count_nc, stride, map, new_missing, ncerr) |
---|
323 | |
---|
324 | use typesizes, only: eightByteReal |
---|
325 | |
---|
326 | integer, intent(in):: ncid, varid |
---|
327 | real(kind = EightByteReal), intent(out):: values(:, :, :) |
---|
328 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
329 | real(kind = EightByteReal), optional, intent(in):: new_missing |
---|
330 | integer, intent(out), optional:: ncerr |
---|
331 | |
---|
332 | ! Local: |
---|
333 | character(len=*), parameter:: procedure_name = & |
---|
334 | "nf95_get_var_3D_EightByteReal" |
---|
335 | integer, parameter:: rank_values = 3 |
---|
336 | real(kind = EightByteReal) missing |
---|
337 | |
---|
338 | !------------------- |
---|
339 | |
---|
340 | include "nf95_get_var_array.h" |
---|
341 | |
---|
342 | end subroutine nf95_get_var_3D_EightByteReal |
---|
343 | |
---|
344 | !*********************** |
---|
345 | |
---|
346 | subroutine nf95_get_var_4D_FourByteReal(ncid, varid, values, start, & |
---|
347 | count_nc, stride, map, new_missing, ncerr) |
---|
348 | |
---|
349 | use typesizes, only: FourByteReal |
---|
350 | |
---|
351 | integer, intent(in):: ncid, varid |
---|
352 | real(kind = FourByteReal), intent(out):: values(:, :, :, :) |
---|
353 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
354 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
355 | integer, intent(out), optional:: ncerr |
---|
356 | |
---|
357 | ! Local: |
---|
358 | character(len=*), parameter:: procedure_name = & |
---|
359 | "nf95_get_var_4D_FourByteReal" |
---|
360 | integer, parameter:: rank_values = 4 |
---|
361 | real(kind = FourByteReal) missing |
---|
362 | |
---|
363 | !------------------- |
---|
364 | |
---|
365 | include "nf95_get_var_array.h" |
---|
366 | |
---|
367 | end subroutine nf95_get_var_4D_FourByteReal |
---|
368 | |
---|
369 | !*********************** |
---|
370 | |
---|
371 | subroutine nf95_get_var_4D_EightByteReal(ncid, varid, values, start, & |
---|
372 | count_nc, stride, map, new_missing, ncerr) |
---|
373 | |
---|
374 | use typesizes, only: EightByteReal |
---|
375 | |
---|
376 | integer, intent(in):: ncid, varid |
---|
377 | real(kind = EightByteReal), intent(out):: values(:, :, :, :) |
---|
378 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
379 | real(kind = EightByteReal), optional, intent(in):: new_missing |
---|
380 | integer, intent(out), optional:: ncerr |
---|
381 | |
---|
382 | ! Local: |
---|
383 | character(len=*), parameter:: procedure_name = & |
---|
384 | "nf95_get_var_4D_EightByteReal" |
---|
385 | integer, parameter:: rank_values = 4 |
---|
386 | real(kind = EightByteReal) missing |
---|
387 | |
---|
388 | !------------------- |
---|
389 | |
---|
390 | include "nf95_get_var_array.h" |
---|
391 | |
---|
392 | end subroutine nf95_get_var_4D_EightByteReal |
---|
393 | |
---|
394 | !*********************** |
---|
395 | |
---|
396 | subroutine nf95_get_var_5D_FourByteReal(ncid, varid, values, start, & |
---|
397 | count_nc, stride, map, new_missing, ncerr) |
---|
398 | |
---|
399 | use typesizes, only: FourByteReal |
---|
400 | |
---|
401 | integer, intent(in):: ncid, varid |
---|
402 | real(kind = FourByteReal), intent(out):: values(:, :, :, :, :) |
---|
403 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
404 | real(kind = FourByteReal), optional, intent(in):: new_missing |
---|
405 | integer, intent(out), optional:: ncerr |
---|
406 | |
---|
407 | ! Local: |
---|
408 | character(len=*), parameter:: procedure_name = & |
---|
409 | "nf95_get_var_5D_FourByteReal" |
---|
410 | integer, parameter:: rank_values = 5 |
---|
411 | real(kind = FourByteReal) missing |
---|
412 | |
---|
413 | !------------------- |
---|
414 | |
---|
415 | include "nf95_get_var_array.h" |
---|
416 | |
---|
417 | end subroutine nf95_get_var_5D_FourByteReal |
---|
418 | |
---|
419 | !*********************** |
---|
420 | |
---|
421 | subroutine nf95_get_var_5D_EightByteReal(ncid, varid, values, start, & |
---|
422 | count_nc, stride, map, new_missing, ncerr) |
---|
423 | |
---|
424 | use typesizes, only: EightByteReal |
---|
425 | |
---|
426 | integer, intent(in):: ncid, varid |
---|
427 | real(kind = EightByteReal), intent(out):: values(:, :, :, :, :) |
---|
428 | integer, dimension(:), optional, intent(in):: start, count_nc, stride, map |
---|
429 | real(kind = EightByteReal), optional, intent(in):: new_missing |
---|
430 | integer, intent(out), optional:: ncerr |
---|
431 | |
---|
432 | ! Local: |
---|
433 | character(len=*), parameter:: procedure_name = & |
---|
434 | "nf95_get_var_5D_EightByteReal" |
---|
435 | integer, parameter:: rank_values = 5 |
---|
436 | real(kind = EightByteReal) missing |
---|
437 | |
---|
438 | !------------------- |
---|
439 | |
---|
440 | include "nf95_get_var_array.h" |
---|
441 | |
---|
442 | end subroutine nf95_get_var_5D_EightByteReal |
---|
443 | |
---|
444 | end module nf95_get_var_m |
---|