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