Complex Numbers#

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",
)

Operations on Complex Data#

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.

Reading and Writing Complex Data#

Writing complex data to NetCDF files (see netCDF) 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.

Alternative: Manual Handling#

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" with auto_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.

See also#