xarray.Dataset.astropy.to¶
- Dataset.astropy.to(units=None, equivalencies=None, **unit_kwargs)¶
convert the quantities in a Dataset
- Parameters:
units (unit-like or mapping of hashable to unit-like, optional) – The units to convert to. If a unit name or
astropy.units.Unitobject, convert all the object’s data variables. If a dict-like, it maps variable names to unit names orastropy.units.Unitobjects.equivalencies (
list) – A list of equivalence pairs to try if the units are not directly convertible. See Equivalencies. This list is in addition to possible global defaults set by, e.g.,astropy.units.set_enabled_equivalencies(). Use None to turn off all equivalencies.**unit_kwargs – The kwargs form of
units. Can only be used for variable names that are strings and valid python identifiers.
- Returns:
object – A new object with converted units.
- Return type:
Examples
>>> ds = xr.Dataset( ... data_vars={ ... "a": ("x", np.linspace(0, 1, 5) * u.m), ... "b": ("x", np.linspace(-1, 0, 5) * u.kg), ... }, ... coords={"u": ("x", np.arange(5) * u.s)}, ... ) >>> ds <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [s] 0.0 1.0 2.0 3.0 4.0 Dimensions without coordinates: x Data variables: a (x) float64 40B [m] 0.0 0.25 0.5 0.75 1.0 b (x) float64 40B [kg] -1.0 -0.75 -0.5 -0.25 0.0
Convert the data
>>> ds.astropy.to({"a": "mm", "b": u.g}) <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [s] 0.0 1.0 2.0 3.0 4.0 Dimensions without coordinates: x Data variables: a (x) float64 40B [mm] 0.0 250.0 500.0 750.0 1e+03 b (x) float64 40B [g] -1e+03 -750.0 -500.0 -250.0 0.0 >>> ds.astropy.to(a=u.mm, b="g") <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [s] 0.0 1.0 2.0 3.0 4.0 Dimensions without coordinates: x Data variables: a (x) float64 40B [mm] 0.0 250.0 500.0 750.0 1e+03 b (x) float64 40B [g] -1e+03 -750.0 -500.0 -250.0 0.0
Convert coordinates
>>> ds.astropy.to({"u": u.ms}) <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [ms] 0.0 1e+03 2e+03 3e+03 4e+03 Dimensions without coordinates: x Data variables: a (x) float64 40B [m] 0.0 0.25 0.5 0.75 1.0 b (x) float64 40B [kg] -1.0 -0.75 -0.5 -0.25 0.0 >>> ds.astropy.to(u="ms") <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [ms] 0.0 1e+03 2e+03 3e+03 4e+03 Dimensions without coordinates: x Data variables: a (x) float64 40B [m] 0.0 0.25 0.5 0.75 1.0 b (x) float64 40B [kg] -1.0 -0.75 -0.5 -0.25 0.0
Convert both simultaneously
>>> ds.astropy.to(a=u.mm, b=u.g, u="ms") <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [ms] 0.0 1e+03 2e+03 3e+03 4e+03 Dimensions without coordinates: x Data variables: a (x) float64 40B [mm] 0.0 250.0 500.0 750.0 1e+03 b (x) float64 40B [g] -1e+03 -750.0 -500.0 -250.0 0.0 >>> ds.astropy.to({"a": "mm", "b": "g", "u": u.ms}) <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [ms] 0.0 1e+03 2e+03 3e+03 4e+03 Dimensions without coordinates: x Data variables: a (x) float64 40B [mm] 0.0 250.0 500.0 750.0 1e+03 b (x) float64 40B [g] -1e+03 -750.0 -500.0 -250.0 0.0
Convert homogeneous data
>>> ds = xr.Dataset( ... data_vars={ ... "a": ("x", np.linspace(0, 1, 5) * u.kg), ... "b": ("x", np.linspace(-1, 0, 5) * u.mg), ... }, ... coords={"u": ("x", np.arange(5) * u.s)}, ... ) >>> ds <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [s] 0.0 1.0 2.0 3.0 4.0 Dimensions without coordinates: x Data variables: a (x) float64 40B [kg] 0.0 0.25 0.5 0.75 1.0 b (x) float64 40B [mg] -1.0 -0.75 -0.5 -0.25 0.0 >>> ds.astropy.to("g") <xarray.Dataset> Size: 120B Dimensions: (x: 5) Coordinates: u (x) float64 40B [s] 0.0 1.0 2.0 3.0 4.0 Dimensions without coordinates: x Data variables: a (x) float64 40B [g] 0.0 250.0 500.0 750.0 1e+03 b (x) float64 40B [g] -0.001 -0.00075 -0.0005 -0.00025 0.0