source: LMDZ6/branches/Amaury_dev/tools/netcdf95/Variables/nf95_gw_var.f90 @ 5080

Last change on this file since 5080 was 4918, checked in by Laurent Fairhead, 7 months ago

Reintegrated NetCDF95 in LMDZ so that it is compiled and made available by the makelmdz_fcm script.
The makelmdz_fcm creates the libnetcdf95 library and copies it in the tools/netcdf/lib directory, copying
the mod files in the tools/netcdf/include library.

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