xarray.Dataset.pint.quantify
- Dataset.pint.quantify(units=<object object>, unit_registry=None, **unit_kwargs)
Attach units to the variables of the Dataset.
Units can be specified as a
pint.Unit
or as a string, which will be parsed by the given unit registry. If no units are specified then the units will be parsed from the"units"
entry of the Dataset variable’s.attrs
. Will raise a ValueError if any of the variables already contain a unit-aware array with a different unit.Note
Be aware that unless you’re using
dask
this will load the data into memory. To avoid that, consider converting todask
first (e.g. usingchunk
).Warning
As units in dimension coordinates are not supported until
xarray
changes the way it implements indexes, these units will be set as attributes.- Parameters
units (mapping of hashable to unit-like, optional) – Physical units to use for particular DataArrays in this Dataset. It should map variable names to units (unit names or
pint.Unit
objects). If not provided,quantify
will try to read them fromDataset[var].attrs['units']
using pint’s parser. The"units"
attribute will be removed from all variables except from dimension coordinates.unit_registry (
pint.UnitRegistry
, optional) – Unit registry to be used for the units attached to each DataArray in this Dataset. If not given then a default registry will be created.**unit_kwargs – Keyword argument form of
units
.
- Returns
quantified – Dataset whose variables will now contain Quantity arrays with units.
- Return type
Notes
"none"
andNone
can be used to mark variables that should not be quantified.Examples
>>> ds = xr.Dataset( ... {"a": ("x", [0, 3, 2], {"units": "m"}), "b": ("x", [5, -2, 1])}, ... coords={"x": [0, 1, 2], "u": ("x", [-1, 0, 1], {"units": "s"})}, ... )
>>> ds.pint.quantify() <xarray.Dataset> Dimensions: (x: 3) Coordinates: * x (x) int64 0 1 2 u (x) int64 [s] -1 0 1 Data variables: a (x) int64 [m] 0 3 2 b (x) int64 5 -2 1 >>> ds.pint.quantify({"b": "dm"}) <xarray.Dataset> Dimensions: (x: 3) Coordinates: * x (x) int64 0 1 2 u (x) int64 [s] -1 0 1 Data variables: a (x) int64 [m] 0 3 2 b (x) int64 [dm] 5 -2 1
Don’t quantify specific variables:
>>> ds.pint.quantify({"a": None}) <xarray.Dataset> Dimensions: (x: 3) Coordinates: * x (x) int64 0 1 2 u (x) int64 [s] -1 0 1 Data variables: a (x) int64 0 3 2 b (x) int64 5 -2 1
Quantify with the same unit:
>>> q = ds.pint.quantify() >>> q <xarray.Dataset> Dimensions: (x: 3) Coordinates: * x (x) int64 0 1 2 u (x) int64 [s] -1 0 1 Data variables: a (x) int64 [m] 0 3 2 b (x) int64 5 -2 1 >>> q.pint.quantify({"a": "m"}) <xarray.Dataset> Dimensions: (x: 3) Coordinates: * x (x) int64 0 1 2 u (x) int64 [s] -1 0 1 Data variables: a (x) int64 [m] 0 3 2 b (x) int64 5 -2 1