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