FMS  2026.01
Flexible Modeling System
netcdf_read_data.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @file
19 !> @brief Routines for retrieving variable values from netcdf files with different
20 !! dimension sizes for the @ref netcdf_read_data interface
21 
22 !> @addtogroup netcdf_io_mod
23 !> @{
24 
25 !> @brief Character read 0d function.
26 subroutine char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
27  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
28  character(len=*), intent(in) :: variable_name !< Variable name.
29  character(len=*), intent(inout) :: buf !< Array that the data will be read into
30  integer, intent(in), optional :: corner !< Index of the string to read if the variable
31  !! contains a 1D array of strings.
32  character(len=200), intent(in):: append_error_msg !< Msg to be appended to FATAL error message
33  integer, intent(inout) :: err
34  integer, intent(in) :: varid
35 
36  integer, dimension(2) :: start
37  integer :: ndims
38  character, dimension(:), allocatable :: charbuf
39  integer, dimension(:), allocatable :: dimsizes
40  integer :: i
41 
42  start(:) = 1
43  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
44  allocate(dimsizes(ndims))
45  call get_variable_size(fileobj, variable_name, dimsizes, broadcast=.false.)
46  allocate(charbuf(dimsizes(1)))
47  charbuf(:) = ""
48  if (ndims .eq. 2) then
49  if (present(corner)) then
50  start(2) = corner
51  endif
52  dimsizes(2) = 1
53  elseif (ndims .gt. 2) then
54  call error("Only scalar and 1d string values are currently supported: "//trim(append_error_msg))
55  endif
56  err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
57  if (len(buf) .lt. dimsizes(1)) then
58  call error("character buffer is too small; increase length: "//trim(append_error_msg))
59  endif
60  buf = ""
61  do i = 1, dimsizes(1)
62  if (charbuf(i) .eq. char(0)) then
63  exit
64  endif
65  buf(i:i) = charbuf(i)
66  enddo
67  deallocate(charbuf)
68  deallocate(dimsizes)
69 end subroutine char_read_0d
70 
71 !> @brief Character read 1d function.
72 subroutine char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
73 
74  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
75  character(len=*), intent(in) :: variable_name !< Variable name.
76  character(len=*), dimension(:), intent(inout) :: buf !< Array that the data
77  !! will be read into.
78  integer, dimension(2), intent(in) :: c
79  character(len=200), intent(in) :: append_error_msg !< Msg to be appended to FATAL error message
80  integer, intent(inout) :: err
81  integer, intent(in) :: varid
82 
83  integer :: ndims
84  integer, dimension(2) :: start
85  integer, dimension(2) :: dimsizes
86  character, dimension(:,:), allocatable :: charbuf
87  character(len=1024) :: sbuf
88  integer :: i
89  integer :: j
90 
91  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
92  if (ndims .ne. 2) then
93  call error(trim(variable_name)//"must be 2d dimensional in netcdf file.")
94  endif
95  start(1) = 1
96  start(2) = c(1)
97  call get_variable_size(fileobj, variable_name, dimsizes, .false.)
98  dimsizes(2) = dimsizes(2) - start(2) + 1
99  call allocate_array(charbuf, dimsizes, initialize=.true.)
100  err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
101  if (len(buf(1)) .lt. dimsizes(1)) then
102  call error("character buffer is too small; increase length: "//trim(append_error_msg))
103  endif
104  if (size(buf) .lt. dimsizes(2)) then
105  call error("incorrect buffer array size:: "//trim(append_error_msg))
106  endif
107  do i = start(2), start(2)+dimsizes(2)-1
108  sbuf = ""
109  do j = 1, dimsizes(1)
110  if (charbuf(j, i-start(2)+1) .eq. char(0)) then
111  exit
112  endif
113  sbuf(j:j) = charbuf(j, i-start(2)+1)
114  enddo
115  call string_copy(buf(i-start(2)+1), sbuf)
116  enddo
117  deallocate(charbuf)
118 
119 end subroutine char_read_1d
120 
121 !> @brief Read in data from a variable in a netcdf file.
122 subroutine netcdf_read_data_0d(fileobj, variable_name, buf, unlim_dim_level, &
123  corner, broadcast)
124 
125  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
126  character(len=*), intent(in) :: variable_name !< Variable name.
127  class(*), intent(inout) :: buf !< Array that the data will be read into.
128  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
129  integer, intent(in), optional :: corner !< Index of the string to read if the variable
130  !! contains a 1D array of strings.
131  logical, intent(in), optional :: broadcast !< Flag controlling whether or
132  !! not the data will be
133  !! broadcasted to non
134  !! "I/O root" ranks.
135  !! The broadcast will be done
136  !! by default.
137 
138  logical :: bcast
139  integer :: err
140  integer :: varid
141  integer :: unlim_dim_index
142  integer, dimension(1) :: c
143  character(len=1024), dimension(1) :: buf1d
144  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
145 
146  append_error_msg = "netcdf_read_data_0d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
147 
148  if (present(broadcast)) then
149  bcast = broadcast
150  else
151  bcast = .true.
152  endif
153  if (fileobj%use_netcdf_mpi) bcast = .false.
154 
155  c(1) = 1
156  if (present(unlim_dim_level)) then
157  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
158  broadcast=bcast)
159  if (unlim_dim_index .ne. 1) then
160  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
161  endif
162  c(unlim_dim_index) = unlim_dim_level
163  endif
164 
165  if (fileobj%is_root) then
166  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
167 
168  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
169  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
170  call check_netcdf_code(err, append_error_msg)
171  endif
172 
173  select type(buf)
174  type is (integer(kind=i4_kind))
175  err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
176  type is (integer(kind=i8_kind))
177  err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
178  type is (real(kind=r4_kind))
179  err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
180  type is (real(kind=r8_kind))
181  err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
182  type is (character(len=*))
183  call char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
184  class default
185  call error("Unsupported variable type: "//trim(append_error_msg))
186  end select
187  call check_netcdf_code(err, append_error_msg)
188  call unpack_data_0d(fileobj, varid, variable_name, buf)
189  endif
190 
191  if (bcast) then
192  select type(buf)
193  type is (integer(kind=i4_kind))
194  call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
195  type is (integer(kind=i8_kind))
196  call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
197  type is (real(kind=r4_kind))
198  call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
199  type is (real(kind=r8_kind))
200  call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
201  type is (character(len=*))
202  call string_copy(buf1d(1), buf)
203  call mpp_broadcast(buf1d, len(buf1d(1)), fileobj%io_root, pelist=fileobj%pelist)
204  call string_copy(buf, buf1d(1))
205  class default
206  call error("Unsupported variable type: "//trim(append_error_msg))
207  end select
208  endif
209 end subroutine netcdf_read_data_0d
210 
211 
212 !> @brief Read in data from a variable in a netcdf file.
213 subroutine netcdf_read_data_1d(fileobj, variable_name, buf, unlim_dim_level, &
214  corner, edge_lengths, broadcast)
215 
216  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
217  character(len=*), intent(in) :: variable_name !< Variable name.
218  class(*), dimension(:), intent(inout) :: buf !< Array that the data
219  !! will be read into.
220  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
221  !! level.
222  integer, dimension(:), intent(in), optional :: corner !< Array of starting
223  !! indices describing
224  !! where the data
225  !! will be read from.
226  integer, dimension(:), intent(in), optional :: edge_lengths !< The number of
227  !! elements that
228  !! will be read
229  !! in each dimension.
230  logical, intent(in), optional :: broadcast !< Flag controlling whether or
231  !! not the data will be
232  !! broadcasted to non
233  !! "I/O root" ranks.
234  !! The broadcast will be done
235  !! by default.
236 
237  logical :: bcast
238  integer :: err
239  integer :: varid
240  integer :: unlim_dim_index
241  integer, dimension(2) :: c
242  integer, dimension(2) :: e
243  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
244 
245  append_error_msg = "netcdf_read_data_1d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
246 
247  if (present(broadcast)) then
248  bcast = broadcast
249  else
250  bcast = .true.
251  endif
252  if (fileobj%use_netcdf_mpi) bcast = .false.
253 
254  c(:) = 1
255  if (present(corner)) then
256  c(1) = corner(1)
257  endif
258  e(:) = 1
259  if (present(edge_lengths)) then
260  e(1) = edge_lengths(1)
261  else
262  e(1) = size(buf)
263  endif
264  if (present(unlim_dim_level)) then
265  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
266  broadcast=bcast)
267  if (unlim_dim_index .ne. 2) then
268  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
269  endif
270  c(unlim_dim_index) = unlim_dim_level
271  endif
272 
273  if (fileobj%is_root) then
274  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
275 
276  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
277  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
278  call check_netcdf_code(err, append_error_msg)
279  endif
280 
281  select type(buf)
282  type is (integer(kind=i4_kind))
283  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
284  type is (integer(kind=i8_kind))
285  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
286  type is (real(kind=r4_kind))
287  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
288  type is (real(kind=r8_kind))
289  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
290  type is (character(len=*))
291  call char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
292  class default
293  call error("Unsupported variable type: "//trim(append_error_msg))
294  end select
295  call check_netcdf_code(err, append_error_msg)
296  call unpack_data_1d(fileobj, varid, variable_name, buf)
297  endif
298 
299  if (bcast) then
300  select type(buf)
301  type is (integer(kind=i4_kind))
302  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
303  type is (integer(kind=i8_kind))
304  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
305  type is (real(kind=r4_kind))
306  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
307  type is (real(kind=r8_kind))
308  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
309  type is (character(len=*))
310  call mpp_broadcast(buf, len(buf(1)), fileobj%io_root, pelist=fileobj%pelist)
311  class default
312  call error("Unsupported variable type: "//trim(append_error_msg))
313  end select
314  endif
315 end subroutine netcdf_read_data_1d
316 
317 
318 !> @brief Read in data from a variable in a netcdf file.
319 subroutine netcdf_read_data_2d(fileobj, variable_name, buf, unlim_dim_level, &
320  corner, edge_lengths, broadcast)
321 
322  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
323  character(len=*), intent(in) :: variable_name !< Variable name.
324  class(*), dimension(:,:), intent(inout) :: buf !< Array that the data
325  !! will be read into.
326  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
327  !! level.
328  integer, dimension(:), intent(in), optional :: corner !< Array of starting
329  !! indices describing
330  !! where the data
331  !! will be read from.
332  integer, dimension(:), intent(in), optional :: edge_lengths !< The number of
333  !! elements that
334  !! will be read
335  !! in each dimension.
336  logical, intent(in), optional :: broadcast !< Flag controlling whether or
337  !! not the data will be
338  !! broadcasted to non
339  !! "I/O root" ranks.
340  !! The broadcast will be done
341  !! by default.
342 
343  logical :: bcast
344  integer :: err
345  integer :: varid
346  integer :: unlim_dim_index
347  integer, dimension(3) :: c
348  integer, dimension(3) :: e
349  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
350 
351  append_error_msg = "netcdf_read_data_2d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
352 
353  if (present(broadcast)) then
354  bcast = broadcast
355  else
356  bcast = .true.
357  endif
358  if (fileobj%use_netcdf_mpi) bcast = .false.
359 
360  c(:) = 1
361  if (present(corner)) then
362  c(1:2) = corner(1:2)
363  endif
364  e(:) = 1
365  if (present(edge_lengths)) then
366  e(1:2) = edge_lengths(1:2)
367  else
368  e(1:2) = shape(buf)
369  endif
370  if (present(unlim_dim_level)) then
371  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
372  broadcast=bcast)
373  if (unlim_dim_index .ne. 3) then
374  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
375  endif
376  c(unlim_dim_index) = unlim_dim_level
377  endif
378 
379  if (fileobj%is_root) then
380  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
381 
382  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
383  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
384  call check_netcdf_code(err, append_error_msg)
385  endif
386 
387  select type(buf)
388  type is (integer(kind=i4_kind))
389  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
390  type is (integer(kind=i8_kind))
391  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
392  type is (real(kind=r4_kind))
393  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
394  type is (real(kind=r8_kind))
395  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
396  class default
397  call error("Unsupported variable type: "//trim(append_error_msg))
398  end select
399  call check_netcdf_code(err, append_error_msg)
400  call unpack_data_2d(fileobj, varid, variable_name, buf)
401  endif
402 
403  if (bcast) then
404  select type(buf)
405  type is (integer(kind=i4_kind))
406  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
407  type is (integer(kind=i8_kind))
408  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
409  type is (real(kind=r4_kind))
410  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
411  type is (real(kind=r8_kind))
412  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
413  class default
414  call error("Unsupported variable type: "//trim(append_error_msg))
415  end select
416  endif
417 end subroutine netcdf_read_data_2d
418 
419 
420 !> @brief Read in data from a variable in a netcdf file.
421 subroutine netcdf_read_data_3d(fileobj, variable_name, buf, unlim_dim_level, &
422  corner, edge_lengths, broadcast)
423 
424  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
425  character(len=*), intent(in) :: variable_name !< Variable name.
426  class(*), dimension(:,:,:), intent(inout) :: buf !< Array that the data
427  !! will be read into.
428  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
429  !! level.
430  integer, dimension(:), intent(in), optional :: corner !< Array of starting
431  !! indices describing
432  !! where the data
433  !! will be read from.
434  integer, dimension(:), intent(in), optional :: edge_lengths !< The number of
435  !! elements that
436  !! will be read
437  !! in each dimension.
438  logical, intent(in), optional :: broadcast !< Flag controlling whether or
439  !! not the data will be
440  !! broadcasted to non
441  !! "I/O root" ranks.
442  !! The broadcast will be done
443  !! by default.
444 
445  logical :: bcast
446  integer :: err
447  integer :: varid
448  integer :: unlim_dim_index
449  integer, dimension(4) :: c
450  integer, dimension(4) :: e
451  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
452 
453  append_error_msg = "netcdf_read_data_3d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
454 
455  if (present(broadcast)) then
456  bcast = broadcast
457  else
458  bcast = .true.
459  endif
460  if (fileobj%use_netcdf_mpi) bcast = .false.
461 
462  c(:) = 1
463  if (present(corner)) then
464  c(1:3) = corner(1:3)
465  endif
466  e(:) = 1
467  if (present(edge_lengths)) then
468  e(1:3) = edge_lengths(1:3)
469  else
470  e(1:3) = shape(buf)
471  endif
472  if (present(unlim_dim_level)) then
473  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
474  broadcast=bcast)
475  if (unlim_dim_index .ne. 4) then
476  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
477  endif
478  c(unlim_dim_index) = unlim_dim_level
479  endif
480 
481  if (fileobj%is_root) then
482  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
483 
484  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
485  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
486  call check_netcdf_code(err, append_error_msg)
487  endif
488 
489  select type(buf)
490  type is (integer(kind=i4_kind))
491  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
492  type is (integer(kind=i8_kind))
493  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
494  type is (real(kind=r4_kind))
495  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
496  type is (real(kind=r8_kind))
497  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
498  class default
499  call error("Unsupported variable type: "//trim(append_error_msg))
500  end select
501  call check_netcdf_code(err, append_error_msg)
502  call unpack_data_3d(fileobj, varid, variable_name, buf)
503  endif
504 
505  if (bcast) then
506  select type(buf)
507  type is (integer(kind=i4_kind))
508  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
509  type is (integer(kind=i8_kind))
510  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
511  type is (real(kind=r4_kind))
512  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
513  type is (real(kind=r8_kind))
514  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
515  class default
516  call error("Unsupported variable type: "//trim(append_error_msg))
517  end select
518  endif
519 end subroutine netcdf_read_data_3d
520 
521 
522 !> @brief Read in data from a variable in a netcdf file.
523 subroutine netcdf_read_data_4d(fileobj, variable_name, buf, unlim_dim_level, &
524  corner, edge_lengths, broadcast)
525 
526  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
527  character(len=*), intent(in) :: variable_name !< Variable name.
528  class(*), dimension(:,:,:,:), intent(inout) :: buf !< Array that the data
529  !! will be read into.
530  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
531  !! level.
532  integer, dimension(:), intent(in), optional :: corner !< Array of starting
533  !! indices describing
534  !! where the data
535  !! will be read from.
536  integer, dimension(:), intent(in), optional :: edge_lengths !< The number of
537  !! elements that
538  !! will be read
539  !! in each dimension.
540  logical, intent(in), optional :: broadcast !< Flag controlling whether or
541  !! not the data will be
542  !! broadcasted to non
543  !! "I/O root" ranks.
544  !! The broadcast will be done
545  !! by default.
546 
547  logical :: bcast
548  integer :: err
549  integer :: varid
550  integer :: unlim_dim_index
551  integer, dimension(5) :: c
552  integer, dimension(5) :: e
553  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
554 
555  append_error_msg = "netcdf_read_data_4d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
556 
557  if (present(broadcast)) then
558  bcast = broadcast
559  else
560  bcast = .true.
561  endif
562  if (fileobj%use_netcdf_mpi) bcast = .false.
563 
564  c(:) = 1
565  if (present(corner)) then
566  c(1:4) = corner(1:4)
567  endif
568  e(:) = 1
569  if (present(edge_lengths)) then
570  e(1:4) = edge_lengths(1:4)
571  else
572  e(1:4) = shape(buf)
573  endif
574  if (present(unlim_dim_level)) then
575  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
576  broadcast=bcast)
577  if (unlim_dim_index .ne. 5) then
578  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
579  endif
580  c(unlim_dim_index) = unlim_dim_level
581  endif
582 
583  if (fileobj%is_root) then
584  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
585 
586  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
587  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
588  call check_netcdf_code(err, append_error_msg)
589  endif
590 
591  select type(buf)
592  type is (integer(kind=i4_kind))
593  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
594  type is (integer(kind=i8_kind))
595  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
596  type is (real(kind=r4_kind))
597  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
598  type is (real(kind=r8_kind))
599  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
600  class default
601  call error("Unsupported variable type: "//trim(append_error_msg))
602  end select
603  call check_netcdf_code(err, append_error_msg)
604  call unpack_data_4d(fileobj, varid, variable_name, buf)
605  endif
606 
607  if (bcast) then
608  select type(buf)
609  type is (integer(kind=i4_kind))
610  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
611  type is (integer(kind=i8_kind))
612  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
613  type is (real(kind=r4_kind))
614  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
615  type is (real(kind=r8_kind))
616  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
617  class default
618  call error("Unsupported variable type: "//trim(append_error_msg))
619  end select
620  endif
621 end subroutine netcdf_read_data_4d
622 
623 
624 !> @brief Read in data from a variable in a netcdf file.
625 subroutine netcdf_read_data_5d(fileobj, variable_name, buf, unlim_dim_level, &
626  corner, edge_lengths, broadcast)
627 
628  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
629  character(len=*), intent(in) :: variable_name !< Variable name.
630  class(*), dimension(:,:,:,:,:), intent(inout) :: buf !< Array that the data
631  !! will be read into.
632  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
633  !! level.
634  integer, dimension(:), intent(in), optional :: corner !< Array of starting
635  !! indices describing
636  !! where the data
637  !! will be read from.
638  integer, dimension(:), intent(in), optional :: edge_lengths !< The number of
639  !! elements that
640  !! will be read
641  !! in each dimension.
642  logical, intent(in), optional :: broadcast !< Flag controlling whether or
643  !! not the data will be
644  !! broadcasted to non
645  !! "I/O root" ranks.
646  !! The broadcast will be done
647  !! by default.
648 
649  logical :: bcast
650  integer :: err
651  integer :: varid
652  integer :: unlim_dim_index
653  integer, dimension(6) :: c
654  integer, dimension(6) :: e
655  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
656 
657  append_error_msg = "netcdf_read_data_5d: file:"//trim(fileobj%path)//"- variable:"//trim(variable_name)
658 
659  if (present(broadcast)) then
660  bcast = broadcast
661  else
662  bcast = .true.
663  endif
664  if (fileobj%use_netcdf_mpi) bcast = .false.
665 
666  c(:) = 1
667  if (present(corner)) then
668  c(1:5) = corner(1:5)
669  endif
670  e(:) = 1
671  if (present(edge_lengths)) then
672  e(1:5) = edge_lengths(1:5)
673  else
674  e(1:5) = shape(buf)
675  endif
676  if (present(unlim_dim_level)) then
677  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
678  broadcast=bcast)
679  if (unlim_dim_index .ne. 6) then
680  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
681  endif
682  c(unlim_dim_index) = unlim_dim_level
683  endif
684 
685  if (fileobj%is_root) then
686  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
687 
688  if (fileobj%use_netcdf_mpi.and.fileobj%use_collective) then
689  err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
690  call check_netcdf_code(err, append_error_msg)
691  endif
692 
693  select type(buf)
694  type is (integer(kind=i4_kind))
695  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
696  type is (integer(kind=i8_kind))
697  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
698  type is (real(kind=r4_kind))
699  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
700  type is (real(kind=r8_kind))
701  err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
702  class default
703  call error("Unsupported variable type: "//trim(append_error_msg))
704  end select
705  call check_netcdf_code(err, append_error_msg)
706  call unpack_data_5d(fileobj, varid, variable_name, buf)
707  endif
708 
709  if (bcast) then
710  select type(buf)
711  type is (integer(kind=i4_kind))
712  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
713  type is (integer(kind=i8_kind))
714  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
715  type is (real(kind=r4_kind))
716  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
717  type is (real(kind=r8_kind))
718  call mpp_broadcast(buf, size(buf), fileobj%io_root, pelist=fileobj%pelist)
719  class default
720  call error("Unsupported variable type: "//trim(append_error_msg))
721  end select
722  endif
723 end subroutine netcdf_read_data_5d
724 !> @}
subroutine unpack_data_2d(fileobj, varid, varname, var_data)
subroutine unpack_data_1d(fileobj, varid, varname, var_data)
Definition: unpack_data.inc:68
subroutine netcdf_read_data_1d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine unpack_data_5d(fileobj, varid, varname, var_data)
subroutine unpack_data_0d(fileobj, varid, varname, var_data)
Definition: unpack_data.inc:24
subroutine netcdf_read_data_5d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine netcdf_read_data_3d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
Character read 1d function.
subroutine netcdf_read_data_2d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine unpack_data_3d(fileobj, varid, varname, var_data)
subroutine char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
Character read 0d function.
subroutine unpack_data_4d(fileobj, varid, varname, var_data)
subroutine netcdf_read_data_0d(fileobj, variable_name, buf, unlim_dim_level, corner, broadcast)
Read in data from a variable in a netcdf file.
subroutine netcdf_read_data_4d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.