Specialized Quantities

Beyond astropy.units.Quantity, astropy provides further specialized quantities. astropy-xarray captures this using the units.class attribute.

import xarray as xr
from astropy import units as u

Time and TimeDelta

astropy_xarray patches astropy.time.Time and astropy.time.TimeDelta to be recognized xarray, capturing:

  • format

  • timescale

  • precision

from astropy.time import Time, TimeDelta

import astropy_xarray  # noqa: F401

ds = xr.Dataset(
    data_vars={
        "d_time": (
            "time",
            Time([0, 1, 2], format="unix", scale="utc").copy("datetime64"),
        ),
        "s_time": (
            "time",
            Time([0, 1, 2], format="unix", scale="utc").copy("isot"),
        ),
        "interval": ("time", TimeDelta([0.9, 0.9, 0.9], format="sec", scale="tai")),
    },
    coords=xr.Coordinates(
        {
            "f_time": ("time", Time([0, 1, 2], format="unix", scale="utc")),
        }
    ),
)
ds
<xarray.Dataset> Size: 348B
Dimensions:   (time: 3)
Coordinates:
    f_time    (time) float64 24B [utc unix] 1970-01-01T00:00:00.000000000 ......
Dimensions without coordinates: time
Data variables:
    d_time    (time) datetime64[ns] 24B [utc datetime64] 1970-01-01T00:00:00....
    s_time    (time) <U23 276B [utc isot] 1970-01-01T00:00:00.000000000 ... 1...
    interval  (time) float64 24B [tai sec] 0.9 0.9 0.9
# dequantified attributes
ds.astropy.dequantify()
<xarray.Dataset> Size: 348B
Dimensions:   (time: 3)
Coordinates:
    f_time    (time) float64 24B 0.0 1.0 2.0
Dimensions without coordinates: time
Data variables:
    d_time    (time) datetime64[ns] 24B 1970-01-01 ... 1970-01-01T00:00:02
    s_time    (time) <U23 276B '1970-01-01T00:00:00.000' ... '1970-01-01T00:0...
    interval  (time) float64 24B 0.9 0.9 0.9

SkyCoord, Frame and Angle

SkyCoord and Frame fully support conversion to datasets preserving:

  • Representation

  • Differential

  • Frame

  • Angle, Longitude, Latitude containers

Conversion optionally requires coords argument for applying named coordinates to data variable axes.

Frame-specific mapped names for data variables coming soon.

from astropy.coordinates import SkyCoord

from astropy_xarray.coordinates.sky_coord import (
    skycoord_to_dataset,
)

sc = SkyCoord(ra=[2, 6, 7, 4] * u.deg, dec=[4, 7, 4, 3] * u.deg, frame="icrs")


s = skycoord_to_dataset(sc, coords={"field": ["a", "b", "c", "d"]})
display(s)
<xarray.Dataset> Size: 80B
Dimensions:  (field: 4)
Coordinates:
  * field    (field) <U1 16B 'a' 'b' 'c' 'd'
Data variables:
    ra       (field) float64 32B [°] 2.0 6.0 7.0 4.0
    dec      (field) float64 32B [°] 4.0 7.0 4.0 3.0
