source: LMDZ5/branches/testing/tools/Max_diff_nc_with_lib/NetCDF95/nf95_gw_var.f90 @ 1887

Last change on this file since 1887 was 1795, checked in by Ehouarn Millour, 11 years ago

Version testing basee sur la r1794


Testing release based on r1794

File size: 11.0 KB
Line 
1module nf95_gw_var_m
2
3  use nf95_get_var_m, only: NF95_GET_VAR
4  use nf95_inquire_variable_m, only: nf95_inquire_variable
5  use nf95_inquire_dimension_m, only: nf95_inquire_dimension
6
7  implicit none
8
9  interface nf95_gw_var
10     ! "nf95_gw_var" stands for "NetCDF 1995 get whole variable".
11     ! These procedures read a whole NetCDF variable (coordinate or
12     ! primary) into an array.
13     ! The difference between the procedures is the rank and type of
14     ! argument "values".
15     ! The procedures do not check the type of the NetCDF variable.
16
17     module procedure nf95_gw_var_real_1d, nf95_gw_var_real_2d, &
18          nf95_gw_var_real_3d, nf95_gw_var_real_4d, nf95_gw_var_real_5d, &
19          nf95_gw_var_dble_1d, nf95_gw_var_dble_2d, nf95_gw_var_dble_3d, &
20          nf95_gw_var_dble_4d, nf95_gw_var_int_1d, nf95_gw_var_int_3d
21  end interface
22
23  private
24  public nf95_gw_var
25
26contains
27
28  subroutine nf95_gw_var_real_1d(ncid, varid, values)
29
30    ! Real type, the array has rank 1.
31
32    integer, intent(in):: ncid
33    integer, intent(in):: varid
34    real, pointer:: values(:)
35
36    ! Variables local to the procedure:
37    integer nclen
38    integer, pointer:: dimids(:)
39
40    !---------------------
41
42    call nf95_inquire_variable(ncid, varid, dimids=dimids)
43
44    if (size(dimids) /= 1) then
45       print *, "nf95_gw_var_real_1d:"
46       print *, "varid = ", varid
47       print *, "rank of NetCDF variable is ", size(dimids), ", not 1"
48       stop 1
49    end if
50
51    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen)
52    deallocate(dimids) ! pointer
53
54    allocate(values(nclen))
55    if (nclen /= 0) call NF95_GET_VAR(ncid, varid, values)
56
57  end subroutine nf95_gw_var_real_1d
58
59  !************************************
60
61  subroutine nf95_gw_var_real_2d(ncid, varid, values)
62
63    ! Real type, the array has rank 2.
64
65    integer, intent(in):: ncid
66    integer, intent(in):: varid
67    real, pointer:: values(:, :)
68
69    ! Variables local to the procedure:
70    integer nclen1, nclen2
71    integer, pointer:: dimids(:)
72
73    !---------------------
74
75    call nf95_inquire_variable(ncid, varid, dimids=dimids)
76
77    if (size(dimids) /= 2) then
78       print *, "nf95_gw_var_real_2d:"
79       print *, "varid = ", varid
80       print *, "rank of NetCDF variable is ", size(dimids), ", not 2"
81       stop 1
82    end if
83
84    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen1)
85    call nf95_inquire_dimension(ncid, dimids(2), nclen=nclen2)
86    deallocate(dimids) ! pointer
87
88    allocate(values(nclen1, nclen2))
89    if (nclen1 /= 0 .and. nclen2 /= 0) call NF95_GET_VAR(ncid, varid, values)
90
91  end subroutine nf95_gw_var_real_2d
92
93  !************************************
94
95  subroutine nf95_gw_var_real_3d(ncid, varid, values)
96
97    ! Real type, the array has rank 3.
98
99    integer, intent(in):: ncid
100    integer, intent(in):: varid
101    real, pointer:: values(:, :, :)
102
103    ! Variables local to the procedure:
104    integer nclen1, nclen2, nclen3
105    integer, pointer:: dimids(:)
106
107    !---------------------
108
109    call nf95_inquire_variable(ncid, varid, dimids=dimids)
110
111    if (size(dimids) /= 3) then
112       print *, "nf95_gw_var_real_3d:"
113       print *, "varid = ", varid
114       print *, "rank of NetCDF variable is ", size(dimids), ", not 3"
115       stop 1
116    end if
117
118    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen1)
119    call nf95_inquire_dimension(ncid, dimids(2), nclen=nclen2)
120    call nf95_inquire_dimension(ncid, dimids(3), nclen=nclen3)
121    deallocate(dimids) ! pointer
122
123    allocate(values(nclen1, nclen2, nclen3))
124    if (nclen1 * nclen2 * nclen3 /= 0) call NF95_GET_VAR(ncid, varid, values)
125
126  end subroutine nf95_gw_var_real_3d
127
128  !************************************
129
130  subroutine nf95_gw_var_real_4d(ncid, varid, values)
131
132    ! Real type, the array has rank 4.
133
134    integer, intent(in):: ncid
135    integer, intent(in):: varid
136    real, pointer:: values(:, :, :, :)
137
138    ! Variables local to the procedure:
139    integer len_dim(4), i
140    integer, pointer:: dimids(:)
141
142    !---------------------
143
144    call nf95_inquire_variable(ncid, varid, dimids=dimids)
145
146    if (size(dimids) /= 4) then
147       print *, "nf95_gw_var_real_4d:"
148       print *, "varid = ", varid
149       print *, "rank of NetCDF variable is ", size(dimids), ", not 4"
150       stop 1
151    end if
152
153    do i = 1, 4
154       call nf95_inquire_dimension(ncid, dimids(i), nclen=len_dim(i))
155    end do
156    deallocate(dimids) ! pointer
157
158    allocate(values(len_dim(1), len_dim(2), len_dim(3), len_dim(4)))
159    if (all(len_dim /= 0)) call NF95_GET_VAR(ncid, varid, values)
160
161  end subroutine nf95_gw_var_real_4d
162
163  !************************************
164
165  subroutine nf95_gw_var_real_5d(ncid, varid, values)
166
167    ! Real type, the array has rank 5.
168
169    integer, intent(in):: ncid
170    integer, intent(in):: varid
171    real, pointer:: values(:, :, :, :, :)
172
173    ! Variables local to the procedure:
174    integer len_dim(5), i
175    integer, pointer:: dimids(:)
176
177    !---------------------
178
179    call nf95_inquire_variable(ncid, varid, dimids=dimids)
180
181    if (size(dimids) /= 5) then
182       print *, "nf95_gw_var_real_5d:"
183       print *, "varid = ", varid
184       print *, "rank of NetCDF variable is ", size(dimids), ", not 5"
185       stop 1
186    end if
187
188    do i = 1, 5
189       call nf95_inquire_dimension(ncid, dimids(i), nclen=len_dim(i))
190    end do
191    deallocate(dimids) ! pointer
192
193    allocate(values(len_dim(1), len_dim(2), len_dim(3), len_dim(4), len_dim(5)))
194    if (all(len_dim /= 0)) call NF95_GET_VAR(ncid, varid, values)
195
196  end subroutine nf95_gw_var_real_5d
197
198  !************************************
199
200  subroutine nf95_gw_var_dble_1d(ncid, varid, values)
201
202    ! Double precision, the array has rank 1.
203
204    integer, intent(in):: ncid
205    integer, intent(in):: varid
206    double precision, pointer:: values(:)
207
208    ! Variables local to the procedure:
209    integer nclen
210    integer, pointer:: dimids(:)
211
212    !---------------------
213
214    call nf95_inquire_variable(ncid, varid, dimids=dimids)
215
216    if (size(dimids) /= 1) then
217       print *, "nf95_gw_var_dble_1d:"
218       print *, "varid = ", varid
219       print *, "rank of NetCDF variable is ", size(dimids), ", not 1"
220        stop 1
221    end if
222
223    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen)
224    deallocate(dimids) ! pointer
225
226    allocate(values(nclen))
227    if (nclen /= 0) call NF95_GET_VAR(ncid, varid, values)
228
229  end subroutine nf95_gw_var_dble_1d
230
231  !************************************
232
233  subroutine nf95_gw_var_dble_2d(ncid, varid, values)
234
235    ! Double precision, the array has rank 2.
236
237    integer, intent(in):: ncid
238    integer, intent(in):: varid
239    double precision, pointer:: values(:, :)
240
241    ! Variables local to the procedure:
242    integer nclen1, nclen2
243    integer, pointer:: dimids(:)
244
245    !---------------------
246
247    call nf95_inquire_variable(ncid, varid, dimids=dimids)
248
249    if (size(dimids) /= 2) then
250       print *, "nf95_gw_var_dble_2d:"
251       print *, "varid = ", varid
252       print *, "rank of NetCDF variable is ", size(dimids), ", not 2"
253       stop 1
254    end if
255
256    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen1)
257    call nf95_inquire_dimension(ncid, dimids(2), nclen=nclen2)
258    deallocate(dimids) ! pointer
259
260    allocate(values(nclen1, nclen2))
261    if (nclen1 /= 0 .and. nclen2 /= 0) call NF95_GET_VAR(ncid, varid, values)
262
263  end subroutine nf95_gw_var_dble_2d
264
265  !************************************
266
267  subroutine nf95_gw_var_dble_3d(ncid, varid, values)
268
269    ! Double precision, the array has rank 3.
270
271    integer, intent(in):: ncid
272    integer, intent(in):: varid
273    double precision, pointer:: values(:, :, :)
274
275    ! Variables local to the procedure:
276    integer nclen1, nclen2, nclen3
277    integer, pointer:: dimids(:)
278
279    !---------------------
280
281    call nf95_inquire_variable(ncid, varid, dimids=dimids)
282
283    if (size(dimids) /= 3) then
284       print *, "nf95_gw_var_dble_3d:"
285       print *, "varid = ", varid
286       print *, "rank of NetCDF variable is ", size(dimids), ", not 3"
287       stop 1
288    end if
289
290    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen1)
291    call nf95_inquire_dimension(ncid, dimids(2), nclen=nclen2)
292    call nf95_inquire_dimension(ncid, dimids(3), nclen=nclen3)
293    deallocate(dimids) ! pointer
294
295    allocate(values(nclen1, nclen2, nclen3))
296    if (nclen1 * nclen2 * nclen3 /= 0) call NF95_GET_VAR(ncid, varid, values)
297
298  end subroutine nf95_gw_var_dble_3d
299
300  !************************************
301
302  subroutine nf95_gw_var_dble_4d(ncid, varid, values)
303
304    ! Double precision, the array has rank 4.
305
306    integer, intent(in):: ncid
307    integer, intent(in):: varid
308    double precision, pointer:: values(:, :, :, :)
309
310    ! Variables local to the procedure:
311    integer len_dim(4), i
312    integer, pointer:: dimids(:)
313
314    !---------------------
315
316    call nf95_inquire_variable(ncid, varid, dimids=dimids)
317
318    if (size(dimids) /= 4) then
319       print *, "nf95_gw_var_dble_4d:"
320       print *, "varid = ", varid
321       print *, "rank of NetCDF variable is ", size(dimids), ", not 4"
322       stop 1
323    end if
324
325    do i = 1, 4
326       call nf95_inquire_dimension(ncid, dimids(i), nclen=len_dim(i))
327    end do
328    deallocate(dimids) ! pointer
329
330    allocate(values(len_dim(1), len_dim(2), len_dim(3), len_dim(4)))
331    if (all(len_dim /= 0)) call NF95_GET_VAR(ncid, varid, values)
332
333  end subroutine nf95_gw_var_dble_4d
334
335  !************************************
336
337  subroutine nf95_gw_var_int_1d(ncid, varid, values)
338
339    ! Integer type, the array has rank 1.
340
341    integer, intent(in):: ncid
342    integer, intent(in):: varid
343    integer, pointer:: values(:)
344
345    ! Variables local to the procedure:
346    integer nclen
347    integer, pointer:: dimids(:)
348
349    !---------------------
350
351    call nf95_inquire_variable(ncid, varid, dimids=dimids)
352
353    if (size(dimids) /= 1) then
354       print *, "nf95_gw_var_int_1d:"
355       print *, "varid = ", varid
356       print *, "rank of NetCDF variable is ", size(dimids), ", not 1"
357       stop 1
358    end if
359
360    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen)
361    deallocate(dimids) ! pointer
362
363    allocate(values(nclen))
364    if (nclen /= 0) call NF95_GET_VAR(ncid, varid, values)
365
366  end subroutine nf95_gw_var_int_1d
367
368  !************************************
369
370  subroutine nf95_gw_var_int_3d(ncid, varid, values)
371
372    ! Integer type, the array has rank 3.
373
374    integer, intent(in):: ncid
375    integer, intent(in):: varid
376    integer, pointer:: values(:, :, :)
377
378    ! Variables local to the procedure:
379    integer nclen1, nclen2, nclen3
380    integer, pointer:: dimids(:)
381
382    !---------------------
383
384    call nf95_inquire_variable(ncid, varid, dimids=dimids)
385
386    if (size(dimids) /= 3) then
387       print *, "nf95_gw_var_int_3d:"
388       print *, "varid = ", varid
389       print *, "rank of NetCDF variable is ", size(dimids), ", not 3"
390       stop 1
391    end if
392
393    call nf95_inquire_dimension(ncid, dimids(1), nclen=nclen1)
394    call nf95_inquire_dimension(ncid, dimids(2), nclen=nclen2)
395    call nf95_inquire_dimension(ncid, dimids(3), nclen=nclen3)
396    deallocate(dimids) ! pointer
397
398    allocate(values(nclen1, nclen2, nclen3))
399    if (nclen1 * nclen2 * nclen3 /= 0) call NF95_GET_VAR(ncid, varid, values)
400
401  end subroutine nf95_gw_var_int_3d
402
403end module nf95_gw_var_m
Note: See TracBrowser for help on using the repository browser.