netCDF and HDF5#
netCDF#
The recommended way to store xarray data structures is netCDF, which
is a binary file format for self-described datasets that originated
in the geosciences. Xarray is based on the netCDF data model, so netCDF files
on disk directly correspond to Dataset
objects (more accurately,
a group in a netCDF file directly corresponds to a Dataset
object.
See Groups for more.)
NetCDF is supported on almost all platforms, and parsers exist for the vast majority of scientific programming languages. Recent versions of netCDF are based on the even more widely used HDF5 file-format.
Tip
If you aren’t familiar with this data format, the netCDF FAQ is a good place to start.
Reading and writing netCDF files with xarray requires scipy, h5netcdf, or the netCDF4-Python library to be installed. SciPy only supports reading and writing of netCDF V3 files.
We can save a Dataset to disk using the
Dataset.to_netcdf()
method:
ds = xr.Dataset(
{"foo": (("x", "y"), np.random.rand(4, 5))},
coords={
"x": [10, 20, 30, 40],
"y": pd.date_range("2000-01-01", periods=5),
"z": ("x", list("abcd")),
},
)
ds.to_netcdf("saved_on_disk.nc")
By default, the file is saved as netCDF4 (assuming netCDF4-Python is
installed). You can control the format and engine used to write the file with
the format
and engine
arguments.
Tip
Using the h5netcdf package
by passing engine='h5netcdf'
to open_dataset()
can
sometimes be quicker than the default engine='netcdf4'
that uses the
netCDF4 package.
We can load netCDF files to create a new Dataset using
open_dataset()
:
ds_disk = xr.open_dataset("saved_on_disk.nc")
ds_disk
<xarray.Dataset> Size: 248B Dimensions: (x: 4, y: 5) Coordinates: * x (x) int64 32B 10 20 30 40 * y (y) datetime64[ns] 40B 2000-01-01 2000-01-02 ... 2000-01-05 z (x) <U1 16B ... Data variables: foo (x, y) float64 160B ...
Similarly, a DataArray can be saved to disk using the
DataArray.to_netcdf()
method, and loaded
from disk using the open_dataarray()
function. As netCDF files
correspond to Dataset
objects, these functions internally
convert the DataArray
to a Dataset
before saving, and then convert back
when loading, ensuring that the DataArray
that is loaded is always exactly
the same as the one that was saved.
A dataset can also be loaded or written to a specific group within a netCDF
file. To load from a group, pass a group
keyword argument to the
open_dataset
function. The group can be specified as a path-like
string, e.g., to access subgroup ‘bar’ within group ‘foo’ pass
‘/foo/bar’ as the group
argument. When writing multiple groups in one file,
pass mode='a'
to to_netcdf
to ensure that each call does not delete the
file.
Tip
It is recommended to use DataTree
to represent
hierarchical data, and to use the xarray.DataTree.to_netcdf()
method
when writing hierarchical data to a netCDF file.
Data is always loaded lazily from netCDF files. You can manipulate, slice and subset Dataset and DataArray objects, and no array values are loaded into memory until you try to perform some sort of actual computation. For an example of how these lazy arrays work, see the OPeNDAP section below.
There may be minor differences in the Dataset
object returned
when reading a NetCDF file with different engines.
It is important to note that when you modify values of a Dataset, even one linked to files on disk, only the in-memory copy you are manipulating in xarray is modified: the original file on disk is never touched.
Tip
Xarray’s lazy loading of remote or on-disk datasets is often but not always
desirable. Before performing computationally intense operations, it is
often a good idea to load a Dataset (or DataArray) entirely into memory by
invoking the Dataset.load()
method.
Datasets have a Dataset.close()
method to close the associated
netCDF file. However, it’s often cleaner to use a with
statement:
# this automatically closes the dataset after use
with xr.open_dataset("saved_on_disk.nc") as ds:
print(ds.keys())
KeysView(<xarray.Dataset> Size: 248B
Dimensions: (x: 4, y: 5)
Coordinates:
* x (x) int64 32B 10 20 30 40
* y (y) datetime64[ns] 40B 2000-01-01 2000-01-02 ... 2000-01-05
z (x) <U1 16B ...
Data variables:
foo (x, y) float64 160B ...)
Although xarray provides reasonable support for incremental reads of files on disk, it does not support incremental writes, which can be a useful strategy for dealing with datasets too big to fit into memory. Instead, xarray integrates with dask.array (see Parallel Computing with Dask), which provides a fully featured engine for streaming computation.
It is possible to append or overwrite netCDF variables using the mode='a'
argument. When using this option, all variables in the dataset will be written
to the original netCDF file, regardless if they exist in the original dataset.
Groups#
Whilst netCDF groups can only be loaded individually as Dataset
objects, a
whole file of many nested groups can be loaded as a single
xarray.DataTree
object. To open a whole netCDF file as a tree of groups
use the xarray.open_datatree()
function. To save a DataTree object as a
netCDF file containing many groups, use the xarray.DataTree.to_netcdf()
method.
Note
Due to file format specifications the on-disk root group name is always "/"
,
overriding any given DataTree
root node name.
Warning
DataTree
objects do not follow the exact same data model as netCDF
files, which means that perfect round-tripping is not always possible.
In particular in the netCDF data model dimensions are entities that can exist regardless of whether any variable possesses them. This is in contrast to xarray’s data model (and hence DataTree’s data model) in which the dimensions of a (Dataset/Tree) object are simply the set of dimensions present across all variables in that dataset.
This means that if a netCDF file contains dimensions but no variables which possess those dimensions, these dimensions will not be present when that file is opened as a DataTree object. Saving this DataTree object to file will therefore not preserve these “unused” dimensions.
Reading encoded data#
NetCDF files follow some conventions for encoding datetime arrays (as numbers
with a “units” attribute) and for packing and unpacking data (as
described by the “scale_factor” and “add_offset” attributes). If the argument
decode_cf=True
(default) is given to open_dataset()
, xarray will attempt
to automatically decode the values in the netCDF objects according to
CF conventions. Sometimes this will fail, for example, if a variable
has an invalid “units” or “calendar” attribute. For these cases, you can
turn this decoding off manually.
You can view this encoding information (among others) in the
DataArray.encoding
and
DataArray.encoding
attributes:
ds_disk["y"].encoding
{'dtype': dtype('int64'),
'zlib': False,
'szip': False,
'zstd': False,
'bzip2': False,
'blosc': False,
'shuffle': False,
'complevel': 0,
'fletcher32': False,
'contiguous': True,
'chunksizes': None,
'source': '/home/docs/checkouts/readthedocs.org/user_builds/xray/checkouts/10526/doc/saved_on_disk.nc',
'original_shape': (5,),
'units': 'days since 2000-01-01 00:00:00',
'calendar': 'proleptic_gregorian'}
ds_disk.encoding
{'unlimited_dims': set(),
'source': '/home/docs/checkouts/readthedocs.org/user_builds/xray/checkouts/10526/doc/saved_on_disk.nc'}
Note that all operations that manipulate variables other than indexing will remove encoding information.
In some cases it is useful to intentionally reset a dataset’s original encoding values.
This can be done with either the Dataset.drop_encoding()
or
DataArray.drop_encoding()
methods.
ds_no_encoding = ds_disk.drop_encoding()
ds_no_encoding.encoding
{}
Reading multi-file datasets#
NetCDF files are often encountered in collections, e.g., with different files
corresponding to different model runs or one file per timestamp.
Xarray can straightforwardly combine such files into a single Dataset by making use of
concat()
, merge()
, combine_nested()
and
combine_by_coords()
. For details on the difference between these
functions see Combining data.
Xarray includes support for manipulating datasets that don’t fit into memory
with dask. If you have dask installed, you can open multiple files
simultaneously in parallel using open_mfdataset()
:
xr.open_mfdataset('my/files/*.nc', parallel=True)
This function automatically concatenates and merges multiple files into a
single xarray dataset.
It is the recommended way to open multiple files with xarray.
For more details on parallel reading, see Combining along multiple dimensions, Reading and writing data and a
blog post by Stephan Hoyer.
open_mfdataset()
takes many kwargs that allow you to
control its behaviour (for e.g. parallel
, combine
, compat
, join
, concat_dim
).
See its docstring for more details.
Note
A common use-case involves a dataset distributed across a large number of files with
each file containing a large number of variables. Commonly, a few of these variables
need to be concatenated along a dimension (say "time"
), while the rest are equal
across the datasets (ignoring floating point differences). The following command
with suitable modifications (such as parallel=True
) works well with such datasets:
xr.open_mfdataset('my/files/*.nc', concat_dim="time", combine="nested",
data_vars='minimal', coords='minimal', compat='override')
This command concatenates variables along the "time"
dimension, but only those that
already contain the "time"
dimension (data_vars='minimal', coords='minimal'
).
Variables that lack the "time"
dimension are taken from the first dataset
(compat='override'
).
Sometimes multi-file datasets are not conveniently organized for easy use of open_mfdataset()
.
One can use the preprocess
argument to provide a function that takes a dataset
and returns a modified Dataset.
open_mfdataset()
will call preprocess
on every dataset
(corresponding to each file) prior to combining them.
If open_mfdataset()
does not meet your needs, other approaches are possible.
The general pattern for parallel reading of multiple files
using dask, modifying those datasets and then combining into a single Dataset
is:
def modify(ds):
# modify ds here
return ds
# this is basically what open_mfdataset does
open_kwargs = dict(decode_cf=True, decode_times=False)
open_tasks = [dask.delayed(xr.open_dataset)(f, **open_kwargs) for f in file_names]
tasks = [dask.delayed(modify)(task) for task in open_tasks]
datasets = dask.compute(tasks) # get a list of xarray.Datasets
combined = xr.combine_nested(datasets) # or some combination of concat, merge
As an example, here’s how we could approximate MFDataset
from the netCDF4
library:
from glob import glob
import xarray as xr
def read_netcdfs(files, dim):
# glob expands paths with * to a list of files, like the unix shell
paths = sorted(glob(files))
datasets = [xr.open_dataset(p) for p in paths]
combined = xr.concat(datasets, dim)
return combined
combined = read_netcdfs('/all/my/files/*.nc', dim='time')
This function will work in many cases, but it’s not very robust. First, it never closes files, which means it will fail if you need to load more than a few thousand files. Second, it assumes that you want all the data from each file and that it can all fit into memory. In many situations, you only need a small subset or an aggregated summary of the data from each file.
Here’s a slightly more sophisticated example of how to remedy these deficiencies:
def read_netcdfs(files, dim, transform_func=None):
def process_one_path(path):
# use a context manager, to ensure the file gets closed after use
with xr.open_dataset(path) as ds:
# transform_func should do some sort of selection or
# aggregation
if transform_func is not None:
ds = transform_func(ds)
# load all data from the transformed dataset, to ensure we can
# use it after closing each original file
ds.load()
return ds
paths = sorted(glob(files))
datasets = [process_one_path(p) for p in paths]
combined = xr.concat(datasets, dim)
return combined
# here we suppose we only care about the combined mean of each file;
# you might also use indexing operations like .sel to subset datasets
combined = read_netcdfs('/all/my/files/*.nc', dim='time',
transform_func=lambda ds: ds.mean())
This pattern works well and is very robust. We’ve used similar code to process tens of thousands of files constituting 100s of GB of data.
Writing encoded data#
Conversely, you can customize how xarray writes netCDF files on disk by
providing explicit encodings for each dataset variable. The encoding
argument takes a dictionary with variable names as keys and variable specific
encodings as values. These encodings are saved as attributes on the netCDF
variables on disk, which allows xarray to faithfully read encoded data back into
memory.
It is important to note that using encodings is entirely optional: if you do not
supply any of these encoding options, xarray will write data to disk using a
default encoding, or the options in the encoding
attribute, if set.
This works perfectly fine in most cases, but encoding can be useful for
additional control, especially for enabling compression.
In the file on disk, these encodings are saved as attributes on each variable, which allow xarray and other CF-compliant tools for working with netCDF files to correctly read the data.
Scaling and type conversions#
These encoding options (based on CF Conventions on packed data) work on any version of the netCDF file format:
dtype
: Any valid NumPy dtype or string convertible to a dtype, e.g.,'int16'
or'float32'
. This controls the type of the data written on disk._FillValue
: Values ofNaN
in xarray variables are remapped to this value when saved on disk. This is important when converting floating point with missing values to integers on disk, becauseNaN
is not a valid value for integer dtypes. By default, variables with float types are attributed a_FillValue
ofNaN
in the output file, unless explicitly disabled with an encoding{'_FillValue': None}
.scale_factor
andadd_offset
: Used to convert from encoded data on disk to to the decoded data in memory, according to the formuladecoded = scale_factor * encoded + add_offset
. Please note thatscale_factor
andadd_offset
must be of same type and determine the type of the decoded data.
These parameters can be fruitfully combined to compress discretized data on disk. For
example, to save the variable foo
with a precision of 0.1 in 16-bit integers while
converting NaN
to -9999
, we would use
encoding={'foo': {'dtype': 'int16', 'scale_factor': 0.1, '_FillValue': -9999}}
.
Compression and decompression with such discretization is extremely fast.
String encoding#
Xarray can write unicode strings to netCDF files in two ways:
As variable length strings. This is only supported on netCDF4 (HDF5) files.
By encoding strings into bytes, and writing encoded bytes as a character array. The default encoding is UTF-8.
By default, we use variable length strings for compatible files and fall-back
to using encoded character arrays. Character arrays can be selected even for
netCDF4 files by setting the dtype
field in encoding
to S1
(corresponding to NumPy’s single-character bytes dtype).
If character arrays are used:
The string encoding that was used is stored on disk in the
_Encoding
attribute, which matches an ad-hoc convention adopted by the netCDF4-Python library. At the time of this writing (October 2017), a standard convention for indicating string encoding for character arrays in netCDF files was still under discussion. Technically, you can use any string encoding recognized by Python if you feel the need to deviate from UTF-8, by setting the_Encoding
field inencoding
. But we don’t recommend it.The character dimension name can be specified by the
char_dim_name
field of a variable’sencoding
. If the name of the character dimension is not specified, the default isf'string{data.shape[-1]}'
. When decoding character arrays from existing files, thechar_dim_name
is added to the variablesencoding
to preserve if encoding happens, but the field can be edited by the user.
Warning
Missing values in bytes or unicode string arrays (represented by NaN
in
xarray) are currently written to disk as empty strings ''
. This means
missing values will not be restored when data is loaded from disk.
This behavior is likely to change in the future (GH1647).
Unfortunately, explicitly setting a _FillValue
for string arrays to handle
missing values doesn’t work yet either, though we also hope to fix this in the
future.
Chunk based compression#
zlib
, complevel
, fletcher32
, contiguous
and chunksizes
can be used for enabling netCDF4/HDF5’s chunk based compression, as described
in the documentation for createVariable for netCDF4-Python. This only works
for netCDF4 files and thus requires using format='netCDF4'
and either
engine='netcdf4'
or engine='h5netcdf'
.
Chunk based gzip compression can yield impressive space savings, especially for sparse data, but it comes with significant performance overhead. HDF5 libraries can only read complete chunks back into memory, and maximum decompression speed is in the range of 50-100 MB/s. Worse, HDF5’s compression and decompression currently cannot be parallelized with dask. For these reasons, we recommend trying discretization based compression (described above) first.
Time units#
The units
and calendar
attributes control how xarray serializes datetime64
and
timedelta64
arrays to datasets on disk as numeric values. The units
encoding
should be a string like 'days since 1900-01-01'
for datetime64
data or a string
like 'days'
for timedelta64
data. calendar
should be one of the calendar types
supported by netCDF4-python: 'standard'
, 'gregorian'
, 'proleptic_gregorian'
, 'noleap'
,
'365_day'
, '360_day'
, 'julian'
, 'all_leap'
, '366_day'
.
By default, xarray uses the 'proleptic_gregorian'
calendar and units of the smallest time
difference between values, with a reference time of the first time value.
Coordinates#
You can control the coordinates
attribute written to disk by specifying DataArray.encoding["coordinates"]
.
If not specified, xarray automatically sets DataArray.encoding["coordinates"]
to a space-delimited list
of names of coordinate variables that share dimensions with the DataArray
being written.
This allows perfect roundtripping of xarray datasets but may not be desirable.
When an xarray Dataset
contains non-dimensional coordinates that do not share dimensions with any of
the variables, these coordinate variable names are saved under a “global” "coordinates"
attribute.
This is not CF-compliant but again facilitates roundtripping of xarray datasets.
Invalid netCDF files#
The library h5netcdf
allows writing some dtypes that aren’t
allowed in netCDF4 (see
h5netcdf documentation).
This feature is available through DataArray.to_netcdf()
and
Dataset.to_netcdf()
when used with engine="h5netcdf"
and currently raises a warning unless invalid_netcdf=True
is set.
Warning
Note that this produces a file that is likely to be not readable by other netCDF libraries!
HDF5#
HDF5 is both a file format and a data model for storing information. HDF5 stores data hierarchically, using groups to create a nested structure. HDF5 is a more general version of the netCDF4 data model, so the nested structure is one of many similarities between the two data formats.
Reading HDF5 files in xarray requires the h5netcdf
engine, which can be installed
with conda install h5netcdf
. Once installed we can use xarray to open HDF5 files:
xr.open_dataset("/path/to/my/file.h5")
The similarities between HDF5 and netCDF4 mean that HDF5 data can be written with the
same Dataset.to_netcdf()
method as used for netCDF4 data:
ds = xr.Dataset(
{"foo": (("x", "y"), np.random.rand(4, 5))},
coords={
"x": [10, 20, 30, 40],
"y": pd.date_range("2000-01-01", periods=5),
"z": ("x", list("abcd")),
},
)
ds.to_netcdf("saved_on_disk.h5")
Groups#
If you have multiple or highly nested groups, xarray by default may not read the group
that you want. A particular group of an HDF5 file can be specified using the group
argument:
xr.open_dataset("/path/to/my/file.h5", group="/my/group")
While xarray cannot interrogate an HDF5 file to determine which groups are available, the HDF5 Python reader h5py can be used instead.
Natively the xarray data structures can only handle one level of nesting, organized as DataArrays inside of Datasets. If your HDF5 file has additional levels of hierarchy you can only access one group and a time and will need to specify group names.
Complex Data Types#
Xarray leverages NumPy to seamlessly handle complex numbers in DataArray
and Dataset
objects.
In the examples below, we are using a DataArray named da
with complex elements (of \(\mathbb{C}\)):
data = np.array([[1 + 2j, 3 + 4j], [5 + 6j, 7 + 8j]])
da = xr.DataArray(
data,
dims=["x", "y"],
coords={"x": ["a", "b"], "y": [1, 2]},
name="complex_nums",
)
You can access real and imaginary components using the .real
and .imag
attributes. Most NumPy universal functions (ufuncs) like numpy.abs or numpy.angle work directly.
da.real
<xarray.DataArray 'complex_nums' (x: 2, y: 2)> Size: 32B array([[1., 3.], [5., 7.]]) Coordinates: * x (x) <U1 8B 'a' 'b' * y (y) int64 16B 1 2
np.abs(da)
<xarray.DataArray 'complex_nums' (x: 2, y: 2)> Size: 32B array([[ 2.23606798, 5. ], [ 7.81024968, 10.63014581]]) Coordinates: * x (x) <U1 8B 'a' 'b' * y (y) int64 16B 1 2
Note
Like NumPy, .real
and .imag
typically return views, not copies, of the original data.
Writing complex data to NetCDF files is supported via to_netcdf()
using specific backend engines that handle complex types:
This requires the h5netcdf library to be installed.
# write the data to disk
da.to_netcdf("complex_nums_h5.nc", engine="h5netcdf")
# read the file back into memory
ds_h5 = xr.open_dataset("complex_nums_h5.nc", engine="h5netcdf")
# check the dtype
ds_h5[da.name].dtype
dtype('complex128')
Requires the netcdf4-python (>= 1.7.1) library and you have to enable auto_complex=True
.
# write the data to disk
da.to_netcdf("complex_nums_nc4.nc", engine="netcdf4", auto_complex=True)
# read the file back into memory
ds_nc4 = xr.open_dataset(
"complex_nums_nc4.nc", engine="netcdf4", auto_complex=True
)
# check the dtype
ds_nc4[da.name].dtype
dtype('complex128')
Warning
The scipy
engine only supports NetCDF V3 and does not support complex arrays; writing with engine="scipy"
raises a TypeError
.
If direct writing is not supported (e.g., targeting NetCDF3), you can manually split the complex array into separate real and imaginary variables before saving:
# Write data to file
ds_manual = xr.Dataset(
{
f"{da.name}_real": da.real,
f"{da.name}_imag": da.imag,
}
)
ds_manual.to_netcdf("complex_manual.nc", engine="scipy") # Example
# Read data from file
ds = xr.open_dataset("complex_manual.nc", engine="scipy")
reconstructed = ds[f"{da.name}_real"] + 1j * ds[f"{da.name}_imag"]
Recommendations:
Use
engine="netcdf4"
withauto_complex=True
for full compliance and ease.Use
h5netcdf
for HDF5-based storage when interoperability with HDF5 is desired.For maximum legacy support (NetCDF3), manually handle real/imaginary components.