Attributes:
    frame:    {'name': 'icrs', 'representation_type': 'spherical', 'different...

Logarithmic Units and Quantities

Astropy provides a generic LogQuantity and other various log scale specialized classes and units:

\[scalar = A * log_{B}(\frac{V}{V_0})\]

Name

Equiv SI Unit

Scale Factor

Log Base

0 Reference

Magnitude (absolute)

Unitless

-2.5

10

Apparent brightness at 10pc away

Dex

Unitless

1

10

1

Decibels (power)

Unitless

10

10

Variable

Decibels (amplitude)

Unitless

20

10

Variable

Magnitude (apparent)

Unitless

-2.5

10

Usually Vega or AB from Earth

Magnitude (AB)

Jy

-2.5

10

3631 Jy

Custom Unit References

import types

import astropy.units
import astropy.units as u
import numpy as np
from astropy.units import Dex, Equivalency, Magnitude, Quantity


# extend registry
def solar():
    solar = types.ModuleType("solar")

    # Solar relative metalicity
    astropy.units.def_unit(
        ["solH", "Xsun"],
        format={"latex": "X_{\\odot}", "unicode": "X☉"},
        namespace=solar.__dict__,
    )
    astropy.units.def_unit(
        ["solHe", "Ysun"],
        format={"latex": "Y_{\\odot}", "unicode": "Y☉"},
        namespace=solar.__dict__,
    )
    astropy.units.def_unit(
        ["solMetal", "Zsun"],
        format={"latex": "Z_{\\odot}", "unicode": "Z☉"},
        namespace=solar.__dict__,
    )

    return solar


if not hasattr(u, "solar"):
    u.solar = solar()
    u.add_enabled_units(u.solar)
    u.__dict__.update(u.get_current_unit_registry().registry)

asplund et. al equivalencies

def asplund_solar_relative_mass() -> Equivalency:
    """measurements from https://www.aanda.org/articles/aa/full_html/2021/09/aa40445-21/aa40445-21.html"""
    return Equivalency(
        [
            (u.Xsun, u.Msun, lambda x: x * 0.7438),
            (u.Ysun, u.Msun, lambda x: x * 0.2423),
            (u.Zsun, u.Msun, lambda x: x * 0.0139),
        ]
    )


if not hasattr(u, "asplund_solar_relative_mass"):
    u.asplund_solar_relative_mass = asplund_solar_relative_mass
    # uncomment for implicit conversions
    # u.add_enabled_equivalencies(u.asplund_solar_relative_mass())
ds = xr.Dataset(
    coords=xr.Coordinates(
        {
            "timestamp": ("time", Time([1.7e9], format="unix")),
            "body_label": ("body", ["sun", "proxima centauri"]),
            "wavelength": ("frequency", [565, 535, 445] * u.nm),
        }
    ),
    data_vars={
        "mass": (["time", "body"], Dex([[0.0, 0.02]], unit=u.dex(u.Msun))),
        "distance": (["time", "body"], [[1.007, 12950]] * u.AU),
        # radiant flux
        "luminosity": (
            ["time", "body"],
            Quantity([[1, 0.001]], unit=u.Lsun).to(u.dex(u.Lsun)),
        ),
        "metallicity": (["time", "body"], Dex([[0.0, np.nan]], unit=u.dex(u.Zsun))),
        # relative to Vega
        "apparent_magnitude": (["time", "body"], Magnitude([[-26.74, 10.67]])),
        "absolute_magnitude": (["time", "body"], Magnitude([[4.83, 15.60]])),
    },
)
ds
<xarray.Dataset> Size: 256B
Dimensions:             (time: 1, body: 2, frequency: 3)
Coordinates:
    timestamp           (time) float64 8B [utc unix] 2023-11-14T22:13:20.0000...
    body_label          (body) <U16 128B 'sun' 'proxima centauri'
    wavelength          (frequency) float64 24B [nm] 565.0 535.0 445.0
Dimensions without coordinates: time, body, frequency
Data variables:
    mass                (time, body) float64 16B [dex(M☉)] 0.0 0.02
    distance            (time, body) float64 16B [AU] 1.007 1.295e+04
    luminosity          (time, body) float64 16B [dex(L☉)] 0.0 -3.0
    metallicity         (time, body) float64 16B [dex(Z☉)] 0.0 nan
    apparent_magnitude  (time, body) float64 16B [mag(1)] -26.74 10.67
    absolute_magnitude  (time, body) float64 16B [mag(1)] 4.83 15.6
ds.astropy.sel(body=[0]).astropy.to(
    {
        "mass": u.kg,
        "distance": u.km,
        "illuminance": u.cd * u.sr / u.m**2,
        "luminosity": u.kW,
        "metallicity": u.kg,
    },
    equivalencies=asplund_solar_relative_mass(),
)
<xarray.Dataset> Size: 144B
Dimensions:             (time: 1, body: 1, frequency: 3)
Coordinates:
    timestamp           (time) float64 8B [utc unix] 2023-11-14T22:13:20.0000...
    body_label          (body) <U16 64B 'sun'
    wavelength          (frequency) float64 24B [nm] 565.0 535.0 445.0
Dimensions without coordinates: time, body, frequency
Data variables:
    mass                (time, body) float64 8B [kg] 1.988e+30
    distance            (time, body) float64 8B [km] 1.506e+08
    luminosity          (time, body) float64 8B [kW] 3.828e+23
    metallicity         (time, body) float64 8B [kg] 2.764e+28
    apparent_magnitude  (time, body) float64 8B [mag(1)] -26.74
    absolute_magnitude  (time, body) float64 8B [mag(1)] 4.83