FMS  2024.03
Flexible Modeling System
netcdf_io.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup netcdf_io_mod netcdf_io_mod
20 !> @ingroup fms2_io
21 !> @brief Creates a basic netcdf type and routines to extend for other uses
22 !!
23 !> Handles type definitions and interfaces for netcdf I/O.
24 
25 !> @addtogroup netcdf_io_mod
26 !> @{
27 module netcdf_io_mod
28 #ifndef MAX_NUM_RESTART_VARS_
29 #define MAX_NUM_RESTART_VARS_ 250
30 #endif
31 use netcdf
32 use mpp_mod
33 use fms_io_utils_mod
34 use platform_mod
35 implicit none
36 private
37 
38 
39 !Module constants.
40 integer, parameter, public :: default_deflate_level = 0 !< The default (no compression) deflate level to use
41 integer, parameter :: variable_missing = -1
42 integer, parameter :: dimension_missing = -1
43 integer, parameter, public :: no_unlimited_dimension = -1 !> No unlimited dimension in file.
44 character(len=1), parameter :: missing_path = ""
45 integer, parameter :: missing_ncid = -1
46 integer, parameter :: missing_rank = -1
47 integer, parameter, public :: define_mode = 0
48 integer, parameter, public :: data_mode = 1
49 integer, parameter, public :: max_num_restart_vars = max_num_restart_vars_
50 integer, parameter, public :: unlimited = nf90_unlimited !> Wrapper to specify unlimited dimension.
51 integer, parameter :: dimension_not_found = 0
52 integer, parameter, public :: max_num_compressed_dims = 10 !> Maximum number of compressed
53  !! dimensions allowed.
54 integer, private :: fms2_ncchksz = -1 !< Chunksize (bytes) used in nc_open and nc_create
55 integer, private :: fms2_nc_format_param = -1 !< Netcdf format type param used in nc_create
56 character (len = 10), private :: fms2_nc_format !< Netcdf format type used in netcdf_file_open
57 integer, private :: fms2_header_buffer_val = -1 !< value used in NF__ENDDEF
58 integer, private :: fms2_deflate_level = default_deflate_level !< Netcdf deflate level to use in
59  !! nf90_def_var (integer between 1 to 9)
60 logical, private :: fms2_shuffle = .false. !< Flag indicating whether to use the netcdf shuffle filter
61 logical, private :: fms2_is_netcdf4 = .false. !< Flag indicating whether the default netcdf file format is netcdf4
62 
63 !> @}
64 
65 !> @brief information needed fr regional restart variables
66 !> @ingroup netcdf_io_mod
67 type, private :: bc_information
68  integer, dimension(:), allocatable :: indices !< Indices for the halo region for the variable
69  !! (starting x, ending x, starting y, ending y)
70  integer, dimension(:), allocatable :: global_size !< Size of the variable for each dimension
71  integer, dimension(:), allocatable :: pelist !< List of pelist that have the data for the variable
72  logical :: is_root_pe !< Flag indicating if this is the root_pe from the pelist
73  integer :: x_halo !< Number of halos in x
74  integer :: y_halo !< Number of halos in y
75  integer :: jshift !< Shift in the x axis (from center)
76  integer :: ishift !< Shift in the y axis (from center)
77  real(kind=r4_kind), dimension(:,:), allocatable :: globaldata2d_r4 !< 2d data pointer.
78  real(kind=r4_kind), dimension(:,:,:), allocatable :: globaldata3d_r4 !< 3d data pointer.
79  real(kind=r8_kind), dimension(:,:), allocatable :: globaldata2d_r8 !< 2d data pointer.
80  real(kind=r8_kind), dimension(:,:,:), allocatable :: globaldata3d_r8 !< 3d data pointer.
81  character(len=32) :: chksum !< The variable's checksum
82  logical :: data_on_file_root !< Flag indicating if the file root is part of the pelist that
83  !!contains data
84 endtype bc_information
85 
86 !> @brief Restart variable.
87 !> @ingroup netcdf_io_mod
88 type, private :: restartvariable_t
89  character(len=256) :: varname !< Variable name.
90  class(*), pointer :: data0d => null() !< Scalar data pointer.
91  class(*), dimension(:), pointer :: data1d => null() !< 1d data pointer.
92  class(*), dimension(:,:), pointer :: data2d => null() !< 2d data pointer.
93  class(*), dimension(:,:,:), pointer :: data3d => null() !< 3d data pointer.
94  class(*), dimension(:,:,:,:), pointer :: data4d => null() !< 4d data pointer.
95  class(*), dimension(:,:,:,:,:), pointer :: data5d => null() !< 5d data pointer.
96  logical :: was_read !< Flag to support legacy "query_initialized" feature, which
97  !! keeps track if a file was read.
98  logical :: is_bc_variable !< Flag indicating if variable is a bc_variable
99  type(bc_information) :: bc_info !< information about the boundary condition variable
100 endtype restartvariable_t
101 
102 !> @brief Compressed dimension.
103 !> @ingroup netcdf_io_mod
104 type, private :: compresseddimension_t
105  character(len=256) :: dimname !< Dimension name.
106  integer, dimension(:), allocatable :: npes_corner !< Array of starting
107  !! indices for each rank.
108  integer, dimension(:), allocatable :: npes_nelems !< Number of elements
109  !! associated with each
110  !! rank.
111  integer :: nelems !< Total size of the dimension.
112 endtype compresseddimension_t
113 
114 !> @brief information about the current dimensions for regional restart variables
115 !> @ingroup netcdf_io_mod
116 type, private :: dimension_information
117  integer, dimension(5) :: xlen !> The size of each unique x dimension
118  integer, dimension(5) :: ylen !> The size of each unique y dimension
119  integer, dimension(5) :: zlen !> The size of each unique z dimension
120  integer, dimension(3) :: cur_dim_len !> Number of unique:
121  !! cur_dim_len(1) : x dimensions
122  !! cur_dim_len(2) : y dimensions
123  !! cur_dim_len(3) : z dimensions
124 endtype dimension_information
125 
126 !> @brief Netcdf file type.
127 !> @ingroup netcdf_io_mod
128 type, public :: fmsnetcdffile_t
129  character(len=FMS_PATH_LEN) :: path !< File path.
130  logical :: is_readonly !< Flag telling if the file is readonly.
131  integer :: ncid !< Netcdf file id.
132  character(len=256) :: nc_format !< Netcdf file format.
133  logical :: is_netcdf4 !< Flag indicating if the netcdf file type is netcdf4
134  integer, dimension(:), allocatable :: pelist !< List of ranks who will
135  !! communicate.
136  integer :: io_root !< I/O root rank of the pelist.
137  logical :: is_root !< Flag telling if the current rank is the
138  !! I/O root.
139  logical :: is_restart !< Flag telling if the this file is a restart
140  !! file (that has internal pointers to data).
141  logical :: mode_is_append !! true if file is open in "append" mode
142  logical, allocatable :: is_open !< Allocated and set to true if opened.
143  type(restartvariable_t), dimension(:), allocatable :: restart_vars !< Array of registered
144  !! restart variables.
145  integer :: num_restart_vars !< Number of registered restart variables.
146  type(compresseddimension_t), dimension(:), allocatable :: compressed_dims !< "Compressed" dimension.
147  integer :: num_compressed_dims !< Number of compressed dimensions.
148  logical :: is_diskless !< Flag telling whether this is a diskless file.
149  character (len=20) :: time_name
150  type(dimension_information) :: bc_dimensions !<information about the current dimensions for regional
151  !! restart variables
152  logical :: use_collective = .false. !< Flag indicating if we should open the file for collective input
153  !! this should be set to .true. in the user application if they want
154  !! collective reads (put before open_file())
155  integer :: tile_comm=mpp_comm_null !< MPI communicator used for collective reads.
156  !! To be replaced with a real communicator at user request
157 
158 endtype fmsnetcdffile_t
159 
160 
161 !> @brief Range type for a netcdf variable.
162 !> @ingroup netcdf_io_mod
163 type, public :: valid_t
164  logical :: has_range !< Flag that's true if both min/max exist for a variable.
165  logical :: has_min !< Flag that's true if min exists for a variable.
166  logical :: has_max !< Flag that's true if max exists for a variable.
167  logical :: has_fill !< Flag that's true a user defined fill value.
168  logical :: has_missing !< Flag that's true a user defined missing value.
169  real(kind=r8_kind) :: fill_val !< Unpacked fill value for a variable.
170  real(kind=r8_kind) :: min_val !< Unpacked minimum value allowed for a variable.
171  real(kind=r8_kind) :: max_val !< Unpacked maximum value allowed for a variable.
172  real(kind=r8_kind) :: missing_val !< Unpacked missing value for a variable.
173 endtype valid_t
174 
175 
176 public :: netcdf_io_init
177 public :: netcdf_file_open
178 public :: netcdf_file_close
179 public :: netcdf_add_dimension
180 public :: netcdf_add_variable
182 public :: global_att_exists
183 public :: variable_att_exists
186 public :: get_global_attribute
187 public :: get_variable_attribute
188 public :: get_num_dimensions
189 public :: get_dimension_names
190 public :: dimension_exists
191 public :: is_dimension_unlimited
192 public :: get_dimension_size
193 public :: get_num_variables
194 public :: get_variable_names
195 public :: variable_exists
198 public :: get_variable_size
200 public :: netcdf_read_data
201 public :: netcdf_write_data
202 public :: compressed_write
203 public :: netcdf_save_restart
204 public :: netcdf_restore_state
205 public :: get_valid
206 public :: is_valid
208 public :: netcdf_file_open_wrap
209 public :: netcdf_file_close_wrap
210 public :: netcdf_add_variable_wrap
211 public :: netcdf_save_restart_wrap
212 public :: compressed_write_0d_wrap
213 public :: compressed_write_1d_wrap
214 public :: compressed_write_2d_wrap
215 public :: compressed_write_3d_wrap
216 public :: compressed_write_4d_wrap
217 public :: compressed_write_5d_wrap
218 public :: compressed_read_0d
219 public :: compressed_read_1d
220 public :: compressed_read_2d
221 public :: compressed_read_3d
222 public :: compressed_read_4d
223 public :: compressed_read_5d
235 public :: get_fill_value
236 public :: get_variable_sense
237 public :: get_variable_missing
238 public :: get_variable_units
239 public :: get_time_calendar
240 public :: is_registered_to_restart
241 public :: set_netcdf_mode
242 public :: check_netcdf_code
243 public :: check_if_open
244 public :: set_fileobj_time_name
245 public :: write_restart_bc
246 public :: read_restart_bc
247 public :: flush_file
248 
249 !> @ingroup netcdf_io_mod
251  module procedure netcdf_add_restart_variable_0d
252  module procedure netcdf_add_restart_variable_1d
253  module procedure netcdf_add_restart_variable_2d
254  module procedure netcdf_add_restart_variable_3d
255  module procedure netcdf_add_restart_variable_4d
256  module procedure netcdf_add_restart_variable_5d
257 end interface netcdf_add_restart_variable
258 
259 !> @ingroup netcdf_io_mod
261  module procedure netcdf_read_data_0d
262  module procedure netcdf_read_data_1d
263  module procedure netcdf_read_data_2d
264  module procedure netcdf_read_data_3d
265  module procedure netcdf_read_data_4d
266  module procedure netcdf_read_data_5d
267 end interface netcdf_read_data
268 
269 
270 !> @ingroup netcdf_io_mod
272  module procedure netcdf_write_data_0d
273  module procedure netcdf_write_data_1d
274  module procedure netcdf_write_data_2d
275  module procedure netcdf_write_data_3d
276  module procedure netcdf_write_data_4d
277  module procedure netcdf_write_data_5d
278 end interface netcdf_write_data
279 
280 
281 !> @ingroup netcdf_io_mod
283  module procedure compressed_write_0d
284  module procedure compressed_write_1d
285  module procedure compressed_write_2d
286  module procedure compressed_write_3d
287  module procedure compressed_write_4d
288  module procedure compressed_write_5d
289 end interface compressed_write
290 
291 
292 !> @ingroup netcdf_io_mod
294  module procedure register_global_attribute_0d
295  module procedure register_global_attribute_1d
296 end interface register_global_attribute
297 
298 
299 !> @ingroup netcdf_io_mod
301  module procedure register_variable_attribute_0d
302  module procedure register_variable_attribute_1d
303 end interface register_variable_attribute
304 
305 
306 !> @ingroup netcdf_io_mod
308  module procedure get_global_attribute_0d
309  module procedure get_global_attribute_1d
310 end interface get_global_attribute
311 
312 
313 !> @ingroup netcdf_io_mod
315  module procedure get_variable_attribute_0d
316  module procedure get_variable_attribute_1d
317 end interface get_variable_attribute
318 
319 
320 !> @ingroup netcdf_io_mod
322  module procedure scatter_data_bc_2d
323  module procedure scatter_data_bc_3d
324 end interface scatter_data_bc
325 
326 !> @ingroup netcdf_io_mod
327 interface gather_data_bc
328  module procedure gather_data_bc_2d
329  module procedure gather_data_bc_3d
330 end interface gather_data_bc
331 
332 !> The interface is needed to accomodate pgi because it can't handle class * and there was no other way around it
333 !> @ingroup netcdf_io_mod
334 interface is_valid
335  module procedure is_valid_r8
336  module procedure is_valid_r4
337 end interface is_valid
338 
339 !> @addtogroup netcdf_io_mod
340 !> @{
341 
342 contains
343 
344 !> @brief Accepts the namelist fms2_io_nml variables relevant to netcdf_io_mod
345 subroutine netcdf_io_init (chksz, header_buffer_val, netcdf_default_format, deflate_level, shuffle)
346 integer, intent(in) :: chksz !< Chunksize (bytes) used in nc_open and nc_create
347 character (len = 10), intent(in) :: netcdf_default_format !< Netcdf format type param used in nc_create
348 integer, intent(in) :: header_buffer_val !< Value used in NF__ENDDEF
349 integer, intent(in) :: deflate_level !< Netcdf deflate level to use in nf90_def_var
350  !! (integer between 1 to 9)
351 logical, intent(in) :: shuffle !< Flag indicating whether to use the netcdf shuffle filter
352 
353  fms2_ncchksz = chksz
354  fms2_deflate_level = deflate_level
355  fms2_shuffle = shuffle
356  fms2_is_netcdf4 = .false.
357  fms2_header_buffer_val = header_buffer_val
358  if (string_compare(netcdf_default_format, "64bit", .true.)) then
359  fms2_nc_format_param = nf90_64bit_offset
360  call string_copy(fms2_nc_format, "64bit")
361  elseif (string_compare(netcdf_default_format, "classic", .true.)) then
362  fms2_nc_format_param = nf90_classic_model
363  call string_copy(fms2_nc_format, "classic")
364  elseif (string_compare(netcdf_default_format, "netcdf4", .true.)) then
365  fms2_nc_format_param = nf90_netcdf4
366  fms2_is_netcdf4 = .true.
367  call string_copy(fms2_nc_format, "netcdf4")
368  else
369  call error("unrecognized netcdf file format "//trim(netcdf_default_format)// &
370  '. The acceptable values are "64bit", "classic", "netcdf4". Check fms2_io_nml: netcdf_default_format')
371  endif
372 
373 end subroutine netcdf_io_init
374 
375 !> @brief Check for errors returned by netcdf.
376 !! @internal
377 subroutine check_netcdf_code(err, msg)
378 
379  integer, intent(in) :: err !< Code returned by netcdf.
380  character(len=*), intent(in) :: msg !< Error message to be appended to the FATAL
381 
382  character(len=80) :: buf
383 
384  if (err .ne. nf90_noerr) then
385  buf = nf90_strerror(err)
386  call error(trim(buf)//": "//trim(msg))
387  endif
388 end subroutine check_netcdf_code
389 
390 
391 !> @brief Switch to the correct netcdf mode.
392 !! @internal
393 subroutine set_netcdf_mode(ncid, mode)
394 
395  integer, intent(in) :: ncid !< Netcdf file id.
396  integer, intent(in) :: mode !< Netcdf file mode.
397 
398  integer :: err
399 
400  if (mode .eq. define_mode) then
401  err = nf90_redef(ncid)
402  if (err .eq. nf90_eindefine .or. err .eq. nf90_eperm) then
403  return
404  endif
405  elseif (mode .eq. data_mode) then
406  if (fms2_header_buffer_val == -1) call error("set_netcdf_mode: fms2_header_buffer_val not set, call fms2_io_init")
407  err = nf90_enddef(ncid, h_minfree=fms2_header_buffer_val)
408  if (err .eq. nf90_enotindefine .or. err .eq. nf90_eperm) then
409  return
410  endif
411  else
412  call error("mode must be either define_mode or data_mode.")
413  endif
414  call check_netcdf_code(err, "set_netcdf_mode")
415 end subroutine set_netcdf_mode
416 
417 
418 !> @brief Get the id of a dimension from its name.
419 !! @return Dimension id, or dimension_missing if it doesn't exist.
420 !! @internal
421 function get_dimension_id(ncid, dimension_name, msg, allow_failure) &
422  result(dimid)
423 
424  integer, intent(in) :: ncid !< Netcdf file id.
425  character(len=*), intent(in) :: dimension_name !< Dimension name.
426  character(len=*), intent(in) :: msg !< Error message
427  logical, intent(in), optional :: allow_failure !< Flag that prevents
428  !! crash if dimension
429  !! does not exist.
430 
431  integer :: dimid
432 
433  integer :: err
434 
435  err = nf90_inq_dimid(ncid, trim(dimension_name), dimid)
436  if (present(allow_failure)) then
437  if (allow_failure .and. err .eq. nf90_ebaddim) then
438  dimid = dimension_missing
439  return
440  endif
441  endif
442  call check_netcdf_code(err, msg)
443 end function get_dimension_id
444 
445 
446 !> @brief Get the id of a variable from its name.
447 !! @return Variable id, or variable_missing if it doesn't exist.
448 !! @internal
449 function get_variable_id(ncid, variable_name, msg, allow_failure) &
450  result(varid)
451 
452  integer, intent(in) :: ncid !< Netcdf file object.
453  character(len=*), intent(in) :: variable_name !< Variable name.
454  character(len=*), intent(in) :: msg !< Error message
455  logical, intent(in), optional :: allow_failure !< Flag that prevents
456  !! crash if variable does
457  !! not exist.
458 
459  integer :: varid
460 
461  integer :: err
462 
463  err = nf90_inq_varid(ncid, trim(variable_name), varid)
464  if (present(allow_failure)) then
465  if (allow_failure .and. err .eq. nf90_enotvar) then
466  varid = variable_missing
467  return
468  endif
469  endif
470  call check_netcdf_code(err, msg)
471 end function get_variable_id
472 
473 
474 !> @brief Determine if an attribute exists.
475 !! @return Flag telling if the attribute exists.
476 !! @internal
477 function attribute_exists(ncid, varid, attribute_name, msg) &
478  result(att_exists)
479 
480  integer, intent(in) :: ncid !< Netcdf file id.
481  integer, intent(in) :: varid !< Variable id.
482  character(len=*), intent(in) :: attribute_name !< Attribute name.
483  character(len=*), intent(in), optional :: msg !< Error message
484 
485  logical :: att_exists
486 
487  integer :: err
488 
489  err = nf90_inquire_attribute(ncid, varid, trim(attribute_name))
490  if (err .eq. nf90_enotatt) then
491  att_exists = .false.
492  else
493  call check_netcdf_code(err, msg)
494  att_exists = .true.
495  endif
496 end function attribute_exists
497 
498 
499 !> @brief Get the type of a netcdf attribute.
500 !! @return The netcdf type of the attribute.
501 !! @internal
502 function get_attribute_type(ncid, varid, attname, msg) &
503  result(xtype)
504 
505  integer, intent(in) :: ncid !< Netcdf file id.
506  integer, intent(in) :: varid !< Variable id.
507  character(len=*), intent(in) :: attname !< Attribute name.
508  character(len=*), intent(in), optional :: msg !< Error message
509 
510  integer :: xtype
511 
512  integer :: err
513 
514  err = nf90_inquire_attribute(ncid, varid, attname, xtype=xtype)
515  call check_netcdf_code(err, msg)
516 end function get_attribute_type
517 
518 
519 !> @brief Get the type of a netcdf variable.
520 !! @return The netcdf type of the variable.
521 !! @internal
522 function get_variable_type(ncid, varid, msg) &
523  result(xtype)
524 
525  integer, intent(in) :: ncid !< Netcdf file id.
526  integer, intent(in) :: varid !< Variable id.
527  character(len=*), intent(in), optional :: msg !< Error message to append to netcdf error code
528 
529  integer :: xtype
530 
531  integer :: err
532 
533  err = nf90_inquire_variable(ncid, varid, xtype=xtype)
534  call check_netcdf_code(err, msg)
535 end function get_variable_type
536 
537 
538 !> @brief Open a netcdf file.
539 !! @return .true. if open succeeds, or else .false.
540 function netcdf_file_open(fileobj, path, mode, nc_format, pelist, is_restart, dont_add_res_to_filename) &
541  result(success)
542 
543  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
544  character(len=*), intent(in) :: path !< File path.
545  character(len=*), intent(in) :: mode !< File mode. Allowed values are:
546  !! "read", "append", "write", or
547  !! "overwrite".
548  character(len=*), intent(in), optional :: nc_format !< Netcdf format that
549  !! new files are written
550  !! as. Allowed values
551  !! are: "64bit", "classic",
552  !! or "netcdf4". Defaults to
553  !! "64bit". This overwrites
554  !! the value set in the fms2io
555  !! namelist
556  integer, dimension(:), intent(in), optional :: pelist !< List of ranks associated
557  !! with this file. If not
558  !! provided, only the current
559  !! rank will be able to
560  !! act on the file.
561  logical, intent(in), optional :: is_restart !< Flag telling if this file
562  !! is a restart file. Defaults
563  !! to false.
564  logical, intent(in), optional :: dont_add_res_to_filename !< Flag indicating not to add
565  !! ".res" to the filename
566  logical :: success
567 
568  integer :: nc_format_param
569  integer :: err
570  integer :: netcdf4 !< Query the file for the _IsNetcdf4 global attribute in the event
571  !! that the open for collective reads fails
572  character(len=FMS_PATH_LEN) :: buf !< File path with .res in the filename if it is a restart
573  character(len=FMS_PATH_LEN) :: buf2 !< File path with the filename appendix if there is one
574  logical :: is_res
575  logical :: dont_add_res !< flag indicated to not add ".res" to the filename
576 
577  if (allocated(fileobj%is_open)) then
578  if (fileobj%is_open) then
579  success = .true.
580  return
581  endif
582  endif
583  !< Only add ".res" to the file path if is_restart is set to true
584  !! and dont_add_res_to_filename is set to false.
585  is_res = .false.
586  if (present(is_restart)) then
587  is_res = is_restart
588  endif
589  dont_add_res = .false.
590  if (present(dont_add_res_to_filename)) then
591  dont_add_res = dont_add_res_to_filename
592  endif
593 
594  if (is_res .and. .not. dont_add_res) then
595  call restart_filepath_mangle(buf, trim(path))
596  else
597  call string_copy(buf, trim(path))
598  endif
599 
600  !< If it is a restart add the filename_appendix to the filename
601  if (is_res) then
602  call get_instance_filename(trim(buf), buf2)
603  else
604  call string_copy(buf2, trim(buf))
605  endif
606 
607  !Check if the file exists.
608  success = .true.
609  if (string_compare(mode, "read", .true.) .or. string_compare(mode, "append", .true.)) then
610  success = file_exists(buf2)
611  if (.not. success) then
612  return
613  endif
614  endif
615 
616  !Store properties in the derived type.
617  call string_copy(fileobj%path, trim(buf2))
618  if (present(pelist)) then
619  allocate(fileobj%pelist(size(pelist)))
620  fileobj%pelist(:) = pelist(:)
621  else
622  allocate(fileobj%pelist(1))
623  fileobj%pelist(1) = mpp_pe()
624  endif
625  fileobj%io_root = fileobj%pelist(1)
626  fileobj%is_root = mpp_pe() .eq. fileobj%io_root
627 
628  fileobj%is_netcdf4 = .false.
629  if (fms2_ncchksz == -1) call error("netcdf_file_open:: fms2_ncchksz not set, call fms2_io_init")
630  if (fms2_nc_format_param == -1) call error("netcdf_file_open:: fms2_nc_format_param not set, call fms2_io_init")
631 
632  if (present(nc_format)) then
633  if (string_compare(nc_format, "64bit", .true.)) then
634  nc_format_param = nf90_64bit_offset
635  elseif (string_compare(nc_format, "classic", .true.)) then
636  nc_format_param = nf90_classic_model
637  elseif (string_compare(nc_format, "netcdf4", .true.)) then
638  fileobj%is_netcdf4 = .true.
639  nc_format_param = nf90_netcdf4
640  else
641  call error("unrecognized netcdf file format: '"//trim(nc_format)//"' for file:"//trim(fileobj%path)//&
642  &"Check your open_file call, the acceptable values are 64bit, classic, netcdf4")
643  endif
644  call string_copy(fileobj%nc_format, nc_format)
645  else
646  call string_copy(fileobj%nc_format, trim(fms2_nc_format))
647  nc_format_param = fms2_nc_format_param
648  fileobj%is_netcdf4 = fms2_is_netcdf4
649  endif
650 
651  !Open the file with netcdf if this rank is the I/O root.
652  if (fileobj%is_root .and. .not.(fileobj%use_collective)) then
653  if (string_compare(mode, "read", .true.)) then
654  err = nf90_open(trim(fileobj%path), nf90_nowrite, fileobj%ncid, chunksize=fms2_ncchksz)
655  elseif (string_compare(mode, "append", .true.)) then
656  err = nf90_open(trim(fileobj%path), nf90_write, fileobj%ncid, chunksize=fms2_ncchksz)
657  elseif (string_compare(mode, "write", .true.)) then
658  err = nf90_create(trim(fileobj%path), ior(nf90_noclobber, nc_format_param), fileobj%ncid, chunksize=fms2_ncchksz)
659  elseif (string_compare(mode,"overwrite",.true.)) then
660  err = nf90_create(trim(fileobj%path), ior(nf90_clobber, nc_format_param), fileobj%ncid, chunksize=fms2_ncchksz)
661  else
662  call error("unrecognized file mode: '"//trim(mode)//"' for file:"//trim(fileobj%path)//&
663  &"Check your open_file call, the acceptable values are read, append, write, overwrite")
664  endif
665  call check_netcdf_code(err, "netcdf_file_open:"//trim(fileobj%path))
666  elseif(fileobj%use_collective .and. (fileobj%tile_comm /= mpp_comm_null)) then
667  if(string_compare(mode, "read", .true.)) then
668  ! Open the file for collective reads if the user requested that treatment in their application.
669  ! NetCDF does not have the ability to specify collective I/O at the file basis
670  ! so we must activate each variable in netcdf_read_data_2d() and netcdf_read_data_3d()
671  err = nf90_open(trim(fileobj%path), ior(nf90_nowrite, nf90_mpiio), fileobj%ncid, &
672  comm=fileobj%tile_comm, info=mpp_info_null)
673  if(err /= nf90_noerr) then
674  err = nf90_open(trim(fileobj%path), nf90_nowrite, fileobj%ncid)
675  err = nf90_get_att(fileobj%ncid, nf90_global, "_IsNetcdf4", netcdf4)
676  err = nf90_close(fileobj%ncid)
677  if(netcdf4 /= 1) then
678  call mpp_error(note,"netcdf_file_open: Open for collective read failed because the file is not &
679  netCDF-4 format."// &
680  " Falling back to parallel independent for file "// trim(fileobj%path))
681  endif
682  err = nf90_open(trim(fileobj%path), nf90_nowrite, fileobj%ncid, chunksize=fms2_ncchksz)
683  endif
684  elseif (string_compare(mode, "write", .true.)) then
685  call mpp_error(fatal,"netcdf_file_open: Attempt to create a file for collective write"// &
686  " This feature is not implemented"// trim(fileobj%path))
687  !err = nf90_create(trim(fileobj%path), ior(nf90_noclobber, nc_format_param), fileobj%ncid, &
688  ! comm=fileobj%tile_comm, info=MPP_INFO_NULL)
689  elseif (string_compare(mode,"overwrite",.true.)) then
690  call mpp_error(fatal,"netcdf_file_open: Attempt to create a file for collective overwrite"// &
691  " This feature is not implemented"// trim(fileobj%path))
692  !err = nf90_create(trim(fileobj%path), ior(nf90_clobber, nc_format_param), fileobj%ncid, &
693  ! comm=fileobj%tile_comm, info=MPP_INFO_NULL)
694  else
695  call error("unrecognized file mode: '"//trim(mode)//"' for file:"//trim(fileobj%path)//&
696  &"Check your open_file call, the acceptable values are read, append, write, overwrite")
697  endif
698  call check_netcdf_code(err, "netcdf_file_open:"//trim(fileobj%path))
699  else
700  fileobj%ncid = missing_ncid
701  endif
702 
703  fileobj%is_diskless = .false.
704 
705  !Allocate memory.
706  fileobj%is_restart = is_res
707  if (fileobj%is_restart) then
708  allocate(fileobj%restart_vars(max_num_restart_vars))
709  fileobj%num_restart_vars = 0
710  endif
711  fileobj%is_readonly = string_compare(mode, "read", .true.)
712  fileobj%mode_is_append = string_compare(mode, "append", .true.)
713  allocate(fileobj%compressed_dims(max_num_compressed_dims))
714  fileobj%num_compressed_dims = 0
715  ! Set the is_open flag to true for this file object.
716  if (.not.allocated(fileobj%is_open)) allocate(fileobj%is_open)
717  fileobj%is_open = .true.
718 
719  fileobj%bc_dimensions%xlen = 0
720  fileobj%bc_dimensions%ylen = 0
721  fileobj%bc_dimensions%zlen = 0
722  fileobj%bc_dimensions%cur_dim_len = 0
723 
724 end function netcdf_file_open
725 
726 
727 !> @brief Close a netcdf file.
728 subroutine netcdf_file_close(fileobj)
729 
730  class(fmsnetcdffile_t),intent(inout) :: fileobj !< File object.
731 
732  integer :: err
733  integer :: i
734 
735  if (fileobj%is_root) then
736  err = nf90_close(fileobj%ncid)
737  call check_netcdf_code(err, "netcdf_file_close:"//trim(fileobj%path))
738  endif
739  if (allocated(fileobj%is_open)) fileobj%is_open = .false.
740  fileobj%path = missing_path
741  fileobj%ncid = missing_ncid
742  if (allocated(fileobj%pelist)) then
743  deallocate(fileobj%pelist)
744  endif
745  fileobj%io_root = missing_rank
746  fileobj%is_root = .false.
747  if (allocated(fileobj%restart_vars)) then
748  deallocate(fileobj%restart_vars)
749  endif
750  fileobj%is_restart = .false.
751  fileobj%num_restart_vars = 0
752  do i = 1, fileobj%num_compressed_dims
753  if (allocated(fileobj%compressed_dims(i)%npes_corner)) then
754  deallocate(fileobj%compressed_dims(i)%npes_corner)
755  endif
756  if (allocated(fileobj%compressed_dims(i)%npes_nelems)) then
757  deallocate(fileobj%compressed_dims(i)%npes_nelems)
758  endif
759  enddo
760  if (allocated(fileobj%compressed_dims)) then
761  deallocate(fileobj%compressed_dims)
762  endif
763 end subroutine netcdf_file_close
764 
765 
766 !> @brief Get the index of a compressed dimension in a file object.
767 !! @return Index of the compressed dimension.
768 !! @internal
769 function get_compressed_dimension_index(fileobj, dim_name) &
770  result(dindex)
771 
772  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
773  character(len=*), intent(in) :: dim_name !< Dimension name.
774 
775  integer :: dindex
776  integer :: i
777 
778  dindex = dimension_not_found
779  do i = 1, fileobj%num_compressed_dims
780  if (string_compare(fileobj%compressed_dims(i)%dimname, dim_name)) then
781  dindex = i
782  return
783  endif
784  enddo
786 
787 
788 !> @brief Add a compressed dimension to a file object.
789 !! @internal
790 subroutine append_compressed_dimension(fileobj, dim_name, npes_corner, &
791  npes_nelems)
792 
793  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
794  character(len=*), intent(in) :: dim_name !< Dimension name.
795  integer, dimension(:), intent(in) :: npes_corner !< Array of starting
796  !! indices for each rank.
797  integer, dimension(:), intent(in) :: npes_nelems !< Number of elements
798  !! associated with each
799  !! rank.
800 
801  integer :: n
802 
803  if (get_compressed_dimension_index(fileobj, dim_name) .ne. dimension_not_found) then
804  call error("dimension "//trim(dim_name)//" already registered" &
805  //" to file "//trim(fileobj%path)//".")
806  endif
807  fileobj%num_compressed_dims = fileobj%num_compressed_dims + 1
808  n = fileobj%num_compressed_dims
809  if (n .gt. max_num_compressed_dims) then
810  call error("number of compressed dimensions exceeds limit.")
811  endif
812  call string_copy(fileobj%compressed_dims(n)%dimname, dim_name)
813  if (size(npes_corner) .ne. size(fileobj%pelist) .or. &
814  size(npes_nelems) .ne. size(fileobj%pelist)) then
815  call error("incorrect size for input npes_corner or npes_nelems arrays.")
816  endif
817  allocate(fileobj%compressed_dims(n)%npes_corner(size(fileobj%pelist)))
818  fileobj%compressed_dims(n)%npes_corner(:) = npes_corner(:)
819  allocate(fileobj%compressed_dims(n)%npes_nelems(size(fileobj%pelist)))
820  fileobj%compressed_dims(n)%npes_nelems(:) = npes_nelems(:)
821  fileobj%compressed_dims(n)%nelems = sum(fileobj%compressed_dims(n)%npes_nelems)
822 end subroutine append_compressed_dimension
823 
824 !> @brief Add a "compressed" unlimited dimension to a netcdf file.
825 !! @note Here compressed means that every rank has a different dimension_length
826 !! compressed. This was written specifically for the icebergs restarts.
827 subroutine register_unlimited_compressed_axis(fileobj, dimension_name, dimension_length)
828  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
829  character(len=*), intent(in) :: dimension_name !< Dimension name.
830  integer, intent(in) :: dimension_length !< Dimension length for the current rank
831 
832  integer, dimension(:), allocatable :: npes_start !< The starting index of the dimension for each of the PEs
833  integer, dimension(:), allocatable :: npes_count !< The size of the dimension for each of the PEs
834  integer :: i !< For do loops
835  integer :: err !< Netcdf error
836  integer :: dimid !< Netcdf id for the dimension
837 
838  !Gather all local dimension lengths on the I/O root pe.
839  allocate(npes_start(size(fileobj%pelist)))
840  allocate(npes_count(size(fileobj%pelist)))
841 
842  call mpp_gather((/dimension_length/),npes_count,pelist=fileobj%pelist)
843 
844  npes_start(1) = 1
845  do i = 1, size(fileobj%pelist)-1
846  npes_start(i+1) = npes_start(i) + npes_count(i)
847  enddo
848  call append_compressed_dimension(fileobj, dimension_name, npes_start, &
849  npes_count)
850 
851  if (fileobj%is_root .and. .not. fileobj%is_readonly) then
852  call set_netcdf_mode(fileobj%ncid, define_mode)
853  err = nf90_def_dim(fileobj%ncid, trim(dimension_name), unlimited, dimid)
854  call check_netcdf_code(err, "Netcdf_add_dimension: file:"//trim(fileobj%path)//" dimension name:"// &
855  & trim(dimension_name))
856  endif
858 
859 !> @brief Add a dimension to a file.
860 subroutine netcdf_add_dimension(fileobj, dimension_name, dimension_length, &
861  is_compressed)
862 
863  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
864  character(len=*), intent(in) :: dimension_name !< Dimension name.
865  integer, intent(in) :: dimension_length !< Dimension length.
866  logical, intent(in), optional :: is_compressed !< Changes the meaning of dim_len from
867  !! referring to the total size of the
868  !! dimension (when false) to the local
869  !! size for the current rank (when true).
870 
871  integer :: dim_len
872  integer, dimension(:), allocatable :: npes_start
873  integer, dimension(:), allocatable :: npes_count
874  integer :: i
875  integer :: err
876  integer :: dimid
877 
878  dim_len = dimension_length
879  if (present(is_compressed)) then
880  if (is_compressed) then
881  !Gather all local dimension lengths on the I/O root pe.
882  allocate(npes_start(size(fileobj%pelist)))
883  allocate(npes_count(size(fileobj%pelist)))
884  do i = 1, size(fileobj%pelist)
885  if (fileobj%pelist(i) .eq. mpp_pe()) then
886  npes_count(i) = dim_len
887  else
888  call mpp_recv(npes_count(i), fileobj%pelist(i), block=.false.)
889  call mpp_send(dim_len, fileobj%pelist(i))
890  endif
891  enddo
892  call mpp_sync_self(check=event_recv)
893  call mpp_sync_self(check=event_send)
894  npes_start(1) = 1
895  do i = 1, size(fileobj%pelist)-1
896  npes_start(i+1) = npes_start(i) + npes_count(i)
897  enddo
898  call append_compressed_dimension(fileobj, dimension_name, npes_start, &
899  npes_count)
900  dim_len = sum(npes_count)
901  endif
902  endif
903  if (fileobj%is_root .and. .not. fileobj%is_readonly) then
904  call set_netcdf_mode(fileobj%ncid, define_mode)
905  err = nf90_def_dim(fileobj%ncid, trim(dimension_name), dim_len, dimid)
906  call check_netcdf_code(err, "Netcdf_add_dimension: file:"//trim(fileobj%path)//" dimension name:"// &
907  & trim(dimension_name))
908  endif
909 end subroutine netcdf_add_dimension
910 
911 
912 !> @brief Add a compressed dimension.
913 subroutine register_compressed_dimension(fileobj, dimension_name, &
914  npes_corner, npes_nelems)
915 
916  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
917  character(len=*), intent(in) :: dimension_name !< Dimension name.
918  integer, dimension(:), intent(in) :: npes_corner !< Array of starting
919  !! indices for each rank.
920  integer, dimension(:), intent(in) :: npes_nelems !< Number of elements
921  !! associated with each
922  !! rank.
923 
924  integer :: dsize
925  integer :: fdim_size
926 
927  call append_compressed_dimension(fileobj, dimension_name, npes_corner, npes_nelems)
928  dsize = sum(npes_nelems)
929  if (fileobj%is_readonly) then
930  call get_dimension_size(fileobj, dimension_name, fdim_size, broadcast=.true.)
931  if (fdim_size .ne. dsize) then
932  call error("dimension "//trim(dimension_name)//" does not match" &
933  //" the size of the associated compressed axis.")
934  endif
935  else
936  call netcdf_add_dimension(fileobj, dimension_name, dsize)
937  endif
938 end subroutine register_compressed_dimension
939 
940 
941 !> @brief Add a variable to a file.
942 subroutine netcdf_add_variable(fileobj, variable_name, variable_type, dimensions, chunksizes)
943 
944  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
945  character(len=*), intent(in) :: variable_name !< Variable name.
946  character(len=*), intent(in) :: variable_type !< Variable type. Allowed
947  !! values are: "char", "int", "int64",
948  !! "float", or "double".
949  character(len=*), dimension(:), intent(in), optional :: dimensions !< Dimension names.
950  integer, optional, intent(in) :: chunksizes(:) !< netcdf chunksize to use for this variable
951  !! This feature is only
952  !! available for netcdf4 files
953  integer :: err
954  integer, dimension(:), allocatable :: dimids
955  integer :: vtype
956  integer :: varid
957  integer :: i
958  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
959 
960  append_error_msg = "netcdf_add_variable: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
961 
962  if (fileobj%is_root) then
963  call set_netcdf_mode(fileobj%ncid, define_mode)
964  if (string_compare(variable_type, "int", .true.)) then
965  vtype = nf90_int
966  elseif (string_compare(variable_type, "int64", .true.)) then
967  if ( .not. fileobj%is_netcdf4) call error(trim(fileobj%path)//&
968  &": 64 bit integers are only supported with 'netcdf4' file format"//&
969  &". Set netcdf_default_format='netcdf4' in the fms2_io namelist OR "//&
970  &"add nc_format='netcdf4' to your open_file call")
971  vtype = nf90_int64
972  elseif (string_compare(variable_type, "float", .true.)) then
973  vtype = nf90_float
974  elseif (string_compare(variable_type, "double", .true.)) then
975  vtype = nf90_double
976  elseif (string_compare(variable_type, "char", .true.)) then
977  vtype = nf90_char
978  if (.not. present(dimensions)) then
979  call error("String variables require a string length dimension:"//trim(append_error_msg))
980  endif
981  else
982  call error("Unsupported variable type:"//trim(append_error_msg))
983  endif
984  if (present(dimensions)) then
985  allocate(dimids(size(dimensions)))
986  do i = 1, size(dimids)
987  dimids(i) = get_dimension_id(fileobj%ncid, trim(dimensions(i)),msg=append_error_msg)
988  enddo
989  if (fileobj%is_netcdf4) then
990  err = nf90_def_var(fileobj%ncid, trim(variable_name), vtype, dimids, varid, &
991  &deflate_level=fms2_deflate_level, shuffle=fms2_shuffle, chunksizes=chunksizes)
992  else
993  if (fms2_deflate_level .ne. default_deflate_level .or. fms2_shuffle .or. present(chunksizes)) &
994  &call mpp_error(note,"Not able to use deflate_level or chunksizes if not using netcdf4"// &
995  & " ignoring them")
996  err = nf90_def_var(fileobj%ncid, trim(variable_name), vtype, dimids, varid)
997  endif
998  deallocate(dimids)
999  else
1000  err = nf90_def_var(fileobj%ncid, trim(variable_name), vtype, varid)
1001  endif
1002  call check_netcdf_code(err, append_error_msg)
1003  endif
1004 end subroutine netcdf_add_variable
1005 
1006 
1007 !> @brief Given a compressed variable, get the index of the compressed
1008 !! dimension.
1009 !! @return Index of the compressed dimension or dimension_not_found.
1010 function get_variable_compressed_dimension_index(fileobj, variable_name, broadcast) &
1011  result(compressed_dimension_index)
1012 
1013  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1014  character(len=*), intent(in) :: variable_name !< Variable name.
1015  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1016  !! not the data will be
1017  !! broadcasted to non
1018  !! "I/O root" ranks.
1019  !! The broadcast will be done
1020  !! by default.
1021  integer, dimension(2) :: compressed_dimension_index
1022 
1023  integer :: ndims
1024  character(len=nf90_max_name), dimension(:), allocatable :: dim_names
1025  integer :: i
1026  integer :: j
1027 
1028  compressed_dimension_index = dimension_not_found
1029  if (fileobj%is_root) then
1030  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
1031  if (ndims .gt. 0) then
1032  allocate(dim_names(ndims))
1033  call get_variable_dimension_names(fileobj, variable_name, dim_names, broadcast=.false.)
1034  do i = 1, size(dim_names)
1035  j = get_compressed_dimension_index(fileobj,dim_names(i))
1036  if (j .ne. dimension_not_found) then
1037  compressed_dimension_index(1) = i
1038  compressed_dimension_index(2) = j
1039  exit
1040  endif
1041  enddo
1042  deallocate(dim_names)
1043  endif
1044  endif
1045  if (present(broadcast)) then
1046  if (.not. broadcast) then
1047  return
1048  endif
1049  endif
1050  call mpp_broadcast(compressed_dimension_index(1), fileobj%io_root, pelist=fileobj%pelist)
1051  call mpp_broadcast(compressed_dimension_index(2), fileobj%io_root, pelist=fileobj%pelist)
1053 
1054 
1055 !> @brief Add a restart variable to a FmsNetcdfFile_t type.
1056 !! @internal
1057 subroutine add_restart_var_to_array(fileobj, variable_name)
1058 
1059  class(fmsnetcdffile_t), intent(inout) :: fileobj !< Netcdf file object.
1060  character(len=*), intent(in) :: variable_name !< Variable name.
1061 
1062  integer :: i
1063 
1064  if (.not. fileobj%is_restart) then
1065  call error("file "//trim(fileobj%path)//" is not a restart file.")
1066  endif
1067  do i = 1, fileobj%num_restart_vars
1068  if (string_compare(fileobj%restart_vars(i)%varname, variable_name, .true.)) then
1069  call error("variable "//trim(variable_name)//" has already" &
1070  //" been added to restart file "//trim(fileobj%path)//".")
1071  endif
1072  enddo
1073  fileobj%num_restart_vars = fileobj%num_restart_vars + 1
1074  if (fileobj%num_restart_vars .gt. max_num_restart_vars) then
1075  call error("Number of restart variables exceeds limit.")
1076  endif
1077  call string_copy(fileobj%restart_vars(fileobj%num_restart_vars)%varname, &
1078  variable_name)
1079 end subroutine add_restart_var_to_array
1080 
1081 
1082 !> @brief Loop through registered restart variables and write them to
1083 !! a netcdf file.
1084 subroutine netcdf_save_restart(fileobj, unlim_dim_level)
1085 
1086  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1087  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
1088  !! level.
1089 
1090  integer :: i
1091 
1092  if (.not. fileobj%is_restart) then
1093  call error("write_restart:: file "//trim(fileobj%path)//" is not a restart file. &
1094  &Be sure the file was opened with is_restart=.true.")
1095  endif
1096  do i = 1, fileobj%num_restart_vars
1097  if (associated(fileobj%restart_vars(i)%data0d)) then
1098  call compressed_write(fileobj, fileobj%restart_vars(i)%varname, &
1099  fileobj%restart_vars(i)%data0d, &
1100  unlim_dim_level=unlim_dim_level)
1101  elseif (associated(fileobj%restart_vars(i)%data1d)) then
1102  call compressed_write(fileobj, fileobj%restart_vars(i)%varname, &
1103  fileobj%restart_vars(i)%data1d, &
1104  unlim_dim_level=unlim_dim_level)
1105  elseif (associated(fileobj%restart_vars(i)%data2d)) then
1106  call compressed_write(fileobj, fileobj%restart_vars(i)%varname, &
1107  fileobj%restart_vars(i)%data2d, &
1108  unlim_dim_level=unlim_dim_level)
1109  elseif (associated(fileobj%restart_vars(i)%data3d)) then
1110  call compressed_write(fileobj, fileobj%restart_vars(i)%varname, &
1111  fileobj%restart_vars(i)%data3d, &
1112  unlim_dim_level=unlim_dim_level)
1113  elseif (associated(fileobj%restart_vars(i)%data4d)) then
1114  call compressed_write(fileobj, fileobj%restart_vars(i)%varname, &
1115  fileobj%restart_vars(i)%data4d, &
1116  unlim_dim_level=unlim_dim_level)
1117  else
1118  call error("this branch should not be reached.")
1119  endif
1120  enddo
1121 end subroutine netcdf_save_restart
1122 
1123 
1124 !> @brief Loop through registered restart variables and read them from
1125 !! a netcdf file.
1126 subroutine netcdf_restore_state(fileobj, unlim_dim_level)
1127 
1128  type(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
1129  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
1130  !! level.
1131 
1132  integer :: i
1133 
1134  if (.not. fileobj%is_restart) then
1135  call error("read_restart:: file "//trim(fileobj%path)//" is not a restart file. &
1136  &Be sure the file was opened with is_restart=.true.")
1137  endif
1138  do i = 1, fileobj%num_restart_vars
1139  if (associated(fileobj%restart_vars(i)%data0d)) then
1140  call netcdf_read_data(fileobj, fileobj%restart_vars(i)%varname, &
1141  fileobj%restart_vars(i)%data0d, &
1142  unlim_dim_level=unlim_dim_level, &
1143  broadcast=.true.)
1144  elseif (associated(fileobj%restart_vars(i)%data1d)) then
1145  call netcdf_read_data(fileobj, fileobj%restart_vars(i)%varname, &
1146  fileobj%restart_vars(i)%data1d, &
1147  unlim_dim_level=unlim_dim_level, &
1148  broadcast=.true.)
1149  elseif (associated(fileobj%restart_vars(i)%data2d)) then
1150  call netcdf_read_data(fileobj, fileobj%restart_vars(i)%varname, &
1151  fileobj%restart_vars(i)%data2d, &
1152  unlim_dim_level=unlim_dim_level, &
1153  broadcast=.true.)
1154  elseif (associated(fileobj%restart_vars(i)%data3d)) then
1155  call netcdf_read_data(fileobj, fileobj%restart_vars(i)%varname, &
1156  fileobj%restart_vars(i)%data3d, &
1157  unlim_dim_level=unlim_dim_level, &
1158  broadcast=.true.)
1159  elseif (associated(fileobj%restart_vars(i)%data4d)) then
1160  call netcdf_read_data(fileobj, fileobj%restart_vars(i)%varname, &
1161  fileobj%restart_vars(i)%data4d, &
1162  unlim_dim_level=unlim_dim_level, &
1163  broadcast=.true.)
1164  else
1165  call error("this branch should not be reached.")
1166  endif
1167  enddo
1168 end subroutine netcdf_restore_state
1169 
1170 
1171 !> @brief Determine if a global attribute exists.
1172 !! @return Flag telling if a global attribute exists.
1173 function global_att_exists(fileobj, attribute_name, broadcast) &
1174  result(att_exists)
1175 
1176  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1177  character(len=*), intent(in) :: attribute_name !< Attribute name.
1178  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1179  !! not the data will be
1180  !! broadcasted to non
1181  !! "I/O root" ranks.
1182  !! The broadcast will be done
1183  !! by default.
1184  logical :: att_exists
1185 
1186  if (fileobj%is_root) then
1187  att_exists = attribute_exists(fileobj%ncid, nf90_global, trim(attribute_name), &
1188  & msg="global_att_exists: file:"//trim(fileobj%path)//" attribute name:"//trim(attribute_name))
1189  endif
1190  if (present(broadcast)) then
1191  if (.not. broadcast) then
1192  return
1193  endif
1194  endif
1195  call mpp_broadcast(att_exists, fileobj%io_root, pelist=fileobj%pelist)
1196 end function global_att_exists
1197 
1198 
1199 !> @brief Determine if a variable's attribute exists.
1200 !! @return Flag telling if the variable's attribute exists.
1201 function variable_att_exists(fileobj, variable_name, attribute_name, &
1202  broadcast) &
1203  result(att_exists)
1204 
1205  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1206  character(len=*), intent(in) :: variable_name !< Variable name.
1207  character(len=*), intent(in) :: attribute_name !< Attribute name.
1208  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1209  !! not the data will be
1210  !! broadcasted to non
1211  !! "I/O root" ranks.
1212  !! The broadcast will be done
1213  !! by default.
1214  logical :: att_exists
1215 
1216  integer :: varid
1217 
1218  att_exists = .false.
1219  if (fileobj%is_root) then
1220  varid = get_variable_id(fileobj%ncid, trim(variable_name), &
1221  & msg="variable_att_exists: file:"//trim(fileobj%path)//"- variable:"//&
1222  &trim(variable_name))
1223  att_exists = attribute_exists(fileobj%ncid, varid, trim(attribute_name), &
1224  &msg="variable_att_exists: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)//&
1225  &" attribute name:"//trim(attribute_name))
1226  endif
1227  if (present(broadcast)) then
1228  if (.not. broadcast) then
1229  return
1230  endif
1231  endif
1232  call mpp_broadcast(att_exists, fileobj%io_root, pelist=fileobj%pelist)
1233 end function variable_att_exists
1234 
1235 
1236 !> @brief Determine the number of dimensions in a file.
1237 !! @return The number of dimensions in the file.
1238 function get_num_dimensions(fileobj, broadcast) &
1239  result(ndims)
1240 
1241  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1242  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1243  !! not the data will be
1244  !! broadcasted to non
1245  !! "I/O root" ranks.
1246  !! The broadcast will be done
1247  !! by default.
1248  integer :: ndims
1249 
1250  integer :: err
1251 
1252  if (fileobj%is_root) then
1253  err = nf90_inquire(fileobj%ncid, ndimensions=ndims)
1254  call check_netcdf_code(err, "get_num_dimensions: file:"//trim(fileobj%path))
1255  endif
1256  if (present(broadcast)) then
1257  if (.not. broadcast) then
1258  return
1259  endif
1260  endif
1261  call mpp_broadcast(ndims, fileobj%io_root, pelist=fileobj%pelist)
1262 end function get_num_dimensions
1263 
1264 
1265 !> @brief Get the names of the dimensions in a file.
1266 subroutine get_dimension_names(fileobj, names, broadcast)
1267 
1268  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1269  character(len=*), dimension(:), intent(inout) :: names !< Names of the
1270  !! dimensions.
1271  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1272  !! not the data will be
1273  !! broadcasted to non
1274  !! "I/O root" ranks.
1275  !! The broadcast will be done
1276  !! by default.
1277 
1278  integer :: ndims
1279  integer :: i
1280  integer :: err
1281 
1282  if (fileobj%is_root) then
1283  ndims = get_num_dimensions(fileobj, broadcast=.false.)
1284  if (ndims .gt. 0) then
1285  if (size(names) .ne. ndims) then
1286  call error("'names' has to be the same size of the number of dimensions. &
1287  &Check your get_dimension_names call for file "//trim(fileobj%path))
1288  endif
1289  else
1290  call error("get_dimension_names: the file "//trim(fileobj%path)//" does not have any dimensions")
1291  endif
1292  names(:) = ""
1293  do i = 1, ndims
1294  err = nf90_inquire_dimension(fileobj%ncid, i, name=names(i))
1295  call check_netcdf_code(err, "get_dimension_names: file:"//trim(fileobj%path))
1296  enddo
1297  endif
1298  if (present(broadcast)) then
1299  if (.not. broadcast) then
1300  return
1301  endif
1302  endif
1303  call mpp_broadcast(ndims, fileobj%io_root, pelist=fileobj%pelist)
1304  if (.not. fileobj%is_root) then
1305  if (ndims .gt. 0) then
1306  if (size(names) .ne. ndims) then
1307  call error("'names' has to be the same size of the number of dimensions. &
1308  &Check your get_dimension_names call for file "//trim(fileobj%path))
1309  endif
1310  else
1311  call error("get_dimension_names: the file "//trim(fileobj%path)//" does not have any dimensions")
1312  endif
1313  names(:) = ""
1314  endif
1315  call mpp_broadcast(names, len(names(ndims)), fileobj%io_root, &
1316  pelist=fileobj%pelist)
1317 end subroutine get_dimension_names
1318 
1319 
1320 !> @brief Determine if a dimension exists.
1321 !! @return Flag telling if the dimension exists.
1322 function dimension_exists(fileobj, dimension_name, broadcast) &
1323  result(dim_exists)
1324 
1325  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1326  character(len=*), intent(in) :: dimension_name !< Dimension name.
1327  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1328  !! not the data will be
1329  !! broadcasted to non
1330  !! "I/O root" ranks.
1331  !! The broadcast will be done
1332  !! by default.
1333  logical :: dim_exists
1334 
1335  integer :: dimid
1336 
1337  if (fileobj%is_root) then
1338  dimid = get_dimension_id(fileobj%ncid, trim(dimension_name), &
1339  msg="dimension_exists: file:"//trim(fileobj%path)//" dimension:"//trim(dimension_name), &
1340  allow_failure=.true.)
1341  if (dimid .eq. dimension_missing) then
1342  dim_exists = .false.
1343  else
1344  dim_exists = .true.
1345  endif
1346  endif
1347  if (present(broadcast)) then
1348  if (.not. broadcast) then
1349  return
1350  endif
1351  endif
1352  call mpp_broadcast(dim_exists, fileobj%io_root, pelist=fileobj%pelist)
1353 end function dimension_exists
1354 
1355 
1356 !> @brief Determine where or not the dimension is unlimited.
1357 !! @return True if the dimension is unlimited, or else false.
1358 function is_dimension_unlimited(fileobj, dimension_name, broadcast) &
1359  result(is_unlimited)
1360 
1361  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1362  character(len=*), intent(in) :: dimension_name !< Dimension name.
1363  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1364  !! not the data will be
1365  !! broadcasted to non
1366  !! "I/O root" ranks.
1367  !! The broadcast will be done
1368  !! by default.
1369  logical :: is_unlimited
1370 
1371  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1372  integer :: dimid
1373  integer :: err
1374  integer :: ulim_dimid
1375 
1376  if (fileobj%is_root) then
1377  append_error_msg="is_dimension_unlimited: file:"//trim(fileobj%path)//&
1378  & " dimension_name:"//trim(dimension_name)
1379  dimid = get_dimension_id(fileobj%ncid, trim(dimension_name), msg=append_error_msg)
1380  err = nf90_inquire(fileobj%ncid, unlimiteddimid=ulim_dimid)
1381  call check_netcdf_code(err, append_error_msg)
1382  is_unlimited = dimid .eq. ulim_dimid
1383  endif
1384  if (present(broadcast)) then
1385  if (.not. broadcast) then
1386  return
1387  endif
1388  endif
1389  call mpp_broadcast(is_unlimited, fileobj%io_root, pelist=fileobj%pelist)
1390 end function is_dimension_unlimited
1391 
1392 
1393 !> @brief Get the name of the unlimited dimension.
1394 subroutine get_unlimited_dimension_name(fileobj, dimension_name, broadcast)
1395 
1396  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1397  character(len=*), intent(out) :: dimension_name !< Dimension name.
1398  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1399  !! not the data will be
1400  !! broadcasted to non
1401  !! "I/O root" ranks.
1402  !! The broadcast will be done
1403  !! by default.
1404 
1405  integer :: err
1406  integer :: dimid
1407  character(len=nf90_max_name), dimension(1) :: buffer
1408 
1409  dimension_name = ""
1410  if (fileobj%is_root) then
1411  err = nf90_inquire(fileobj%ncid, unlimiteddimid=dimid)
1412  call check_netcdf_code(err, "get_unlimited_dimension_name: file:"//trim(fileobj%path))
1413  err = nf90_inquire_dimension(fileobj%ncid, dimid, dimension_name)
1414  call check_netcdf_code(err, "get_unlimited_dimension_name: file:"//trim(fileobj%path))
1415  call string_copy(buffer(1), dimension_name)
1416  endif
1417  if (present(broadcast)) then
1418  if (.not. broadcast) then
1419  return
1420  endif
1421  endif
1422  call mpp_broadcast(buffer, nf90_max_name, fileobj%io_root, &
1423  pelist=fileobj%pelist)
1424  call string_copy(dimension_name, buffer(1))
1425 end subroutine get_unlimited_dimension_name
1426 
1427 
1428 !> @brief Get the length of a dimension.
1429 subroutine get_dimension_size(fileobj, dimension_name, dim_size, broadcast)
1430 
1431  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1432  character(len=*), intent(in) :: dimension_name !< Dimension name.
1433  integer, intent(inout) :: dim_size !< Dimension size.
1434  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1435  !! not the data will be
1436  !! broadcasted to non
1437  !! "I/O root" ranks.
1438  !! The broadcast will be done
1439  !! by default.
1440 
1441  integer :: dimid
1442  integer :: err
1443  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1444 
1445  if (fileobj%is_root) then
1446  append_error_msg = "get_dimension_size: file:"//trim(fileobj%path)//" dimension_name: "//trim(dimension_name)
1447  dimid = get_dimension_id(fileobj%ncid, trim(dimension_name), msg=append_error_msg)
1448  err = nf90_inquire_dimension(fileobj%ncid, dimid, len=dim_size)
1449  call check_netcdf_code(err, append_error_msg)
1450  endif
1451  if (present(broadcast)) then
1452  if (.not. broadcast) then
1453  return
1454  endif
1455  endif
1456  call mpp_broadcast(dim_size, fileobj%io_root, pelist=fileobj%pelist)
1457 end subroutine get_dimension_size
1458 
1459 
1460 !> @brief Determine the number of variables in a file.
1461 !! @return The number of variables in the file.
1462 function get_num_variables(fileobj, broadcast) &
1463  result(nvars)
1464 
1465  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1466  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1467  !! not the data will be
1468  !! broadcasted to non
1469  !! "I/O root" ranks.
1470  !! The broadcast will be done
1471  !! by default.
1472  integer :: nvars
1473 
1474  integer :: err
1475 
1476  if (fileobj%is_root) then
1477  err = nf90_inquire(fileobj%ncid, nvariables=nvars)
1478  call check_netcdf_code(err, "get_num_variables: file: "//trim(fileobj%path))
1479  endif
1480  if (present(broadcast)) then
1481  if (.not. broadcast) then
1482  return
1483  endif
1484  endif
1485  call mpp_broadcast(nvars, fileobj%io_root, pelist=fileobj%pelist)
1486 end function get_num_variables
1487 
1488 
1489 !> @brief Get the names of the variables in a file.
1490 subroutine get_variable_names(fileobj, names, broadcast)
1491 
1492  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1493  character(len=*), dimension(:), intent(inout) :: names !< Names of the
1494  !! variables.
1495  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1496  !! not the data will be
1497  !! broadcasted to non
1498  !! "I/O root" ranks.
1499  !! The broadcast will be done
1500  !! by default.
1501 
1502  integer :: nvars
1503  integer :: i
1504  integer :: err
1505 
1506  if (fileobj%is_root) then
1507  nvars = get_num_variables(fileobj, broadcast=.false.)
1508  if (nvars .gt. 0) then
1509  if (size(names) .ne. nvars) then
1510  call error("'names' has to be the same size of the number of variables. &
1511  &Check your get_variable_names call for file "//trim(fileobj%path))
1512  endif
1513  else
1514  call error("get_variable_names: the file "//trim(fileobj%path)//" does not have any variables")
1515  endif
1516  names(:) = ""
1517  do i = 1, nvars
1518  err = nf90_inquire_variable(fileobj%ncid, i, name=names(i))
1519  call check_netcdf_code(err, "get_variable_names: "//trim(fileobj%path))
1520  enddo
1521  endif
1522  if (present(broadcast)) then
1523  if (.not. broadcast) then
1524  return
1525  endif
1526  endif
1527  call mpp_broadcast(nvars, fileobj%io_root, pelist=fileobj%pelist)
1528  if (.not. fileobj%is_root) then
1529  if (nvars .gt. 0) then
1530  if (size(names) .ne. nvars) then
1531  call error("'names' has to be the same size of the number of variables. &
1532  &Check your get_variable_names call for file "//trim(fileobj%path))
1533  endif
1534  else
1535  call error("get_variable_names: the file "//trim(fileobj%path)//" does not have any variables")
1536  endif
1537  names(:) = ""
1538  endif
1539  call mpp_broadcast(names, len(names(nvars)), fileobj%io_root, &
1540  pelist=fileobj%pelist)
1541 end subroutine get_variable_names
1542 
1543 
1544 !> @brief Determine if a variable exists.
1545 !! @return Flag telling if the variable exists.
1546 function variable_exists(fileobj, variable_name, broadcast) &
1547  result(var_exists)
1548 
1549  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1550  character(len=*), intent(in) :: variable_name !< Variable name.
1551  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1552  !! not the data will be
1553  !! broadcasted to non
1554  !! "I/O root" ranks.
1555  !! The broadcast will be done
1556  !! by default.
1557  logical :: var_exists
1558 
1559  integer :: varid
1560 
1561  if (fileobj%is_root) then
1562  varid = get_variable_id(fileobj%ncid, trim(variable_name), &
1563  msg="variable_exists: file:"//trim(fileobj%path)//" variable:"//trim(variable_name), &
1564  allow_failure=.true.)
1565  var_exists = varid .ne. variable_missing
1566  endif
1567  if (present(broadcast)) then
1568  if (.not. broadcast) then
1569  return
1570  endif
1571  endif
1572  call mpp_broadcast(var_exists, fileobj%io_root, pelist=fileobj%pelist)
1573 end function variable_exists
1574 
1575 
1576 !> @brief Get the number of dimensions a variable depends on.
1577 !! @return Number of dimensions that the variable depends on.
1578 function get_variable_num_dimensions(fileobj, variable_name, broadcast) &
1579  result(ndims)
1580 
1581  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1582  character(len=*), intent(in) :: variable_name !< Variable name.
1583  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1584  !! not the data will be
1585  !! broadcasted to non
1586  !! "I/O root" ranks.
1587  !! The broadcast will be done
1588  !! by default.
1589  integer :: ndims
1590 
1591  integer :: varid
1592  integer :: err
1593  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1594 
1595 
1596  if (fileobj%is_root) then
1597  append_error_msg = "get_variable_num_dimension: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
1598  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
1599  err = nf90_inquire_variable(fileobj%ncid, varid, ndims=ndims)
1600  call check_netcdf_code(err, append_error_msg)
1601  endif
1602  if (present(broadcast)) then
1603  if (.not. broadcast) then
1604  return
1605  endif
1606  endif
1607  call mpp_broadcast(ndims, fileobj%io_root, pelist=fileobj%pelist)
1608 end function get_variable_num_dimensions
1609 
1610 
1611 !> @brief Get the name of a variable's dimensions.
1612 subroutine get_variable_dimension_names(fileobj, variable_name, dim_names, &
1613  broadcast)
1614 
1615  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1616  character(len=*), intent(in) :: variable_name !< Variable name.
1617  character(len=*), dimension(:), intent(inout) :: dim_names !< Array of
1618  !! dimension
1619  !! names.
1620  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1621  !! not the data will be
1622  !! broadcasted to non
1623  !! "I/O root" ranks.
1624  !! The broadcast will be done
1625  !! by default.
1626 
1627  integer :: varid
1628  integer :: err
1629  integer :: ndims
1630  integer,dimension(nf90_max_var_dims) :: dimids
1631  integer :: i
1632  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1633 
1634 
1635  if (fileobj%is_root) then
1636  append_error_msg = "get_variable_dimension_names: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
1637 
1638  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
1639  err = nf90_inquire_variable(fileobj%ncid, varid, ndims=ndims, &
1640  dimids=dimids)
1641  call check_netcdf_code(err, append_error_msg)
1642  if (ndims .gt. 0) then
1643  if (size(dim_names) .ne. ndims) then
1644  call error("'names' has to be the same size of the number of dimensions for the variable. &
1645  &Check your get_variable_dimension_names call for file "//trim(fileobj%path)// &
1646  " and variable:"//trim(variable_name))
1647  endif
1648  else
1649  call error("get_variable_dimension_names: the variable: "//trim(variable_name)//" in file: "//trim(fileobj%path)&
1650  & //" does not any dimensions. ")
1651  endif
1652  dim_names(:) = ""
1653  do i = 1, ndims
1654  err = nf90_inquire_dimension(fileobj%ncid, dimids(i), name=dim_names(i))
1655  call check_netcdf_code(err, append_error_msg)
1656  enddo
1657  endif
1658  if (present(broadcast)) then
1659  if (.not. broadcast) then
1660  return
1661  endif
1662  endif
1663  call mpp_broadcast(ndims, fileobj%io_root, pelist=fileobj%pelist)
1664  if (.not. fileobj%is_root) then
1665  if (ndims .gt. 0) then
1666  if (size(dim_names) .ne. ndims) then
1667  call error("'names' has to be the same size of the number of dimensions for the variable. &
1668  & Check your get_variable_dimension_names call for file "//trim(fileobj%path)// &
1669  " and variable:"//trim(variable_name))
1670  endif
1671  else
1672  call error("get_variable_dimension_names: the variable: "//trim(variable_name)//" in file: "//trim(fileobj%path)&
1673  & //" does not any dimensions. ")
1674  endif
1675  dim_names(:) = ""
1676  endif
1677  call mpp_broadcast(dim_names, len(dim_names(ndims)), fileobj%io_root, &
1678  pelist=fileobj%pelist)
1679 end subroutine get_variable_dimension_names
1680 
1681 
1682 !> @brief Get the size of a variable's dimensions.
1683 subroutine get_variable_size(fileobj, variable_name, dim_sizes, broadcast)
1684 
1685  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1686  character(len=*), intent(in) :: variable_name !< Variable name.
1687  integer, dimension(:), intent(inout) :: dim_sizes !< Array of dimension
1688  !! sizes.
1689  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1690  !! not the data will be
1691  !! broadcasted to non
1692  !! "I/O root" ranks.
1693  !! The broadcast will be done
1694  !! by default.
1695 
1696  integer :: varid
1697  integer :: err
1698  integer :: ndims
1699  integer,dimension(nf90_max_var_dims) :: dimids
1700  integer :: i
1701  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1702 
1703  if (fileobj%is_root) then
1704  append_error_msg = "get_variable_size: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
1705  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
1706  err = nf90_inquire_variable(fileobj%ncid, varid, ndims=ndims, dimids=dimids)
1707  call check_netcdf_code(err, append_error_msg)
1708  if (ndims .gt. 0) then
1709  if (size(dim_sizes) .ne. ndims) then
1710  call error("'dim_sizes' has to be the same size of the number of dimensions for the variable. &
1711  &Check your get_variable_size call for file "//trim(fileobj%path)// &
1712  " and variable:"//trim(variable_name))
1713  endif
1714  else
1715  call error("get_variable_size: the variable: "//trim(variable_name)//" in file: "//trim(fileobj%path)//&
1716  &" does not any dimensions. ")
1717  endif
1718  do i = 1, ndims
1719  err = nf90_inquire_dimension(fileobj%ncid, dimids(i), len=dim_sizes(i))
1720  call check_netcdf_code(err, append_error_msg)
1721  enddo
1722  endif
1723  if (present(broadcast)) then
1724  if (.not. broadcast) then
1725  return
1726  endif
1727  endif
1728  call mpp_broadcast(ndims, fileobj%io_root, pelist=fileobj%pelist)
1729  if (.not. fileobj%is_root) then
1730  if (ndims .gt. 0) then
1731  if (size(dim_sizes) .ne. ndims) then
1732  call error("'dim_sizes' has to be the same size of the number of dimensions for the variable. &
1733  &Check your get_variable_size call for file "//trim(fileobj%path)// &
1734  " and variable:"//trim(variable_name))
1735  endif
1736  else
1737  call error("get_variable_size: the variable: "//trim(variable_name)//" in file: "//trim(fileobj%path)//&
1738  &" does not any dimensions. ")
1739  endif
1740  endif
1741  call mpp_broadcast(dim_sizes, ndims, fileobj%io_root, pelist=fileobj%pelist)
1742 end subroutine get_variable_size
1743 
1744 
1745 !> @brief Get the index of a variable's unlimited dimensions.
1746 !! @return The index of the unlimited dimension, or else
1747 !! no_unlimited_dimension.
1748 function get_variable_unlimited_dimension_index(fileobj, variable_name, &
1749  broadcast) &
1750  result(unlim_dim_index)
1751 
1752  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1753  character(len=*), intent(in) :: variable_name !< Variable name.
1754  logical, intent(in), optional :: broadcast !< Flag controlling whether or
1755  !! not the data will be
1756  !! broadcasted to non
1757  !! "I/O root" ranks.
1758  !! The broadcast will be done
1759  !! by default.
1760  integer :: unlim_dim_index
1761 
1762  integer :: ndims
1763  character(len=nf90_max_name), dimension(:), allocatable :: dim_names
1764  integer :: i
1765 
1766  unlim_dim_index = no_unlimited_dimension
1767  if (fileobj%is_root) then
1768  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
1769  allocate(dim_names(ndims))
1770  call get_variable_dimension_names(fileobj, variable_name, dim_names, &
1771  broadcast=.false.)
1772  do i = 1, size(dim_names)
1773  if (is_dimension_unlimited(fileobj, dim_names(i), .false.)) then
1774  unlim_dim_index = i
1775  exit
1776  endif
1777  enddo
1778  deallocate(dim_names)
1779  endif
1780  if (present(broadcast)) then
1781  if (.not. broadcast) then
1782  return
1783  endif
1784  endif
1785  call mpp_broadcast(unlim_dim_index, fileobj%io_root, pelist=fileobj%pelist)
1787 
1788 
1789 !> @brief Store the valid range for a variable.
1790 !! @return A ValidType_t object containing data about the valid
1791 !! range data for this variable can take.
1792 function get_valid(fileobj, variable_name) &
1793  result(valid)
1794 
1795  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
1796  character(len=*), intent(in) :: variable_name !< Variable name.
1797  type(valid_t) :: valid
1798 
1799  integer :: varid
1800  real(kind=r8_kind) :: scale_factor
1801  real(kind=r8_kind) :: add_offset
1802  real(kind=r8_kind), dimension(2) :: buffer
1803  integer :: xtype
1804  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
1805 
1806  append_error_msg = "get_valid: file:"//trim(fileobj%path)
1807  if (fileobj%is_root) then
1808  varid = get_variable_id(fileobj%ncid, variable_name, msg=append_error_msg)
1809  valid%has_max = .false.
1810  valid%has_min = .false.
1811  valid%has_fill = .false.
1812  valid%has_missing = .false.
1813  valid%has_range = .false.
1814 
1815  !This routine makes use of netcdf's automatic type conversion to
1816  !store all range information in double precision.
1817  if (attribute_exists(fileobj%ncid, varid, "scale_factor", msg=append_error_msg)) then
1818  call get_variable_attribute(fileobj, variable_name, "scale_factor", scale_factor, &
1819  broadcast=.false.)
1820  else
1821  scale_factor = 1._r8_kind
1822  endif
1823  if (attribute_exists(fileobj%ncid, varid, "add_offset", msg=append_error_msg)) then
1824  call get_variable_attribute(fileobj, variable_name, "add_offset", add_offset, &
1825  broadcast=.false.)
1826  else
1827  add_offset = 0._r8_kind
1828  endif
1829 
1830  !valid%max_val and valid%min_val are defined by the "valid_range", "valid_min", and
1831  !"valid_max" variable attributes if they are present in the file. If either the maximum value
1832  !or minimum value is defined, valid%has_range is set to .true. (i.e. open ended ranges
1833  !are valid and should be tested within the is_valid function).
1834  if (attribute_exists(fileobj%ncid, varid, "valid_range", msg=append_error_msg)) then
1835  call get_variable_attribute(fileobj, variable_name, "valid_range", buffer, &
1836  broadcast=.false.)
1837  valid%max_val = buffer(2)*scale_factor + add_offset
1838  valid%has_max = .true.
1839  valid%min_val = buffer(1)*scale_factor + add_offset
1840  valid%has_min = .true.
1841  else
1842  if (attribute_exists(fileobj%ncid, varid, "valid_max", msg=append_error_msg)) then
1843  call get_variable_attribute(fileobj, variable_name, "valid_max", buffer(1), &
1844  broadcast=.false.)
1845  valid%max_val = buffer(1)*scale_factor + add_offset
1846  valid%has_max = .true.
1847  endif
1848  if (attribute_exists(fileobj%ncid, varid, "valid_min", msg=append_error_msg)) then
1849  call get_variable_attribute(fileobj, variable_name, "valid_min", buffer(1), &
1850  broadcast=.false.)
1851  valid%min_val = buffer(1)*scale_factor + add_offset
1852  valid%has_min = .true.
1853  endif
1854  endif
1855  valid%has_range = valid%has_min .or. valid%has_max
1856 
1857 
1858  !Get the missing value from the file if it exists.
1859  if (attribute_exists(fileobj%ncid, varid, "missing_value", msg=append_error_msg)) then
1860  call get_variable_attribute(fileobj, variable_name, "missing_value", buffer(1), &
1861  broadcast=.false.)
1862  valid%missing_val = buffer(1)*scale_factor + add_offset
1863  valid%has_missing = .true.
1864  endif
1865 
1866  !Get the fill value from the file if it exists.
1867  !If the _FillValue attribute is present and the maximum or minimum value is not defined,
1868  !then the maximum or minimum value will be determined by the _FillValue according to the NUG convention.
1869  !The NUG convention states that a positive fill value will be the exclusive upper
1870  !bound (i.e. valid values are less than the fill value), while a
1871  !non-positive fill value will be the exclusive lower bound (i.e. valis
1872  !values are greater than the fill value). As before, valid%has_range is true
1873  !if either a maximum or minimum value is set.
1874  if (attribute_exists(fileobj%ncid, varid, "_FillValue", msg=append_error_msg)) then
1875  call get_variable_attribute(fileobj, variable_name, "_FillValue", buffer(1), &
1876  broadcast=.false.)
1877  valid%fill_val = buffer(1)*scale_factor + add_offset
1878  valid%has_fill = .true.
1879  xtype = get_variable_type(fileobj%ncid, varid, msg=append_error_msg)
1880 
1881  if (.not. valid%has_range) then
1882  if (xtype .eq. nf90_short .or. xtype .eq. nf90_int) then
1883  if (buffer(1) .gt. 0) then
1884  valid%max_val = (buffer(1) - 1._r8_kind)*scale_factor + add_offset
1885  valid%has_max = .true.
1886  else
1887  valid%min_val = (buffer(1) + 1._r8_kind)*scale_factor + add_offset
1888  valid%has_min = .true.
1889  endif
1890  elseif (xtype .eq. nf90_float .or. xtype .eq. nf90_double) then
1891  if (buffer(1) .gt. 0) then
1892  valid%max_val = (nearest(nearest(buffer(1), -1._r8_kind), -1._r8_kind)) &
1893  *scale_factor + add_offset
1894  valid%has_max = .true.
1895  else
1896  valid%min_val = (nearest(nearest(buffer(1), 1._r8_kind), 1._r8_kind)) &
1897  *scale_factor + add_offset
1898  valid%has_min = .true.
1899  endif
1900  else
1901  call error("Unsupported variable type:"//trim(append_error_msg))
1902  endif
1903  valid%has_range = .true.
1904  endif
1905  endif
1906 
1907  endif
1908 
1909  call mpp_broadcast(valid%has_min, fileobj%io_root, pelist=fileobj%pelist)
1910  if (valid%has_min) then
1911  call mpp_broadcast(valid%min_val, fileobj%io_root, pelist=fileobj%pelist)
1912  endif
1913  call mpp_broadcast(valid%has_max, fileobj%io_root, pelist=fileobj%pelist)
1914  if (valid%has_max) then
1915  call mpp_broadcast(valid%max_val, fileobj%io_root, pelist=fileobj%pelist)
1916  endif
1917  call mpp_broadcast(valid%has_range, fileobj%io_root, pelist=fileobj%pelist)
1918 
1919  call mpp_broadcast(valid%has_fill, fileobj%io_root, pelist=fileobj%pelist)
1920  if (valid%has_fill) then
1921  call mpp_broadcast(valid%fill_val, fileobj%io_root, pelist=fileobj%pelist)
1922  endif
1923 
1924  call mpp_broadcast(valid%has_missing, fileobj%io_root, pelist=fileobj%pelist)
1925  if (valid%has_missing) then
1926  call mpp_broadcast(valid%missing_val, fileobj%io_root, pelist=fileobj%pelist)
1927  endif
1928 
1929 end function get_valid
1930 
1931 !> @brief Determine if a piece of (r4) data is "valid" (in the correct range.)
1932 !! @return A flag telling if the data element is "valid."
1933 elemental function is_valid_r8(datum, validobj) &
1934  result(valid_data)
1935 
1936  real(kind=r8_kind), intent(in) :: datum !< Unpacked data element.
1937  type(valid_t), intent(in) :: validobj !< Valid object.
1938  logical :: valid_data
1939 
1940  valid_data = check_if_valid(datum, validobj)
1941 
1942 end function is_valid_r8
1943 
1944 !> @brief Determine if a piece of (r8) data is "valid" (in the correct range.)
1945 !! @return A flag telling if the data element is "valid."
1946 elemental function is_valid_r4(datum, validobj) &
1947  result(valid_data)
1948 
1949  real(kind=r4_kind), intent(in) :: datum !< Unpacked data element.
1950  type(valid_t), intent(in) :: validobj !< Valid object.
1951  logical :: valid_data
1952 
1953  real(kind=r8_kind) :: rdatum
1954 
1955  !< Convert the data to r8 so it can be compared correctly since validobj%*_val are defined as r8
1956  rdatum = real(datum, kind=r8_kind)
1957  valid_data = check_if_valid(rdatum, validobj)
1958 
1959 end function is_valid_r4
1960 
1961 !> @brief Determine if a piece of data is "valid" (in the correct range.)
1962 !! @return A flag telling if the data element is "valid."
1963 elemental function check_if_valid(rdatum, validobj) &
1964  result(valid_data)
1965 
1966  real(kind=r8_kind), intent(in) :: rdatum !< packed data element.
1967  type(valid_t), intent(in) :: validobj !< Valid object.
1968  logical :: valid_data
1969 
1970  valid_data = .true.
1971  ! If the variable has a range (open or closed), valid values must be in that
1972  ! range.
1973  if (validobj%has_range) then
1974  if (validobj%has_min .and. .not. validobj%has_max) then
1975  valid_data = rdatum .ge. validobj%min_val
1976  elseif (validobj%has_max .and. .not. validobj%has_min) then
1977  valid_data = rdatum .le. validobj%max_val
1978  else
1979  valid_data = .not. (rdatum .lt. validobj%min_val .or. rdatum .gt. validobj%max_val)
1980  endif
1981  endif
1982  ! If the variable has a fill value or missing value, valid values must not be
1983  ! equal to either.
1984  if (validobj%has_fill .or. validobj%has_missing) then
1985  if (validobj%has_fill .and. .not. validobj%has_missing) then
1986  valid_data = rdatum .ne. validobj%fill_val
1987  elseif (validobj%has_missing .and. .not. validobj%has_fill) then
1988  valid_data = rdatum .ne. validobj%missing_val
1989  else
1990  valid_data = .not. (rdatum .eq. validobj%missing_val .or. rdatum .eq. validobj%fill_val)
1991  endif
1992  endif
1993 end function check_if_valid
1994 
1995 
1996 !> @brief Gathers a compressed arrays size and offset for each pe.
1997 subroutine compressed_start_and_count(fileobj, nelems, npes_start, npes_count)
1998 
1999  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2000  integer, intent(in) :: nelems !< Number of elements on the current pe.
2001  integer, dimension(:), allocatable, intent(out) :: npes_start !< Offset for each pe.
2002  integer, dimension(:), allocatable, intent(out) :: npes_count !< Number of elements for
2003  !! each pe.
2004 
2005  integer :: i
2006 
2007  allocate(npes_start(size(fileobj%pelist)))
2008  allocate(npes_count(size(fileobj%pelist)))
2009  do i = 1, size(fileobj%pelist)
2010  if (fileobj%pelist(i) .eq. mpp_pe()) then
2011  npes_count(i) = nelems
2012  else
2013  call mpp_recv(npes_count(i), fileobj%pelist(i), block=.false.)
2014  call mpp_send(nelems, fileobj%pelist(i))
2015  endif
2016  enddo
2017  call mpp_sync_self(check=event_recv)
2018  call mpp_sync_self(check=event_send)
2019  npes_start(1) = 1
2020  do i = 1, size(fileobj%pelist)-1
2021  npes_start(i+1) = npes_start(i) + npes_count(i)
2022  enddo
2023 end subroutine compressed_start_and_count
2024 
2025 
2026 include "netcdf_add_restart_variable.inc"
2027 include "netcdf_read_data.inc"
2028 include "netcdf_write_data.inc"
2029 include "register_global_attribute.inc"
2030 include "register_variable_attribute.inc"
2031 include "get_global_attribute.inc"
2032 include "get_variable_attribute.inc"
2033 include "compressed_write.inc"
2034 include "compressed_read.inc"
2035 include "scatter_data_bc.inc"
2036 include "gather_data_bc.inc"
2037 include "unpack_data.inc"
2038 
2039 !> @brief Wrapper to distinguish interfaces.
2040 function netcdf_file_open_wrap(fileobj, path, mode, nc_format, pelist, is_restart, dont_add_res_to_filename) &
2041  result(success)
2042 
2043  type(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
2044  character(len=*), intent(in) :: path !< File path.
2045  character(len=*), intent(in) :: mode !< File mode. Allowed values are:
2046  !! "read", "append", "write", or
2047  !! "overwrite".
2048  character(len=*), intent(in), optional :: nc_format !< Netcdf format that
2049  !! new files are written
2050  !! as. Allowed values
2051  !! are: "64bit", "classic",
2052  !! or "netcdf4". Defaults to
2053  !! "64bit".
2054  integer, dimension(:), intent(in), optional :: pelist !< List of ranks associated
2055  !! with this file. If not
2056  !! provided, only the current
2057  !! rank will be able to
2058  !! act on the file.
2059  logical, intent(in), optional :: is_restart !< Flag telling if this file
2060  !! is a restart file. Defaults
2061  !! to false.
2062  logical, intent(in), optional :: dont_add_res_to_filename !< Flag indicating not to add
2063  !! ".res" to the filename
2064  logical :: success
2065 
2066  success = netcdf_file_open(fileobj, path, mode, nc_format, pelist, is_restart, dont_add_res_to_filename)
2067 end function netcdf_file_open_wrap
2068 
2069 
2070 !> @brief Wrapper to distinguish interfaces.
2071 subroutine netcdf_file_close_wrap(fileobj)
2072 
2073  type(fmsnetcdffile_t), intent(inout) :: fileobj !< File object.
2074 
2075  call netcdf_file_close(fileobj)
2076 end subroutine netcdf_file_close_wrap
2077 
2078 
2079 !> @brief Wrapper to distinguish interfaces.
2080 subroutine netcdf_add_variable_wrap(fileobj, variable_name, variable_type, dimensions)
2081 
2082  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2083  character(len=*), intent(in) :: variable_name !< Variable name.
2084  character(len=*), intent(in) :: variable_type !< Variable type. Allowed
2085  !! values are: "int", "int64",
2086  !! "float", or "double".
2087  character(len=*), dimension(:), intent(in), optional :: dimensions !< Dimension names.
2088 
2089  call netcdf_add_variable(fileobj, variable_name, variable_type, dimensions)
2090 end subroutine netcdf_add_variable_wrap
2091 
2092 !> @brief Wrapper to distinguish interfaces.
2093 subroutine netcdf_save_restart_wrap(fileobj, unlim_dim_level)
2094 
2095  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2096  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
2097  !! level.
2098 
2099  call netcdf_save_restart(fileobj, unlim_dim_level)
2100 end subroutine netcdf_save_restart_wrap
2101 
2102 
2103 !> @brief Returns a variable's fill value if it exists in the file.
2104 !! @return Flag telling if a fill value exists.
2105 function get_fill_value(fileobj, variable_name, fill_value, broadcast) &
2106  result(fill_exists)
2107 
2108  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2109  character(len=*), intent(in) :: variable_name !< Variable name.
2110  class(*), intent(out) :: fill_value !< Fill value.
2111  logical, intent(in), optional :: broadcast !< Flag controlling whether or
2112  !! not the data will be
2113  !! broadcasted to non
2114  !! "I/O root" ranks.
2115  !! The broadcast will be done
2116  !! by default.
2117  logical :: fill_exists
2118 
2119  character(len=32), dimension(2) :: attribute_names
2120  logical :: bcast
2121  integer :: i
2122 
2123  fill_exists = .false.
2124  call string_copy(attribute_names(1), "_FillValue")
2125  call string_copy(attribute_names(2), "missing_value")
2126  if (present(broadcast)) then
2127  bcast = broadcast
2128  else
2129  bcast = .true.
2130  endif
2131  do i = 1, size(attribute_names)
2132  fill_exists = variable_att_exists(fileobj, variable_name, attribute_names(i), &
2133  broadcast=bcast)
2134  if (fill_exists) then
2135  call get_variable_attribute(fileobj, variable_name, attribute_names(i), &
2136  fill_value, broadcast=bcast)
2137  exit
2138  endif
2139  enddo
2140 end function get_fill_value
2141 
2142 
2143 function get_variable_sense(fileobj, variable_name) &
2144  result(variable_sense)
2145 
2146  class(fmsnetcdffile_t), intent(in) :: fileobj
2147  character(len=*), intent(in) :: variable_name
2148  integer :: variable_sense
2149 
2150  character(len=256) :: buf
2151 
2152  variable_sense = 0
2153  if (variable_att_exists(fileobj, variable_name, "positive")) then
2154  call get_variable_attribute(fileobj, variable_name, "positive", buf)
2155  if (string_compare(buf, "down")) then
2156  variable_sense = -1
2157  elseif (string_compare(buf, "up")) then
2158  variable_sense = 1
2159  endif
2160  endif
2161 end function get_variable_sense
2162 
2163 
2164 function get_variable_missing(fileobj, variable_name) &
2165  result(variable_missing)
2166 
2167  type(fmsnetcdffile_t), intent(in) :: fileobj
2168  character(len=*), intent(in) :: variable_name
2169  real(kind=r8_kind) :: variable_missing
2170 
2171  real(kind=r8_kind) :: variable_missing_1d(1) !< Workaround for pgi
2172 
2173  if (variable_att_exists(fileobj, variable_name, "_FillValue")) then
2174  call get_variable_attribute(fileobj, variable_name, "_FillValue", variable_missing_1d)
2175  elseif (variable_att_exists(fileobj, variable_name, "missing_value")) then
2176  call get_variable_attribute(fileobj, variable_name, "missing_value", variable_missing_1d)
2177  elseif (variable_att_exists(fileobj, variable_name, "missing")) then
2178  call get_variable_attribute(fileobj, variable_name, "missing", variable_missing_1d)
2179  else
2180  variable_missing_1d = mpp_fill_double
2181  endif
2182 
2183  variable_missing = variable_missing_1d(1)
2184 
2185 end function get_variable_missing
2186 
2187 
2188 subroutine get_variable_units(fileobj, variable_name, units)
2189 
2190  class(fmsnetcdffile_t), intent(in) :: fileobj
2191  character(len=*), intent(in) :: variable_name
2192  character(len=*), intent(out) :: units
2193 
2194  if (variable_att_exists(fileobj, variable_name, "units")) then
2195  call get_variable_attribute(fileobj, variable_name, "units", units)
2196  else
2197  units = "nounits"
2198  endif
2199 end subroutine get_variable_units
2200 
2201 
2202 subroutine get_time_calendar(fileobj, time_name, calendar_type)
2203 
2204  class(fmsnetcdffile_t), intent(in) :: fileobj
2205  character(len=*), intent(in) :: time_name
2206  character(len=*), intent(out) :: calendar_type
2207 
2208  if (variable_att_exists(fileobj, time_name, "calendar")) then
2209  call get_variable_attribute(fileobj, time_name, "calendar", calendar_type)
2210  elseif (variable_att_exists(fileobj, time_name, "calendar_type")) then
2211  call get_variable_attribute(fileobj, time_name, "calendar_type", calendar_type)
2212  else
2213  calendar_type = "unspecified"
2214  endif
2215 end subroutine get_time_calendar
2216 
2217 
2218 !> @brief Determine if a variable has been registered to a restart file..
2219 !! @return Flag telling if the variable has been registered to a restart file.
2220 function is_registered_to_restart(fileobj, variable_name) &
2221  result(is_registered)
2222 
2223  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2224  character(len=*), intent(in) :: variable_name !< Variable name.
2225  logical :: is_registered
2226 
2227  integer :: i
2228 
2229  if (.not. fileobj%is_restart) then
2230  call error("file "//trim(fileobj%path)//" is not a restart file. &
2231  &Add is_restart=.true. to your open_file call")
2232  endif
2233  is_registered = .false.
2234  do i = 1, fileobj%num_restart_vars
2235  if (string_compare(fileobj%restart_vars(i)%varname, variable_name, .true.)) then
2236  is_registered = .true.
2237  exit
2238  endif
2239  enddo
2240 end function is_registered_to_restart
2241 
2242 
2243 function check_if_open(fileobj, fname) result(is_open)
2244  logical :: is_open !< True if the file in the file object is opened
2245  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
2246  character(len=*), intent(in), optional :: fname !< Optional filename for checking
2247 
2248  !Check if the is_open variable in the object has been allocated
2249  if (allocated(fileobj%is_open)) then
2250  is_open = fileobj%is_open !Return the value of the fileobj%is_open
2251  else
2252  is_open = .false. !If fileobj%is_open is not allocated, that the file has not been opened
2253  endif
2254 
2255  if (present(fname)) then
2256  !If the filename does not match the name in path,
2257  !then this is considered not open
2258  if (is_open .AND. trim(fname) .ne. trim(fileobj%path)) is_open = .false.
2259  endif
2260 end function check_if_open
2261 
2262 subroutine set_fileobj_time_name (fileobj,time_name)
2263  class(fmsnetcdffile_t), intent(inout) :: fileobj
2264  character(*),intent(in) :: time_name
2265  integer :: len_of_name
2266  len_of_name = len(trim(time_name))
2267  fileobj%time_name = ' '
2268  fileobj%time_name = time_name(1:len_of_name)
2269 ! if (.not. allocated(fileobj%time_name)) then
2270 ! allocate(character(len=len_of_name) :: fileobj%time_name)
2271 ! fileobj%time_name = time_name(1:len_of_name)
2272 ! else
2273 ! call error ("set_fileobj_time_name :: The time_name has already been set")
2274 ! endif
2275 end subroutine set_fileobj_time_name
2276 
2277 !> @brief Loop through the registered restart variables (including regional
2278 !! variables) and read them from the netcdf file
2279 subroutine read_restart_bc(fileobj, unlim_dim_level, ignore_checksum)
2280  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object
2281  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
2282  !! level.
2283  logical, intent(in), optional :: ignore_checksum !< Checksum data integrity flag.
2284 
2285  integer :: i !< No description
2286 
2287  if (.not. fileobj%is_restart) then
2288  call error("file "//trim(fileobj%path)//" is not a restart file.")
2289  endif
2290 
2291  do i = 1, fileobj%num_restart_vars
2292  !> Go away if you are not in the pelist!
2293  if (.not.any(mpp_pe().eq.fileobj%restart_vars(i)%bc_info%pelist(:))) cycle
2294 
2295  !> The file's root pe reads the file and scatters it to the rest of the pes
2296  if (associated(fileobj%restart_vars(i)%data2d)) then
2297  call scatter_data_bc (fileobj, fileobj%restart_vars(i)%varname, &
2298  fileobj%restart_vars(i)%data2d, &
2299  fileobj%restart_vars(i)%bc_info, &
2300  unlim_dim_level = unlim_dim_level, &
2301  ignore_checksum=ignore_checksum)
2302  else if (associated(fileobj%restart_vars(i)%data3d)) then
2303  call scatter_data_bc (fileobj, fileobj%restart_vars(i)%varname, &
2304  fileobj%restart_vars(i)%data3d, &
2305  fileobj%restart_vars(i)%bc_info, &
2306  unlim_dim_level = unlim_dim_level, &
2307  ignore_checksum=ignore_checksum)
2308  endif
2309  end do
2310 
2311 
2312 end subroutine read_restart_bc
2313 
2314 !> @brief Loop through the registered restart variables (including regional
2315 !! variables) and write them to the netcdf file
2316 subroutine write_restart_bc(fileobj, unlim_dim_level)
2317  class(fmsnetcdffile_t), intent(inout) :: fileobj !< File object
2318  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
2319  !! level.
2320  integer :: i !< No description
2321 
2322  if (.not. fileobj%is_restart) then
2323  call error("file "//trim(fileobj%path)//" is not a restart file. &
2324  &Add is_restart=.true. to your open_file call")
2325  endif
2326 
2327  !> Loop through the variables, root pe gathers the data from the other pes and writes out the checksum.
2328  !! Then loop through the variables again to write the data to the netcdf file.
2329  !! All the metadata should be written before the data to prevent netcdf from rewriting the file
2330  !! if the header is not big enough. That is why the two do loops are needed.
2331 
2332  do i = 1, fileobj%num_restart_vars
2333  !> Go away if you are not in the pelist!
2334  if (.not.any(mpp_pe().eq.fileobj%restart_vars(i)%bc_info%pelist(:))) cycle
2335 
2336  !> Go away if this is not a BC variable
2337  if (.not. fileobj%restart_vars(i)%is_bc_variable) cycle
2338 
2339  !> Root pe gathers the data from the other ranks, saves it in a buffer, and writes out the checksum.
2340  if (associated(fileobj%restart_vars(i)%data2d)) then
2341  call gather_data_bc(fileobj, fileobj%restart_vars(i)%data2d, fileobj%restart_vars(i)%bc_info)
2342  call register_variable_attribute(fileobj, fileobj%restart_vars(i)%varname, "checksum", &
2343  fileobj%restart_vars(i)%bc_info%chksum(1:len(fileobj%restart_vars(i)%bc_info%chksum)),&
2344  str_len=len(fileobj%restart_vars(i)%bc_info%chksum))
2345  else if (associated(fileobj%restart_vars(i)%data3d)) then
2346  call gather_data_bc(fileobj, fileobj%restart_vars(i)%data3d, fileobj%restart_vars(i)%bc_info)
2347  call register_variable_attribute(fileobj, fileobj%restart_vars(i)%varname, "checksum", &
2348  fileobj%restart_vars(i)%bc_info%chksum(1:len(fileobj%restart_vars(i)%bc_info%chksum)),&
2349  str_len=len(fileobj%restart_vars(i)%bc_info%chksum))
2350  endif
2351  enddo
2352 
2353  !> Write the data to the netcdf file
2354  do i = 1, fileobj%num_restart_vars
2355  if (allocated(fileobj%restart_vars(i)%bc_info%globaldata2d_r8 )) then
2356  call netcdf_write_data(fileobj, fileobj%restart_vars(i)%varname, &
2357  fileobj%restart_vars(i)%bc_info%globaldata2d_r8 , &
2358  unlim_dim_level=unlim_dim_level)
2359  deallocate(fileobj%restart_vars(i)%bc_info%globaldata2d_r8)
2360  else if (allocated(fileobj%restart_vars(i)%bc_info%globaldata2d_r4 )) then
2361  call netcdf_write_data(fileobj, fileobj%restart_vars(i)%varname, &
2362  fileobj%restart_vars(i)%bc_info%globaldata2d_r4 , &
2363  unlim_dim_level=unlim_dim_level)
2364  deallocate(fileobj%restart_vars(i)%bc_info%globaldata2d_r4)
2365  else if (allocated(fileobj%restart_vars(i)%bc_info%globaldata3d_r8 )) then
2366  call netcdf_write_data(fileobj, fileobj%restart_vars(i)%varname, &
2367  fileobj%restart_vars(i)%bc_info%globaldata3d_r8 , &
2368  unlim_dim_level=unlim_dim_level)
2369  deallocate(fileobj%restart_vars(i)%bc_info%globaldata3d_r8)
2370  else if (allocated(fileobj%restart_vars(i)%bc_info%globaldata3d_r4 )) then
2371  call netcdf_write_data(fileobj, fileobj%restart_vars(i)%varname, &
2372  fileobj%restart_vars(i)%bc_info%globaldata3d_r4 , &
2373  unlim_dim_level=unlim_dim_level)
2374  deallocate(fileobj%restart_vars(i)%bc_info%globaldata3d_r4 )
2375  endif
2376 
2377  enddo
2378 
2379 end subroutine write_restart_bc
2380 
2381 !> @brief flushes the netcdf file into disk
2382 subroutine flush_file(fileobj)
2383  class(fmsnetcdffile_t), intent(inout) :: fileobj !< FMS2_io fileobj
2384 
2385  integer :: err !< Netcdf error code
2386 
2387  if (fileobj%is_root) then
2388  err = nf90_sync(fileobj%ncid)
2389  call check_netcdf_code(err, "Flush_file: File:"//trim(fileobj%path))
2390  endif
2391 end subroutine flush_file
2392 
2393 end module netcdf_io_mod
2394 !> @}
2395 ! close documentation grouping
subroutine compressed_write_1d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine register_variable_attribute_1d(fileobj, variable_name, attribute_name, attribute_value, str_len)
Add an attribute to a variable.
subroutine compressed_write_3d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_5d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_0d(fileobj, variable_name, cdata, unlim_dim_level, corner)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_read_0d(fileobj, variable_name, cdata, unlim_dim_level, corner)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_write_5d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_read_5d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_read_3d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_write_0d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner)
Wrapper to distinguish interfaces.
subroutine compressed_read_1d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_read_2d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_write_4d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_2d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_1d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_4d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_read_4d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks....
subroutine compressed_write_2d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine register_variable_attribute_0d(fileobj, variable_name, attribute_name, attribute_value, str_len)
Add an attribute to a variable.
subroutine compressed_write_3d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine, public get_instance_filename(name_in, name_out)
Adds the filename_appendix to name_in and sets it as name_out.
subroutine, public restart_filepath_mangle(dest, source)
Add ".res" to an input file path.
logical function, public file_exists(path)
Determine if a file exists.
logical function, public string_compare(string1, string2, ignore_case)
Compare strings.
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
integer, parameter, public mpp_comm_null
MPP_COMM_NULL acts as an analagous mpp-macro for MPI_COMM_NULL to share with fms2_io NetCDF4 mpi-io....
Definition: mpp.F90:1345
integer, parameter, public mpp_info_null
MPP_INFO_NULL acts as an analagous mpp-macro for MPI_INFO_NULL to share with fms2_io NetCDF4 mpi-io....
Definition: mpp.F90:1336
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Perform parallel broadcasts.
Definition: mpp.F90:1091
Error handler.
Definition: mpp.F90:382
Gather data sent from pelist onto the root pe Wrapper for MPI_gather, can be used with and without in...
Definition: mpp.F90:698
Recieve data from another PE.
Definition: mpp.F90:937
Send data to a receiving PE.
Definition: mpp.F90:1004
integer, private fms2_ncchksz
Chunksize (bytes) used in nc_open and nc_create.
Definition: netcdf_io.F90:54
subroutine, public netcdf_restore_state(fileobj, unlim_dim_level)
Loop through registered restart variables and read them from a netcdf file.
Definition: netcdf_io.F90:1127
subroutine append_compressed_dimension(fileobj, dim_name, npes_corner, npes_nelems)
Add a compressed dimension to a file object.
Definition: netcdf_io.F90:792
subroutine, public netcdf_file_close(fileobj)
Close a netcdf file.
Definition: netcdf_io.F90:729
subroutine netcdf_add_restart_variable_4d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
integer function, dimension(2) get_variable_compressed_dimension_index(fileobj, variable_name, broadcast)
Given a compressed variable, get the index of the compressed dimension.
Definition: netcdf_io.F90:1012
logical function, public netcdf_file_open(fileobj, path, mode, nc_format, pelist, is_restart, dont_add_res_to_filename)
Open a netcdf file.
Definition: netcdf_io.F90:542
elemental logical function is_valid_r8(datum, validobj)
Determine if a piece of (r4) data is "valid" (in the correct range.)
Definition: netcdf_io.F90:1935
logical function attribute_exists(ncid, varid, attribute_name, msg)
Determine if an attribute exists.
Definition: netcdf_io.F90:479
subroutine, public get_variable_size(fileobj, variable_name, dim_sizes, broadcast)
Get the size of a variable's dimensions.
Definition: netcdf_io.F90:1684
subroutine, public netcdf_add_dimension(fileobj, dimension_name, dimension_length, is_compressed)
Add a dimension to a file.
Definition: netcdf_io.F90:862
subroutine netcdf_add_restart_variable_0d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
type(valid_t) function, public get_valid(fileobj, variable_name)
Store the valid range for a variable.
Definition: netcdf_io.F90:1794
subroutine, public read_restart_bc(fileobj, unlim_dim_level, ignore_checksum)
Loop through the registered restart variables (including regional variables) and read them from the n...
Definition: netcdf_io.F90:2280
logical function, public is_dimension_unlimited(fileobj, dimension_name, broadcast)
Determine where or not the dimension is unlimited.
Definition: netcdf_io.F90:1360
character(len=10), private fms2_nc_format
Netcdf format type used in netcdf_file_open.
Definition: netcdf_io.F90:56
integer function, public get_variable_unlimited_dimension_index(fileobj, variable_name, broadcast)
Get the index of a variable's unlimited dimensions.
Definition: netcdf_io.F90:1751
logical function, public netcdf_file_open_wrap(fileobj, path, mode, nc_format, pelist, is_restart, dont_add_res_to_filename)
Wrapper to distinguish interfaces.
Definition: netcdf_io.F90:2042
integer function get_compressed_dimension_index(fileobj, dim_name)
Get the index of a compressed dimension in a file object.
Definition: netcdf_io.F90:771
subroutine, public compressed_start_and_count(fileobj, nelems, npes_start, npes_count)
Gathers a compressed arrays size and offset for each pe.
Definition: netcdf_io.F90:1998
subroutine scatter_data_bc_2d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
integer, private fms2_header_buffer_val
value used in NF__ENDDEF
Definition: netcdf_io.F90:57
integer, private fms2_nc_format_param
Netcdf format type param used in nc_create.
Definition: netcdf_io.F90:55
elemental logical function is_valid_r4(datum, validobj)
Determine if a piece of (r8) data is "valid" (in the correct range.)
Definition: netcdf_io.F90:1948
integer function, public get_variable_num_dimensions(fileobj, variable_name, broadcast)
Get the number of dimensions a variable depends on.
Definition: netcdf_io.F90:1580
subroutine, public get_dimension_size(fileobj, dimension_name, dim_size, broadcast)
Get the length of a dimension.
Definition: netcdf_io.F90:1430
subroutine, public netcdf_add_variable_wrap(fileobj, variable_name, variable_type, dimensions)
Wrapper to distinguish interfaces.
Definition: netcdf_io.F90:2081
subroutine netcdf_add_restart_variable_3d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
subroutine netcdf_write_data_4d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
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 netcdf_add_restart_variable_5d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
elemental logical function check_if_valid(rdatum, validobj)
Determine if a piece of data is "valid" (in the correct range.)
Definition: netcdf_io.F90:1965
subroutine netcdf_add_restart_variable_1d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
logical function, public global_att_exists(fileobj, attribute_name, broadcast)
Determine if a global attribute exists.
Definition: netcdf_io.F90:1175
subroutine, public set_netcdf_mode(ncid, mode)
Switch to the correct netcdf mode.
Definition: netcdf_io.F90:394
subroutine netcdf_write_data_0d(fileobj, variable_name, variable_data, unlim_dim_level, corner)
Write data to a variable in a netcdf file.
subroutine register_global_attribute_0d(fileobj, attribute_name, attribute_value, str_len)
Add a global attribute.
subroutine get_variable_attribute_1d(fileobj, variable_name, attribute_name, attribute_value, broadcast)
Get the value of a variable's attribute.
subroutine netcdf_add_restart_variable_0d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
logical, private fms2_shuffle
Flag indicating whether to use the netcdf shuffle filter.
Definition: netcdf_io.F90:60
subroutine get_global_attribute_1d(fileobj, attribute_name, attribute_value, broadcast)
Get the value of a global attribute.
subroutine netcdf_add_restart_variable_3d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
subroutine netcdf_add_restart_variable_1d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
subroutine get_variable_attribute_0d(fileobj, variable_name, attribute_name, attribute_value, broadcast, reproduce_null_char_bug_flag)
Get the value of a variable's attribute.
integer function get_variable_id(ncid, variable_name, msg, allow_failure)
Get the id of a variable from its name.
Definition: netcdf_io.F90:451
subroutine register_restart_region_3d(fileobj, variable_name, vdata, indices, global_size, pelist, is_root_pe, x_halo, y_halo, jshift, ishift, is_optional)
Registers a regional 3D variable and stores the information needed.
subroutine register_global_attribute_1d(fileobj, attribute_name, attribute_value, str_len)
Add a global attribute.
logical function, public dimension_exists(fileobj, dimension_name, broadcast)
Determine if a dimension exists.
Definition: netcdf_io.F90:1324
logical function, public is_registered_to_restart(fileobj, variable_name)
Determine if a variable has been registered to a restart file..
Definition: netcdf_io.F90:2222
integer function get_dimension_id(ncid, dimension_name, msg, allow_failure)
Get the id of a dimension from its name.
Definition: netcdf_io.F90:423
subroutine, public register_compressed_dimension(fileobj, dimension_name, npes_corner, npes_nelems)
Add a compressed dimension.
Definition: netcdf_io.F90:915
logical function, public variable_att_exists(fileobj, variable_name, attribute_name, broadcast)
Determine if a variable's attribute exists.
Definition: netcdf_io.F90:1204
integer, private fms2_deflate_level
Netcdf deflate level to use in nf90_def_var (integer between 1 to 9)
Definition: netcdf_io.F90:58
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_add_restart_variable_2d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
subroutine netcdf_write_data_3d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
logical function, public get_fill_value(fileobj, variable_name, fill_value, broadcast)
Returns a variable's fill value if it exists in the file.
Definition: netcdf_io.F90:2107
subroutine gather_data_bc_3d(fileobj, vdata, bc_info)
gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a buffer.
subroutine, public get_variable_dimension_names(fileobj, variable_name, dim_names, broadcast)
Get the name of a variable's dimensions.
Definition: netcdf_io.F90:1614
subroutine, public check_netcdf_code(err, msg)
Check for errors returned by netcdf.
Definition: netcdf_io.F90:378
subroutine netcdf_write_data_1d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine, public flush_file(fileobj)
flushes the netcdf file into disk
Definition: netcdf_io.F90:2383
subroutine netcdf_add_restart_variable_4d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
subroutine get_global_attribute_0d(fileobj, attribute_name, attribute_value, broadcast)
Get the value of a global attribute.
subroutine netcdf_add_restart_variable_5d(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Add a restart variable to a netcdf file.
subroutine, public register_unlimited_compressed_axis(fileobj, dimension_name, dimension_length)
Add a "compressed" unlimited dimension to a netcdf file.
Definition: netcdf_io.F90:828
subroutine, public netcdf_save_restart(fileobj, unlim_dim_level)
Loop through registered restart variables and write them to a netcdf file.
Definition: netcdf_io.F90:1085
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 netcdf_write_data_2d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine, public write_restart_bc(fileobj, unlim_dim_level)
Loop through the registered restart variables (including regional variables) and write them to the ne...
Definition: netcdf_io.F90:2317
integer function get_variable_type(ncid, varid, msg)
Get the type of a netcdf variable.
Definition: netcdf_io.F90:524
logical, private fms2_is_netcdf4
Flag indicating whether the default netcdf file format is netcdf4.
Definition: netcdf_io.F90:61
subroutine, public netcdf_save_restart_wrap(fileobj, unlim_dim_level)
Wrapper to distinguish interfaces.
Definition: netcdf_io.F90:2094
subroutine, public get_variable_names(fileobj, names, broadcast)
Get the names of the variables in a file.
Definition: netcdf_io.F90:1491
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 netcdf_add_restart_variable_2d_wrap(fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
Wrapper to distinguish interfaces.
subroutine, public get_dimension_names(fileobj, names, broadcast)
Get the names of the dimensions in a file.
Definition: netcdf_io.F90:1267
subroutine gather_data_bc_2d(fileobj, vdata, bc_info)
gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a buffer.
logical function, public check_if_open(fileobj, fname)
Definition: netcdf_io.F90:2244
subroutine, public netcdf_io_init(chksz, header_buffer_val, netcdf_default_format, deflate_level, shuffle)
Accepts the namelist fms2_io_nml variables relevant to netcdf_io_mod.
Definition: netcdf_io.F90:346
integer function, public get_num_variables(fileobj, broadcast)
Determine the number of variables in a file.
Definition: netcdf_io.F90:1464
subroutine, public get_unlimited_dimension_name(fileobj, dimension_name, broadcast)
Get the name of the unlimited dimension.
Definition: netcdf_io.F90:1395
integer function get_attribute_type(ncid, varid, attname, msg)
Get the type of a netcdf attribute.
Definition: netcdf_io.F90:504
subroutine, public netcdf_add_variable(fileobj, variable_name, variable_type, dimensions, chunksizes)
Add a variable to a file.
Definition: netcdf_io.F90:943
subroutine scatter_data_bc_3d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
subroutine netcdf_write_data_5d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine register_restart_region_2d(fileobj, variable_name, vdata, indices, global_size, pelist, is_root_pe, x_halo, y_halo, jshift, ishift, is_optional)
Registers a regional 2D variable and stores the information needed.
subroutine, public netcdf_file_close_wrap(fileobj)
Wrapper to distinguish interfaces.
Definition: netcdf_io.F90:2072
subroutine netcdf_read_data_0d(fileobj, variable_name, buf, unlim_dim_level, corner, broadcast)
Read in data from a variable in a netcdf file.
logical function, public variable_exists(fileobj, variable_name, broadcast)
Determine if a variable exists.
Definition: netcdf_io.F90:1548
subroutine add_restart_var_to_array(fileobj, variable_name)
Add a restart variable to a FmsNetcdfFile_t type.
Definition: netcdf_io.F90:1058
integer function, public get_num_dimensions(fileobj, broadcast)
Determine the number of dimensions in a file.
Definition: netcdf_io.F90:1240
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.
The interface is needed to accomodate pgi because it can't handle class * and there was no other way ...
Definition: netcdf_io.F90:334
information needed fr regional restart variables
Definition: netcdf_io.F90:67
information about the current dimensions for regional restart variables
Definition: netcdf_io.F90:116
Range type for a netcdf variable.
Definition: netcdf_io.F90:163