FMS  2024.03
Flexible Modeling System
fms2_io_mod

An updated library for parallel IO to replace mpp_io_mod. This module contains the public API for fms2 I/O interfaces and routines defined throughout this directory. A conversion guide for replacing mpp/fms_io code with this module is available below. More...

Data Types

interface  close_file
 Close a netcdf or domain file opened with open_file or open_virtual_file. More...
 
interface  open_file
 Opens a given netcdf or domain file. More...
 
interface  open_virtual_file
 Creates a diskless netcdf or domain file. More...
 
interface  read_data
 Read data from a defined field in a file. More...
 
interface  read_new_restart
 Read registered restarts from a new file Optionally takes directory to write to, model time and filename
Example usage: call read_new_restart(fileobj, unlimted_dimension_level) More...
 
interface  read_restart
 Reads in restart variables from a given file
Example usage: call read_restart(fileobj) Reads registered restart variables from fileobj. More...
 
interface  register_axis
 Add a dimension to a given file. More...
 
interface  register_field
 Defines a new field within the given file
Example usage: More...
 
interface  register_restart_field
 Similar to register_field, but occupies the field with data for restarts
Example usage: More...
 
interface  write_data
 Write data to a defined field within a file
Example usage: More...
 
interface  write_new_restart
 Writes all restart fields in a given restart file to a new restart file
Example usage: More...
 
interface  write_restart
 Writes all restart fields registered within a given restart file
Example usage: More...
 

Functions/Subroutines

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. This routine may only be used with variables that are "compressed". More...
 
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. This routine may only be used with variables that are "compressed". More...
 
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. This routine may only be used with variables that are "compressed". More...
 
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. This routine may only be used with variables that are "compressed". More...
 
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. This routine may only be used with variables that are "compressed". More...
 
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. This routine may only be used with variables that are "compressed". More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_0d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner)
 Wrapper to distinguish interfaces. More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_1d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_2d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_3d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_4d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
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. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file. More...
 
