source: LMDZ4/trunk/libf/bibio/nf95_gw_var_m.F90 @ 2265

Last change on this file since 2265 was 1279, checked in by Laurent Fairhead, 15 years ago

Merged LMDZ4-dev branch changes r1241:1278 into the trunk
Running trunk and LMDZ4-dev in LMDZOR configuration on local
machine (sequential) and SX8 (4-proc) yields identical results
(restart and restartphy are identical binarily)
Log history from r1241 to r1278 is available by switching to
source:LMDZ4/branches/LMDZ4-dev-20091210

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