subroutine compressed_write_5d_wrap (fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
subroutine domain_read_0d (fileobj, variable_name, vdata, unlim_dim_level, corner)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_read_1d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_read_2d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_read_3d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_read_4d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_read_5d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_0d (fileobj, variable_name, vdata, unlim_dim_level, corner)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_1d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_2d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_3d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_4d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine domain_write_5d (fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
 Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed". More...
 
subroutine, public fms2_io_init ()
 Reads the fms2_io_nml. More...
 
subroutine register_domain_restart_variable_0d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_domain_restart_variable_1d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_domain_restart_variable_2d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_domain_restart_variable_3d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_domain_restart_variable_4d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_domain_restart_variable_5d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_0d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_1d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_2d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_3d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_4d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_unstructured_domain_restart_variable_5d (fileobj, variable_name, vdata, dimensions, is_optional, chunksizes)
 Add a domain decomposed variable. More...
 
subroutine register_variable_attribute_0d (fileobj, variable_name, attribute_name, attribute_value, str_len)
 Add an attribute to a variable. More...
 
subroutine register_variable_attribute_1d (fileobj, variable_name, attribute_name, attribute_value, str_len)
 Add an attribute to a variable. More...
 
subroutine unstructured_domain_read_0d (fileobj, variable_name, buf, unlim_dim_level, corner, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_read_1d (fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_read_2d (fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_read_3d (fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_read_4d (fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_read_5d (fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_0d (fileobj, variable_name, variable_data, unlim_dim_level, corner)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_1d (fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_2d (fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_3d (fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_4d (fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 
subroutine unstructured_domain_write_5d (fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
 Wrapper to distinguish interfaces. More...
 

Variables

integer deflate_level = default_deflate_level
 Netcdf deflate level to use in nf90_def_var (integer between 1 to 9)
 
logical, private fms2_io_is_initialized = .false.
 True after fms2_io_init is run.
 
integer header_buffer_val = 16384
 Use defined netCDF header buffer size(in bytes) used in NF__ENDDEF.
 
integer ncchksz = 64*1024
 User defined chunksize (in bytes) argument in netcdf file creation calls. Replaces setting the NC_CHKSZ environment variable.
 
character(len=10) netcdf_default_format = "64bit"
 User defined netcdf file format, acceptable values are: "64bit", "classic", "netcdf4". This can be overwritten if you specify "nc_format" in the open_file call.
 
logical shuffle = .false.
 Flag indicating whether to use the netcdf shuffle filter.
 

Detailed Description

An updated library for parallel IO to replace mpp_io_mod. This module contains the public API for fms2 I/O interfaces and routines defined throughout this directory. A conversion guide for replacing mpp/fms_io code with this module is available below.

readme

FMS_io/mpp_io to FMS2_io conversion guide

Before introducing the FMS2_io module, subroutines and functions in fms_io and mpp_io modules were used to mainly handle read/writes to NetCDF files. However, there were duplicate routines present in both modules that led to redundancy, and the blackbox-like I/O processes restricted user flexibility. The FMS2_io module has thus been implemented for a cleaner set of I/O tools and to give users more control over the information being written/read to NetCDF files. This guide helps convert fms_io/mpp_io code to FMS2_io

Contents

A. FMS2_io Fileobjs

FMS2_io provides three new derived types, which target the different I/O paradigms used in GFDL models.1. FmsNetcdfFile_t: This type provides a thin wrapper over the netCDF4 library, but allows the user to assign a “pelist” to the file. If a pelist is assigned, only the first rank on the list directly interacts with the netCDF library, and performs broadcasts to relay the information read to the rest of the ranks on the list. When writing netcdf files, only the first rank in the pelist will perform the writes.2. FmsNetcdfDomainFile_t: This type does everything that the FmsNetcdfFile_t type does and it adds support for “domain-decomposed” reads/writes. Here, "domain decomposed" refers to data that is on a user-defined mpp_domain and is decomposed in two dimensions, in which each MPI rank has its own section of the global data. This requires a domain to be associated with the fileobj.3. FmsNetcdfUnstructuredDomainFile_t: This type does everything that the FmsNetcdfFile_t type does and it adds support for “domain-decomposed” reads/writes on a user defined mpp_domains unstructured grid. This requires a unstructured domain to be associated with the fileobj.The FMS_io equivalent to these derived types is restart_file_type

B. Writing Restarts

1. Domain Decomposed Restarts

When writing domain decomposed restarts, FMS sets the filename for domain decomposed restarts in the following way:

  • If the domain is a cubesphere (6 tiles) and the io_layout is (/1,1/), it will create restart files: "RESTART/filename.tile{tile_number}.res.nc"
  • If the domain is a cubesphere (6 tiles) and the io_layout is (/x_layout,y_layout/), it write x_layout*y_layout restart files: "RESTART/filename.tile{tile_number}.res.nc.????"
  • If the domain has only 1 tile and the tile number is 1, "tile{tile_number}" will not be added
Additionally, the user can append an appendix (i.e "nest{nest_number}" or ens{ensemble_number}") to the filename by calling `set_filename_appendix (appendix)` - If the domain has more than 1 tile, the appendix will be added as "RESTART/filename.appendix.tile{tile_number}.res.nc" - If the domain has 1 tile, the appendix will added as "RESTART/filename.res.appendix.nc"This is how a restart file can be written with FMS_io:

use fms_io_mod, only: restart_file_type, register_restart_field, save_restart
use mpp_domains_mod, only: domain
integer :: id_restart !< Id for restart variable
type(restart_file_type) :: fileobj !< Fms_io fileobj
real, dimension(:,:,:) :: variable_data !< Variable data in the compute or data domain
type (domain2d) :: domain !< 2d mpp domain
id_restart = register_restart_field(fileobj, "filename", "variable_name", variable_data, domain=domain)
call save_restart(atm_restart)

Metadata:

  • FMS_io named the dimensions: "xaxis_1", "yaxis_1", zaxis_1" and "Time". With FMS2_io the user can name the axis whatever they like.
  • FMS_io wrote the dimensions as variables as well. FMS2_io does not do this by default, the user can do this if they like.
  • FMS_io wrote variable attribute: "longname = {same as variable} and "units = {"none"} to all variables by default. FMS2_io does not do this by default, the user can add real meta data if they like.
This is how a restart file can be written with FMS2-io:

use fms2_io_mod, only: fmsnetcdfdomainfile_t, register_restart_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_restart
use mpp_domains_mod, only: domain2d, center
type(FmsNetcdfDomainFile_t) :: fileobj !< Fms2_io domain decomposed fileobj
real, dimension(:,:,:,:) :: variable_data !< Variable data in the compute or data domain
type (domain2d) :: domain !< 2d mpp domain
character(len=8) :: dim_names(4) !< Array of dimension names
dim_names(1) = "xaxis_1"
dim_names(2) = "yaxis_1"
dim_names(3) = "zaxis_1"
dim_names(4) = "Time"
!> Create domain !<
!> Create an io_domain !<
if (open_file(fileobj, "filename", "overwrite", domain, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1), "x", position=center)
call register_axis(fileobj, dim_names(2), "y", position=center)
call register_axis(fileobj, dim_names(3), dimsize)
call register_axis(fileobj, dim_names(4), unlimited)
call register_restart_field(fileobj, 'variable_name', variable_data, dim_names)
call write_restart(fileobj)
call close_file(fileobj)
endif
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:166
Opens a given netcdf or domain file.
Definition: fms2_io.F90:122
Add a dimension to a given file.
Definition: fms2_io.F90:190
Similar to register_field, but occupies the field with data for restarts Example usage:
Definition: fms2_io.F90:231
Writes all restart fields registered within a given restart file Example usage:
Definition: fms2_io.F90:325
The domain2D type contains all the necessary information to define the global, compute and data domai...
  • open_file
    • a logical function, outputs .true. if the file was opened successfully and .false. if it failed.
    • Mangles the filename in the same manner as fms_io
    • Opens the netcdf file to write (a nf90_create call)
    • Set ups the pelist for each io_domain
    • is_restart indicates that this is a restart file, so it adds ".res" to the filename and it allows user to use the write_restart and register_restart_field functionality
    • NOTE: The filename needs to include the full path to the file, including the directory.
  • register_axis
    • writes the dimension metadata in the netcdf file (a nf90_def_dim call)
    • The "x" and "y" argument indicate that that dimension is domain decomposed in x/y. The only acceptable values are "x" and "y".
    • The position=center indicates the position of the axis (this is the default). The other acceptable values are position=east for "x" andposition=north for "y", in this cases the data is staggered, which may be on the volume face or corner instead of at the centroid.
    • The "unlimited" indicates that the dimension is unlimited (nf90_unlimited)
    • The integer "dimsize" indicates that this is not a domain decomposed dimension with a length equal to dimsize
  • register_restart_field
    • Writes the variable metadata to the file (a nf90_def_var call)
    • Saves the data as pointers, which will be written to the netcdf file later
  • write_restart
    • Loops through the restart variables that were registered
    • Calculates and writes a global checksum for each variables
    • Writes the data to the file (a nf90_put_var call)
  • close_file
    • Cleans up the fileobj (deallocates any arrays)
    • Closes the netcdf file (a nf90_close call)

2. Unstructured Domain Restarts

Restart files with domain decomposed variables in the unstructured domain can be written using the FmsNetcdfUnstructuredDomainFile_t fileobj.
use fms2_io_mod, only: fmsnetcdfunstructureddomainfile_t, register_restart_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_restart
use mpp_domains_mod, only: domainug, center
type(FmsNetcdfUnstructuredDomainFile_t) :: fileobj !< Fms2_io domain decomposed fileobj
real, dimension(:,:) :: variable_data !< Variable data in the unstructured domain
type (domainug) :: domain !< 2d mpp domain
character(len=8) :: dim_names(3) !< Array of dimension names
dim_names(1) = "ucdim"
dim_names(2) = "zaxis_1"
dim_names(3) = "Time"
!> Create domain !<
!> Create an io_domain !<
!> Create an unstructured domain !<
if (open_file(fileobj, "filename", "overwrite", domain, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1))
call register_axis(fileobj, dim_names(2), dimsize)
call register_axis(fileobj, dim_names(3), unlimited)
call register_restart_field(fileobj, 'variable_name', variable_data, dim_names)
call write_restart(fileobj)
call close_file(fileobj)
endif
Domain information for managing data on unstructured grids.
  • open_file
    • a logical function, outputs .true. if the file was opened successfully and .false. if it failed.
    • Mangles the filename: (1) Adds ".tileXX" if you are on multiple tiles (2) Adds .XXXX to the end of the file if uncombined.
    • Opens the netcdf file to write (a nf90_create call)
    • Set ups the pelist for each io_domain
    • is_restart indicates that this is a restart file, so it adds ".res" to the filename and it allows user to use the write_restart and register_restart_field functionality
    • NOTE: The filename needs to include the full path to the file, including the directory.
  • register_axis
    • writes the dimension metadata in the netcdf file (a nf90_def_dim call)
    • The lack of the third argument in the register_axis calls indicates that this dimension is in the unstructured domain
    • The "unlimited" indicates that the dimension is unlimited (nf90_unlimited)
    • The integer "dimsize" indicates that this is a normal dimension of length equal to dimsize
  • register_restart_field
    • Writes the variable metadata to the file (a nf90_def_var call)
    • Saves the data as pointers, which will be written to the netcdf file later
  • write_restart
    • Loops through the restart variables that were registered
    • Calculates and writes a global checksum for each variables
    • Writes the data to the file (a nf90_put_var call)
  • close_file
    • Cleans up the fileobj and deallocated any variables
    • Closes the netcdf file (a nf90_close call)

3. Non-domain Decomposed Restarts

Restart files without domain decomposed variables can be written using the FmsNetcdfFile_t fileobj.
use fms2_io_mod, only: fmsnetcdffile_t, register_restart_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_restart
use mpp_mod, only: mpp_npes, mpp_get_current_pelist
type(FmsNetcdfFile_t) :: fileobj !< Fms2_io fileobj
real, dimension(:,:,:,:) :: variable_data !< Variable data
character(len=8) :: dim_names(4) !< Array of dimension names
integer, allocatable, :: pes(:) !< Array of the pes in the current pelist
!< Get the current pelist
allocate(pes(mpp_npes()))
call mpp_get_current_pelist(pes)
dim_names(1) = "xaxis_1"
dim_names(2) = "yaxis_1"
dim_names(3) = "zaxis_1"
dim_names(4) = "Time"
if (open_file(fileobj, "filename", "overwrite", pelist=pes, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1), 96)
call register_axis(fileobj, dim_names(2), 96)
call register_axis(fileobj, dim_names(3), dimsize)
call register_axis(fileobj, dim_names(4), unlimited)
call register_restart_field(fileobj, 'variable_name', variable_data, dim_names)
call write_restart(fileobj)
call close_file(fileobj)
endif
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:421
  • open_file
    • a logical function, outputs .true. if the file was opened successfully and .false. if it failed.
    • Mangles the filename: ".res" is the only thing that is added to the filename
    • Opens the netcdf file to write (a nf90_create call)
    • is_restart indicates that this is a restart file, so it adds ".res" to the filename and it allows user to use the write_restart and register_restart_field functionality
    • With the optional pelist argument, only the first rank interacts with the file (opens, writes, closes) and broadcasts the information to the rest of the ranks on the list.
    • NOTE: The filename needs to include the full path to the file, including the directory.
  • register_axis
    • writes the dimension metadata in the netcdf file (a nf90_def_dim call)
    • The "unlimited" indicates that the dimension is unlimited (nf90_unlimited)
    • The integer "dimsize" indicates that this is a normal dimension of length equal to dimsize
    • "x" and "y" cannot be added to the register_axis calls
  • register_restart_field
    • Writes the variable metadata to the file (a nf90_def_var call)
    • Saves the data as pointers, which will be written to the netcdf file later
  • write_restart
    • Loops through the restart variables that were registered
    • Calculates and writes a global checksum for each variables
    • Writes the data to the file (a nf90_put_var call)
  • close_file
    • Cleans up the fileobj
    • Closes the netcdf file (a nf90_close call)

4. Other helpful information:

  • By default, FMS2_io writes the file as nf90_64bit_offset, the user can change the netcdf file type by adding nc_format="64bit", "classic", or "netcdf4" to the open_file call (which will change the type for that file only) or by adding netcdf_default_format="64bit", "classic", or "netcdf4" to the fms2_io_nml (which will change the type for all files that don't have the nc_format in the open_file call)
  • If the user wishes to not add ".res" to filename, the user can add dont_add_res_to_filename=.true. to the open_file call
  • Variable attributes can be written by calling register_variable_attribute. Scalar and 1d real and integers (32 and 64 bit) and string values are supported
    call register_variable_attribute(fileobj, "varname", "attribute_name", value)
    This interface can be used with any FMS2_io fileobj, but the open_file needs to be called before using it.
A similar thing can be accomplished in fms_io with mpp_write_meta
  • Global attributes can be written by calling register_global_attribute. Scalar and 1d real and integers (32 and 64 bit) and scalar string values are supported
    call register_global_attribute(fileobj, "global_attribute_name", value)
    This interface can be used with any FMS2_io fileobj, but the open_file needs to be called before using it.
A similar thing can be accomplished in fms_io with mpp_write_meta

C. Reading Restarts

The restarts can be read the same way as the writes. The only difference is that "read" is used in the open_file call and read_restart is used instead of write restart.Additionally, the register_restart_field has the optional argument, is_optional, where is .true. the code will not crash if the variable does not exist in the file when reading. *This argument is the same as mandatory in fms_io.

Other helpful information:

The following subroutines can be used with any of any of the FMS2_io fileobj.

  • Variable attributes can be read by calling get_variable_attribute. Scalar and 1d real and integers (32 and 64 bit) and string values are supported. To check if a variable attribute exists before reading, the logical function variable_att_exists can be used
    if (variable_att_exists(fileobj, "varname", "attribute_name")) then
    call get_variable_attribute(fileobj, "varname", "attribute_name", value)
    endif
    A similar thing can be accomplished in fms_io with get_vart_att_value
  • Global attributes can be read by calling get_global_attribute. Scalar and 1d real and integers (32 and 64 bit) and scalar string values are supported. To check if a global attribute exist before reading, the logical function global_att_exists can be used.
    if (global_att_exists(fileobj, "global_attribute_name)) then
    call get_global_attribute(fileobj, "global_attribute_name", value)
    A similar thing can be accomplished in fms_io with get_global_att_value
  • Reading dimension metadata

    Sample usage:

    ndims_file = get_num_dimensions(fileobj)
    allocate(dim_names_file(ndims_file))
    call get_dimension_names(fileobj, dim_names_file)
    call get_unlimited_dimension_name(fileobj, dimension_name)
    if (is_dimension_unlimited(fileobj, "dimension_name")) call mpp_error(fatal, "dimension_name is not an unlimited dimension")
    if (dimension_exists(fileobj, "dimension_name")) then
    call get_dimension_size(fileobj, "dimension_name", dimsize)
    endif
    Error handler.
    Definition: mpp.F90:382

    A similar thing can be accomplished with fms_io with mpp_get_info

  • Reading variable metadata
Sample usage:

nvar = get_num_variables(fileobj)
allocate(var_names(nvar))
call get_variable_names(fileobj, var_names)
if (variable_exists(fileobj, "varname")) then
ndims = get_variable_num_dimensions(fileobj, "varname")
allocate(dim_sizes(ndims))
call get_variable_size(fileobj, "varname", dim_sizes)
allocate(dim_names(ndims))
call get_variable_dimension_names(fileobj, "varname", dim_names)
endif

A similar thing can be accomplished with fms_io with mpp_get_info and/or field_size

D. Reading/Writting Non-restarts

Reading and writing netcdf files that are not restarts can be done using read_data and write_data calls.

1. Domain decomposed read/write:

use fms2_io_mod, only: fmsnetcdfdomainfile_t, register_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_data
use mpp_domains_mod, only: domain2d, center
type(FmsNetcdfDomainFile_t) :: fileobj !< Fms2_io domain decomposed fileobj
real, dimension(:,:,:,:) :: variable_data !< Variable data in the compute or data domain
type (domain2d) :: domain !< 2d mpp domain
character(len=8) :: dim_names(4) !< Array of dimension names
dim_names(1) = "xaxis_1"
dim_names(2) = "yaxis_1"
dim_names(3) = "zaxis_1"
dim_names(4) = "Time"
if (open_file(fileobj, "filename", "overwrite", domain)) then
call register_axis(fileobj, dim_names(1), "x", position=center)
call register_axis(fileobj, dim_names(2), "y", position=center)
call register_axis(fileobj, dim_names(3), dimsize)
call register_axis(fileobj, dim_names(4), unlimited)
call register_field(fileobj, 'variable_name', 'variable_type', dim_names)
call write_data(fileobj, 'variable_name', variable_data)
call close_file(fileobj)
endif
Defines a new field within the given file Example usage:
Definition: fms2_io.F90:212
Write data to a defined field within a file Example usage:
Definition: fms2_io.F90:262

Difference from writing restarts:

  • open_file does not have is_restart=.true., so the write_restart functionality cannot be used'
  • register_field is basically a wrapper for nf90_def_var which just adds the variable metadata to the file.
    • "variable_type" is a string which indicates the type you want the variable to be written as. The acceptable values are "int", "int64", "double", "float", and "char"
  • write_data
    • The ranks send their data to the io_root pe. The io_root pe receives the data and writes it to the file.
Similarly for domain reads:

if (open_file(fileobj, "filename", "read", domain)) then
call register_axis(fileobj, dim_names(1), "x", position=center)
call register_axis(fileobj, dim_names(2), "y", position=center)
call read_data(fileobj, 'variable_name', variable_data)
call close_file(fileobj)
endif
Read data from a defined field in a file.
Definition: fms2_io.F90:292
  • register_axis is required for the domain decomposed dimensions so the code knows which dimensions and therefore variables are domain decomposed.
  • register_field is not required because the code can get the dimensions information from the file.
  • read_data
    • The io_root pe reads the data and sends it to the other pes.

2. Unstructured Domain non-restart read/writes

use fms2_io_mod, only: fmsnetcdfunstructureddomainfile_t, register_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_data
use mpp_domains_mod, only: domainug, center
type(FmsNetcdfUnstructuredDomainFile_t) :: fileobj !< Fms2_io domain decomposed fileobj
real, dimension(:,:) :: variable_data !< Variable data in the unstructured domain
type (domainug) :: domain !< 2d mpp domain
character(len=8) :: dim_names(3) !< Array of dimension names
dim_names(1) = "ucdim"
dim_names(2) = "zaxis_1"
dim_names(3) = "Time"
!> Create domain !<
!> Create an io_domain !<
!> Create an unstructured domain !<
if (open_file(fileobj, "filename", "overwrite", domain, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1))
call register_axis(fileobj, dim_names(2), dimsize)
call register_axis(fileobj, dim_names(3), unlimited)
call register_field(fileobj, 'variable_name', 'variable_type', dim_names)
call write_data(fileobj, 'variable_name', variable_data)
call close_file(fileobj)
endif

Difference from writing restarts:

  • open_file does not have is_restart=.true., so the write_restart functionality cannot be used'
  • register_field is basically a wrapper for nf90_def_var which just adds the variable metadata to the file.
    • "variable_type" is a string which indicates the type you want the variable to be written as. The acceptable values are "int", "int64", "double", "float", and "char"
  • write_data
    • The ranks send their data to the io_root pe. The io_root pe receives the data and writes it to the file.
Similarly for reads:

if (open_file(fileobj, "filename", "read", domain, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1))
call write_data(fileobj, 'variable_name', variable_data)
call close_file(fileobj)
endif
  • register_axis is required for the dimensions that are on the unstructured domain, so the code knows which dimensions and therefore variables are in the unstructured domain
  • register_field is not required because the code can get the dimensions information from the file.
  • read_data
    • The io_root pe reads the data and sends it to the other pes.

3. Non-domain decomposed read/writes

Writing non-domain decomposed files is similar to writing domain decomposed files:
use fms2_io_mod, only: fmsnetcdffile_t, register_restart_field, register_axis, unlimited
use fms2_io_mod, only: open_file, close_file, write_restart
use mpp_mod, only: mpp_npes, mpp_get_current_pelist
type(FmsNetcdfFile_t) :: fileobj !< Fms2_io fileobj
real, dimension(:,:,:,:) :: variable_data !< Variable data
character(len=8) :: dim_names(4) !< Array of dimension names
integer, allocatable, :: pes(:) !< Array of the pes in the current pelist
!< Get the current pelist
allocate(pes(mpp_npes()))
call mpp_get_current_pelist(pes)
dim_names(1) = "xaxis_1"
dim_names(2) = "yaxis_1"
dim_names(3) = "zaxis_1"
dim_names(4) = "Time"
if (open_file(fileobj, "filename", "overwrite", pelist=pes, is_restart=.true.)) then
call register_axis(fileobj, dim_names(1), 96)
call register_axis(fileobj, dim_names(2), 96)
call register_axis(fileobj, dim_names(3), dimsize)
call register_axis(fileobj, dim_names(4), unlimited)
call register_field(fileobj, 'variable_name', 'variable_type', dim_names)
call write_data(fileobj)
call close_file(fileobj)
endif
  • open_file - the pelist argument is needed so that only the root pe writes the file
  • register_axis is needed so that the code knows the size of each dimension
  • register_field is needed so that the code knows what dimensions the variable is a function of
  • write_data: the root pe in the pelist writes the data
To read non-domain decomposed restarts:
if (open_file(fileobj, "filename", "read", pelist=pes, is_restart=.true.)) then
call read_data(fileobj, 'variable_name', variable_data)
call close_file(fileobj)
endif
  • Here the register_axis and register_field calls are not needed because the code can get the dimensions from the file
Additionally, the corner and edge_lengths optional arguments in the read_data calls can be used to read only a section of a file.

if (open_file(fileobj, "filename", "read", pelist=pes, is_restart=.true.)) then
call read_data(fileobj, 'variable_name', variable_data, corner=(/10,2,2,3/), edge_lengths=(/2, 3, 4, 1/))
call close_file(fileobj)
endif

Here the code is going to start reading at x=10, y=2, z=2, t=3 and read 2, 3, 4, 1 points in each direction. Note that variable_data must be the correct size or the code will fail.Additionally, it is possible to do "compressed" reads and writes. Here "compressed" is when different ranks have different numbers of points for a specific axis. For example:

if (open_file(fileobj, "filename", "overwrite", pelist=pes, is_restart=.true.)) then
call register_axis(fileobj, "xaxis_1", mpp_pe()+1, is_compressed=.true.)
call register_axis(fileobj, "yaxis_1", mpp_pe()+1, is_compressed=.true.)
call register_field(fileobj, 'variable_name', 'variable_type', dim_names)
call write_data(fileobj)
call close_file(fileobj)
endif
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
  • Because is_compressed=.true., the root pe is going to gather the dimension size of each rank and add them to get the total length of the dimension. In the example above, rank 0 has a x/y axis of 1, rank 1 has a x/y of 2. The total dimension size of x/y is equal to 3 (the sum for all pes).

E. Coupler Type Restarts

coupler_type_register_restarts permits interating a list of variables in the coupler_type and register each variable to the corresponding restart file and then read/write the file.

1. Reading coupler_type restarts

In FMS_io, this was accomplished as:
type(coupler_2d_bc_type) :: bc_type
type(restart_file_type), pointer :: bc_rest_files(:)=> null() !< Array of fms2_io fileobjs
call coupler_type_register_restarts(bc_type, bc_rest_files, num_rest_files, domain)
do l = 1, num_rest_files
call restore_state(bc_rest_files(l), directory='INPUT', &
nonfatal_missing_files=.true.)
enddo
In FMS2_io, this can be accomplished as:
type(coupler_2d_bc_type) :: bc_type
type(FmsNetcdfDomainFile_t), pointer :: bc_rest_files(:)=> null() !< Array of fms2_io fileobjs
call coupler_type_register_restarts(bc_type, bc_rest_files, num_rest_files, domain, to_read=.true.)
do l = 1, num_rest_files
call read_restart(bc_rest_files(l))
call close_file(bc_rest_files(l))
enddo
Reads in restart variables from a given file Example usage: call read_restart(fileobj) Reads registe...
Definition: fms2_io.F90:356

Note: num_rest_files is set inside coupler_type_register_restartsDifference from fms_io

  • to_read was added to the coupler_type_register_restarts subroutine
  • coupler_type_register_restarts loops through the variables in the coupler_type, opens the file needed, register the axis, and registers the restart_variables.
  • The file should be closed after writing it by calling close_file
  • It is required that the domain has an io_domain or the code will fail

2. Writting coupler type restarts

This is done the same way as the reads expect write_restart is used instead of read_restart.

F. Boundary Conditions Restarts

Both FMS_io and FMS2_io have the functionality where one can read or write data where only some pes have a section of the data (i.e halos).In FMS_io this was accomplished as:

type(restart_file_type) :: fileobj !< Fms_io fileobj
id_restart = register_restart_field(fileobj, "filename", "variable_name", variable_data, indices, global_size, pelist, is_root_pe)
call save_restart_border(fileobj)
call free_restart_type(fileobj)
In FMS2_io this can be accomplished as:

type(FmsNetcdfFile_t) :: fileobj !< Fms2_io fileobj
if (open_file(fileobj, "filename", "overwrite", is_restart=.true., pelist=all_pelist)) then
call register_restart_field(fileobj, "variable_name", variable_data, indices, global_size, pelist, is_root_pe)
call write_restart_bc(fileobj)
call close_file(fileobj)
endif
  • indices are the starting and indices of the region for that pe (starting x, ending x, starting y, ending y)
  • global_size is the size of the variable in (x,y,z)
  • pelist is the list of pelist that have the data for the variable
  • is_root_pe is a flag indicating if this is the root_pe from the pelist
Difference between FMS2_io and FMS_io:

G. Ascii_io

1. Writing ascii files

With fms_io, ascii writes can be accomplished with a combination of mpp_open, mpp_close from mpp_io and write.

integer :: unit_number
call mpp_open( unit_number, 'filename', nohdrs=.true. )
if ( mpp_pe().EQ.mpp_root_pe() ) write (unit_number, *) 'FMS2_io is flawless and works as inteded'
call mpp_close(unit_number)
With FMS2_io, there is no real support for ascii writes. The user can use fortran's open, write and close. For example,

integer :: unit_number
if ( mpp_pe().EQ.mpp_root_pe() ) then
open(newunit = unit_number, file='filename', status='replace', form='formatted')
write (unit_number, *) 'FMS2_io is flawless and works as inteded'
close (unit_number)
endif
  • The if ( mpp_pe().EQ.mpp_root_pe() ) ensures that the file is only written once.
  • The newunit ensures that an unused unit number is used.

2. Reading ascii files

With fms_io:

  • Opening an ascii file was accomplished with mpp_open, open_namelist_file or open calls.
  • Reading an ascii file was accomplished with read_data, read_distributed or read calls.
  • Closing an ascii file was sometimes done with mpp_close, close_file or close calls.
With FMS2_io, ascii reads can be accomplished using ascii_reads. In ascii_reads, the root_pe opens and reads the file into a a string buffer and broadcasts it to the other ranks. Each rank can then read from the buffer. For example:

character(len=:), dimension(:), allocatable :: restart_file !< Restart file saved as a string
call ascii_read('INPUT/coupler.res', restart_file)
read(restart_file(1), *) calendar_type
read(restart_file(2), *) date_init
read(restart_file(3), *) date
deallocate(restart_file)

H. FMS2_io namelist

  • ncchksz: Sets chunksize (in bytes) argument in netcdf file creation calls. The default is 64*1024.
  • netcdf_default_format: Sets the netcdf file type. The acceptable values are "64bit", "classic", "netcdf4". This can be overwritten per file if you specify nc_format in the open_file call. The default is 64bit.
  • header_buffer_val: Sets the netCDF header buffer size(in bytes). The default is 16384 bytes.
  • deflate_level: Determines how much to compress the variable. Chosen by an integer of 1 through 9. The higher the number the more compression will take place, but will take longer to write the file. The default is no compression. This is loseless compression so that every bit of the original data can be recovered. See the NETCDF user guide for more information: https://www.unidata.ucar.edu/blogs/developer/en/entry/netcdf_compression
  • shuffle: Flag indicating whether to use the netcdf shuffle filter.

I. Chunking

In release 2022.04, "chunksize" was added as an optional argument to the register_diag_field:
call register_restart_field(fileobj, 'variable_name', variable_data,
dim_names, chunsizes=chunksizes)
  • chunksizes: Is an array defining the chunksize of each dimension of the variable. Chunksizes can be used to improve performance. Default chunks are chosen by the library.
  • NOTE: This argument is only valid with "NETCDF4" file formats, otherwise it will be ignored. You can set the netcdf file format for all files using the 'netcdf_default_format' namelist or in a per file basis by using the 'nc_format' argument in the 'open_file' call
  • See the NETCDF user guide for more information: (https://cluster.earlham.edu/bccd-ng/testing/mobeen/GALAXSEEHPC/netcdf-4.1.3/man4/netcdf.html#Chunking)

Data Type Documentation

◆ fms2_io_mod::close_file

interface fms2_io_mod::close_file

Close a netcdf or domain file opened with open_file or open_virtual_file.


Example usage:

         call close_file(fileobj)

Closes any given fileobj opened via open_file or open_virtual_file

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 166 of file fms2_io.F90.

Private Member Functions

 close_domain_file
 
 close_unstructured_domain_file
 
 netcdf_file_close_wrap
 

◆ fms2_io_mod::open_file

interface fms2_io_mod::open_file

Opens a given netcdf or domain file.


Example usage:

         io_success = open_file(fileobj, "filename", "write")

Opens a netcdf file of type fmsnetcdffile_t at the given file path string. File mode is set to one of "read"/"write"/"overwrite"/"append"

         io_success = open_file(fileobj, "filename", "overwrite", domain)

Opens a domain netcdf file of type fmsnetcdfdomainfile_t or fmsnetcdfunstructureddomainfile_t at the given file path name and 2D or unstructured domain.

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 122 of file fms2_io.F90.

Private Member Functions

 netcdf_file_open_wrap
 
 open_domain_file
 
 open_unstructured_domain_file
 

◆ fms2_io_mod::open_virtual_file

interface fms2_io_mod::open_virtual_file

Creates a diskless netcdf or domain file.

Returns
true if successful, false otherwise


Example usage:

         io_success = open_virtual_file(fileobj, "filename", pelist)

Opens a virtual file through fmsnetcdffile_t at an optional file path and pelist

         io_success = open_virtual_file(fileobj, domain, "filename")

Opens a virtual domain file through fmsnetcdfdomainfile_t or fmsnetcdfunstructureddomainfile_t for a given 2D domain at an optional path

Note
For individual documentation on the listed routines, please see the appropriate helper module: blackboxio

Definition at line 146 of file fms2_io.F90.

Private Member Functions

 create_diskless_domain_file
 
 create_diskless_netcdf_file_wrap
 
 create_diskless_unstructured_domain_file
 

◆ fms2_io_mod::read_data

interface fms2_io_mod::read_data

Read data from a defined field in a file.


Example usage:

         call read_data(fileobj, "lat", data)

Read the values for the field "lat" from the file and write them onto data

Definition at line 292 of file fms2_io.F90.

Private Member Functions

 compressed_read_0d
 
 compressed_read_1d
 
 compressed_read_2d
 
 compressed_read_3d
 
 compressed_read_4d
 
 compressed_read_5d
 
 domain_read_0d
 
 domain_read_1d
 
 domain_read_2d
 
 domain_read_3d
 
 domain_read_4d
 
 domain_read_5d
 
 unstructured_domain_read_0d
 
 unstructured_domain_read_1d
 
 unstructured_domain_read_2d
 
 unstructured_domain_read_3d
 
 unstructured_domain_read_4d
 
 unstructured_domain_read_5d
 

◆ fms2_io_mod::read_new_restart

interface fms2_io_mod::read_new_restart

Read registered restarts from a new file Optionally takes directory to write to, model time and filename
Example usage: call read_new_restart(fileobj, unlimted_dimension_level)

call read_new_restart(fileobj, unlimited_dimension_level, directory, timestamp, filename)

Note
For individual documentation on the listed routines, please see the appropriate helper module: blackboxio

Definition at line 369 of file fms2_io.F90.

Private Member Functions

 netcdf_restore_state_wrap
 
 restore_domain_state_wrap
 

◆ fms2_io_mod::read_restart

interface fms2_io_mod::read_restart

Reads in restart variables from a given file
Example usage: call read_restart(fileobj) Reads registered restart variables from fileobj.

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 356 of file fms2_io.F90.

Private Member Functions

 netcdf_restore_state
 
 restore_domain_state
 

◆ fms2_io_mod::register_axis

interface fms2_io_mod::register_axis

Add a dimension to a given file.


Example usage:

         call register_axis(fileobj, "lon", "x")

Adds a dimension named "lon" associated with the x axis of the 2D domain file. For unstructured domains no x or y axis character is provided.

         call register_axis(fileobj, "lon", n)

Adds a dimension named "lon" with length n to a given netcdf file.

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 190 of file fms2_io.F90.

Private Member Functions

 netcdf_add_dimension
 
 register_compressed_dimension
 
 register_domain_decomposed_dimension
 
 register_unstructured_dimension
 

◆ fms2_io_mod::register_field

interface fms2_io_mod::register_field

Defines a new field within the given file
Example usage:

         call register_field(fileobj, "lon", "double", (/"lon"/) )

Adds a double variable named "lon" to the given file, corresponding to the list of dimension names (which must be previously defined in the fileobj). The size of dimension name list provided is the amount of ranks for the created field, scalar if list not provided.

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 212 of file fms2_io.F90.

Private Member Functions

 netcdf_add_variable_wrap
 
 register_domain_variable
 
 register_unstructured_domain_variable
 

◆ fms2_io_mod::register_restart_field

interface fms2_io_mod::register_restart_field

Similar to register_field, but occupies the field with data for restarts
Example usage:

         call register_restart_field(fileobj, "temperature", data_ptr, (/"lon", "time"/) )

Creates a restart variable and sets it to the values from data_ptr, corresponding to the list of dimension names. Rank of data_ptr must equal the amount of corresponding dimensions.

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 231 of file fms2_io.F90.

Private Member Functions

 netcdf_add_restart_variable_0d_wrap
 
 netcdf_add_restart_variable_1d_wrap
 
 netcdf_add_restart_variable_2d_wrap
 
 netcdf_add_restart_variable_3d_wrap
 
 netcdf_add_restart_variable_4d_wrap
 
 netcdf_add_restart_variable_5d_wrap
 
 register_domain_restart_variable_0d
 
 register_domain_restart_variable_1d
 
 register_domain_restart_variable_2d
 
 register_domain_restart_variable_3d
 
 register_domain_restart_variable_4d
 
 register_domain_restart_variable_5d
 
 register_restart_region_2d
 
 register_restart_region_3d
 
 register_unstructured_domain_restart_variable_0d
 
 register_unstructured_domain_restart_variable_1d
 
 register_unstructured_domain_restart_variable_2d
 
 register_unstructured_domain_restart_variable_3d
 
 register_unstructured_domain_restart_variable_4d
 
 register_unstructured_domain_restart_variable_5d
 

◆ fms2_io_mod::write_data

interface fms2_io_mod::write_data

Write data to a defined field within a file
Example usage:

         call write_data(fileobj, "lon", data)

Write the value(s) in data to the field named "lon"

Definition at line 262 of file fms2_io.F90.

Private Member Functions

 compressed_write_0d_wrap
 
 compressed_write_1d_wrap
 
 compressed_write_2d_wrap
 
 compressed_write_3d_wrap
 
 compressed_write_4d_wrap
 
 compressed_write_5d_wrap
 
 domain_write_0d
 
 domain_write_1d
 
 domain_write_2d
 
 domain_write_3d
 
 domain_write_4d
 
 domain_write_5d
 
 unstructured_domain_write_0d
 
 unstructured_domain_write_1d
 
 unstructured_domain_write_2d
 
 unstructured_domain_write_3d
 
 unstructured_domain_write_4d
 
 unstructured_domain_write_5d
 

◆ fms2_io_mod::write_new_restart

interface fms2_io_mod::write_new_restart

Writes all restart fields in a given restart file to a new restart file
Example usage:

         call write_new_restart(fileobj, timestamp="tstring", filename="new_restartfilename")

Creates a new restart file, with the provided timestamp and filename, out of the registered restart fields in the given restart file.

Note
For individual documentation on the listed routines, please see the appropriate helper module: blackboxio

Definition at line 341 of file fms2_io.F90.

Private Member Functions

 netcdf_save_restart_wrap2
 
 save_domain_restart_wrap
 
 unstructured_write_restart_wrap
 

◆ fms2_io_mod::write_restart

interface fms2_io_mod::write_restart

Writes all restart fields registered within a given restart file
Example usage:

         call write_restart(fileobj)

Writes previously registered restart fields to the given restart file

Note
For individual documentation on the listed routines, please see the appropriate helper module. For netcdf files with a structured domain: fms_netcdf_domain_io_mod. For netcdf files with an unstructured domain: fms_netcdf_unstructured_domain_io_mod. For generic netcdf: netcdf_io_mod.

Definition at line 325 of file fms2_io.F90.

Private Member Functions

 netcdf_save_restart_wrap
 
 save_domain_restart
 
 unstructured_write_restart
 

Function/Subroutine Documentation

◆ compressed_read_0d()

subroutine compressed_read_0d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.

Definition at line 28 of file compressed_read.inc.

◆ compressed_read_1d()

subroutine compressed_read_1d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 53 of file compressed_read.inc.

◆ compressed_read_2d()

subroutine compressed_read_2d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 81 of file compressed_read.inc.

◆ compressed_read_3d()

subroutine compressed_read_3d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 109 of file compressed_read.inc.

◆ compressed_read_4d()

subroutine compressed_read_4d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 137 of file compressed_read.inc.

◆ compressed_read_5d()

subroutine compressed_read_5d ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(inout)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths 
)

I/O domain reads in data from the netcdf file and broadcasts the data to the rest of the ranks. This routine may only be used with variables that are "compressed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]cdataBuffer where the read in data will be stored.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 165 of file compressed_read.inc.

◆ compressed_write_0d()

subroutine compressed_write_0d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.

Definition at line 29 of file compressed_write.inc.

◆ compressed_write_0d_wrap()

subroutine compressed_write_0d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.

Definition at line 57 of file compressed_write.inc.

◆ compressed_write_1d()

subroutine compressed_write_1d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(1), intent(in), optional  corner,
integer, dimension(1), intent(in), optional  edge_lengths 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 81 of file compressed_write.inc.

◆ compressed_write_1d_wrap()

subroutine compressed_write_1d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(1), intent(in), optional  corner,
integer, dimension(1), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 641 of file compressed_write.inc.

◆ compressed_write_2d()

subroutine compressed_write_2d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(2), intent(in), optional  corner,
integer, dimension(2), intent(in), optional  edge_lengths 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 180 of file compressed_write.inc.

◆ compressed_write_2d_wrap()

subroutine compressed_write_2d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(2), intent(in), optional  corner,
integer, dimension(2), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 668 of file compressed_write.inc.

◆ compressed_write_3d()

subroutine compressed_write_3d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in), target, contiguous  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(3), intent(in), optional  corner,
integer, dimension(3), intent(in), optional  edge_lengths 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 300 of file compressed_write.inc.

◆ compressed_write_3d_wrap()

subroutine compressed_write_3d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(3), intent(in), optional  corner,
integer, dimension(3), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 694 of file compressed_write.inc.

◆ compressed_write_4d()

subroutine compressed_write_4d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(4), intent(in), optional  corner,
integer, dimension(4), intent(in), optional  edge_lengths 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 428 of file compressed_write.inc.

◆ compressed_write_4d_wrap()

subroutine compressed_write_4d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(4), intent(in), optional  corner,
integer, dimension(4), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 720 of file compressed_write.inc.

◆ compressed_write_5d()

subroutine compressed_write_5d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(5), intent(in), optional  corner,
integer, dimension(5), intent(in), optional  edge_lengths 
)

For variables without a compressed dimension, this routine simply wraps netcdf_write data. If the variable does have a compressed axis, the I/O root gathers the data from the rest of the ranks and then writes the combined data to the netcdf file.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 536 of file compressed_write.inc.

◆ compressed_write_5d_wrap()

subroutine compressed_write_5d_wrap ( type(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in)  cdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(5), intent(in), optional  corner,
integer, dimension(5), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]cdataCompressed data that will be gathered and written to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 746 of file compressed_write.inc.

◆ domain_read_0d()

subroutine domain_read_0d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(inout)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.

Definition at line 30 of file domain_read.inc.

◆ domain_read_1d()

subroutine domain_read_1d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(inout)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(1), intent(in), optional  corner,
integer, dimension(1), intent(in), optional  edge_lengths 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.
[in]edge_lengthsThe number of elements that will be read in each dimension.

Definition at line 56 of file domain_read.inc.

◆ domain_read_2d()

subroutine domain_read_2d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(inout), target, contiguous  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(2), intent(in), optional  corner,
integer, dimension(2), intent(in), optional  edge_lengths 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.
[in]edge_lengthsThe number of elements that will be read in each dimension.

Definition at line 87 of file domain_read.inc.

◆ domain_read_3d()

subroutine domain_read_3d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(inout), target, contiguous  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(3), intent(in), optional  corner,
integer, dimension(3), intent(in), optional  edge_lengths 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.
[in]edge_lengthsThe number of elements that will be read in each dimension.

Definition at line 267 of file domain_read.inc.

◆ domain_read_4d()

subroutine domain_read_4d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(inout)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(4), intent(in), optional  corner,
integer, dimension(4), intent(in), optional  edge_lengths 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.
[in]edge_lengthsThe number of elements that will be read in each dimension.

Definition at line 448 of file domain_read.inc.

◆ domain_read_5d()

subroutine domain_read_5d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(inout)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(5), intent(in), optional  corner,
integer, dimension(5), intent(in), optional  edge_lengths 
)

I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and scatters the data to the rest of the ranks using its I/O compute domain indices. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]vdataData that will be read out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be read to.
[in]edge_lengthsThe number of elements that will be read in each dimension.

Definition at line 669 of file domain_read.inc.

◆ domain_write_0d()

subroutine domain_write_0d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.

Definition at line 29 of file domain_write.inc.

◆ domain_write_1d()

subroutine domain_write_1d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(1), intent(in), optional  corner,
integer, dimension(1), intent(in), optional  edge_lengths 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 53 of file domain_write.inc.

◆ domain_write_2d()

subroutine domain_write_2d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(2), intent(in), optional  corner,
integer, dimension(2), intent(in), optional  edge_lengths 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 83 of file domain_write.inc.

◆ domain_write_3d()

subroutine domain_write_3d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(3), intent(in), optional  corner,
integer, dimension(3), intent(in), optional  edge_lengths 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 379 of file domain_write.inc.

◆ domain_write_4d()

subroutine domain_write_4d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(4), intent(in), optional  corner,
integer, dimension(4), intent(in), optional  edge_lengths 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 675 of file domain_write.inc.

◆ domain_write_5d()

subroutine domain_write_5d ( type(fmsnetcdfdomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in)  vdata,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(5), intent(in), optional  corner,
integer, dimension(5), intent(in), optional  edge_lengths 
)

Gather "compute" domain data on the I/O root rank and then have the I/O root write out the data that spans the "global" domain. This routine may only be used with variables that are "domain decomposed".

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataData that will be written out to the netcdf file.
[in]unlim_dim_levelLevel for the unlimited dimension.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 971 of file domain_write.inc.

◆ fms2_io_init()

subroutine, public fms2_io_mod::fms2_io_init

Reads the fms2_io_nml.

Check if the module has already been initialized

Call initialization routines that this module depends on

Read the namelist

Send the namelist variables to their respective modules

Mark the fms2_io as initialized

Definition at line 395 of file fms2_io.F90.

◆ register_domain_restart_variable_0d()

subroutine register_domain_restart_variable_0d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 27 of file register_domain_restart_variable.inc.

◆ register_domain_restart_variable_1d()

subroutine register_domain_restart_variable_1d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 47 of file register_domain_restart_variable.inc.

◆ register_domain_restart_variable_2d()

subroutine register_domain_restart_variable_2d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 70 of file register_domain_restart_variable.inc.

◆ register_domain_restart_variable_3d()

subroutine register_domain_restart_variable_3d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 90 of file register_domain_restart_variable.inc.

◆ register_domain_restart_variable_4d()

subroutine register_domain_restart_variable_4d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 110 of file register_domain_restart_variable.inc.

◆ register_domain_restart_variable_5d()

subroutine register_domain_restart_variable_5d ( type(fmsnetcdfdomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 130 of file register_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_0d()

subroutine register_unstructured_domain_restart_variable_0d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 27 of file register_unstructured_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_1d()

subroutine register_unstructured_domain_restart_variable_1d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 47 of file register_unstructured_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_2d()

subroutine register_unstructured_domain_restart_variable_2d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 67 of file register_unstructured_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_3d()

subroutine register_unstructured_domain_restart_variable_3d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 87 of file register_unstructured_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_4d()

subroutine register_unstructured_domain_restart_variable_4d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 107 of file register_unstructured_domain_restart_variable.inc.

◆ register_unstructured_domain_restart_variable_5d()

subroutine register_unstructured_domain_restart_variable_5d ( type(fmsnetcdfunstructureddomainfile_t), intent(inout)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in), target  vdata,
character(len=*), dimension(:), intent(in), optional  dimensions,
logical, intent(in), optional  is_optional,
integer, dimension(:), intent(in), optional  chunksizes 
)

Add a domain decomposed variable.

Parameters
[in,out]fileobjFile object.
[in]variable_nameVariable name.
[in]vdataVariable data.
[in]dimensionsDimension names.
[in]is_optionalPrevent errors in read-only files if a variable does not exist.
[in]chunksizesnetcdf chunksize to use for this variable This feature is only available for netcdf4 files

Definition at line 127 of file register_unstructured_domain_restart_variable.inc.

◆ register_variable_attribute_0d()

subroutine register_variable_attribute_0d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
character(len=*), intent(in)  attribute_name,
class(*), intent(in)  attribute_value,
integer, intent(in), optional  str_len 
)

Add an attribute to a variable.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]attribute_nameAttribute name.
[in]attribute_valueAttribute value
[in]str_lenLength of the string

Definition at line 27 of file register_variable_attribute.inc.

◆ register_variable_attribute_1d()

subroutine register_variable_attribute_1d ( class(fmsnetcdffile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
character(len=*), intent(in)  attribute_name,
class(*), dimension(:), intent(in)  attribute_value,
integer, intent(in), optional  str_len 
)

Add an attribute to a variable.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]attribute_nameAttribute name.
[in]attribute_valueAttribute value
[in]str_lenLength of the string

Definition at line 139 of file register_variable_attribute.inc.

◆ unstructured_domain_read_0d()

subroutine unstructured_domain_read_0d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 26 of file unstructured_domain_read.inc.

◆ unstructured_domain_read_1d()

subroutine unstructured_domain_read_1d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]edge_lengthsThe number of elements that will be read in each dimension.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 56 of file unstructured_domain_read.inc.

◆ unstructured_domain_read_2d()

subroutine unstructured_domain_read_2d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]edge_lengthsThe number of elements that will be read in each dimension.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 89 of file unstructured_domain_read.inc.

◆ unstructured_domain_read_3d()

subroutine unstructured_domain_read_3d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]edge_lengthsThe number of elements that will be read in each dimension.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 122 of file unstructured_domain_read.inc.

◆ unstructured_domain_read_4d()

subroutine unstructured_domain_read_4d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]edge_lengthsThe number of elements that will be read in each dimension.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 155 of file unstructured_domain_read.inc.

◆ unstructured_domain_read_5d()

subroutine unstructured_domain_read_5d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(inout)  buf,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(:), intent(in), optional  corner,
integer, dimension(:), intent(in), optional  edge_lengths,
logical, intent(in), optional  broadcast 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in,out]bufArray that the data will be read into.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be read from.
[in]edge_lengthsThe number of elements that will be read in each dimension.
[in]broadcastFlag controlling whether or not the data will be broadcasted to non "I/O root" ranks. The broadcast will be done by default.

Definition at line 188 of file unstructured_domain_read.inc.

◆ unstructured_domain_write_0d()

subroutine unstructured_domain_write_0d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, intent(in), optional  corner 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.

Definition at line 27 of file unstructured_domain_write.inc.

◆ unstructured_domain_write_1d()

subroutine unstructured_domain_write_1d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(1), intent(in), optional  corner,
integer, dimension(1), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 46 of file unstructured_domain_write.inc.

◆ unstructured_domain_write_2d()

subroutine unstructured_domain_write_2d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(2), intent(in), optional  corner,
integer, dimension(2), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 70 of file unstructured_domain_write.inc.

◆ unstructured_domain_write_3d()

subroutine unstructured_domain_write_3d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(3), intent(in), optional  corner,
integer, dimension(3), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 94 of file unstructured_domain_write.inc.

◆ unstructured_domain_write_4d()

subroutine unstructured_domain_write_4d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(4), intent(in), optional  corner,
integer, dimension(4), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 118 of file unstructured_domain_write.inc.

◆ unstructured_domain_write_5d()

subroutine unstructured_domain_write_5d ( type(fmsnetcdfunstructureddomainfile_t), intent(in)  fileobj,
character(len=*), intent(in)  variable_name,
class(*), dimension(:,:,:,:,:), intent(in)  variable_data,
integer, intent(in), optional  unlim_dim_level,
integer, dimension(5), intent(in), optional  corner,
integer, dimension(5), intent(in), optional  edge_lengths 
)

Wrapper to distinguish interfaces.

Parameters
[in]fileobjFile object.
[in]variable_nameVariable name.
[in]variable_dataData that will be written.
[in]unlim_dim_levelUnlimited dimension level.
[in]cornerArray of starting indices describing where the data will be written to.
[in]edge_lengthsThe number of elements that will be written in each dimension.

Definition at line 142 of file unstructured_domain_write.inc.