Note
Click here to download the full example code
Create a new coordinate class (for the Sagittarius stream)¶
This document describes in detail how to subclass and define a custom spherical
coordinate frame, as discussed in Defining a New Frame and
the docstring for BaseCoordinateFrame
. In this example,
we will define a coordinate system defined by the plane of orbit of the
Sagittarius Dwarf Galaxy (hereafter Sgr; as defined in Majewski et al. 2003).
The Sgr coordinate system is often referred to in terms of two angular
coordinates, \(\Lambda,B\).
To do this, we need to define a subclass of
BaseCoordinateFrame
that knows the names and units of the
coordinate system angles in each of the supported representations. In this case
we support SphericalRepresentation
with “Lambda” and
“Beta”. Then we have to define the transformation from this coordinate system to
some other built-in system. Here we will use Galactic coordinates, represented
by the Galactic
class.
See Also¶
The gala package, which defines a number of Astropy coordinate frames for stellar stream coordinate systems.
Majewski et al. 2003, “A Two Micron All Sky Survey View of the Sagittarius Dwarf Galaxy. I. Morphology of the Sagittarius Core and Tidal Arms”, https://arxiv.org/abs/astro-ph/0304198
Law & Majewski 2010, “The Sagittarius Dwarf Galaxy: A Model for Evolution in a Triaxial Milky Way Halo”, https://arxiv.org/abs/1003.1132
David Law’s Sgr info page https://www.stsci.edu/~dlaw/Sgr/
By: Adrian Price-Whelan, Erik Tollerud
License: BSD
Make print
work the same in all versions of Python, set up numpy,
matplotlib, and use a nicer set of plot parameters:
import matplotlib.pyplot as plt
import numpy as np
from astropy.visualization import astropy_mpl_style
plt.style.use(astropy_mpl_style)
Import the packages necessary for coordinates
import astropy.coordinates as coord
import astropy.units as u
from astropy.coordinates import frame_transform_graph
from astropy.coordinates.matrix_utilities import matrix_transpose, rotation_matrix
The first step is to create a new class, which we’ll call
Sagittarius
and make it a subclass of
BaseCoordinateFrame
:
class Sagittarius(coord.BaseCoordinateFrame):
"""
A Heliocentric spherical coordinate system defined by the orbit
of the Sagittarius dwarf galaxy, as described in
https://ui.adsabs.harvard.edu/abs/2003ApJ...599.1082M
and further explained in
https://www.stsci.edu/~dlaw/Sgr/.
Parameters
----------
representation : `~astropy.coordinates.BaseRepresentation` or None
A representation object or None to have no data (or use the other keywords)
Lambda : `~astropy.coordinates.Angle`, optional, must be keyword
The longitude-like angle corresponding to Sagittarius' orbit.
Beta : `~astropy.coordinates.Angle`, optional, must be keyword
The latitude-like angle corresponding to Sagittarius' orbit.
distance : `~astropy.units.Quantity`, optional, must be keyword
The Distance for this object along the line-of-sight.
pm_Lambda_cosBeta : `~astropy.units.Quantity`, optional, must be keyword
The proper motion along the stream in ``Lambda`` (including the
``cos(Beta)`` factor) for this object (``pm_Beta`` must also be given).
pm_Beta : `~astropy.units.Quantity`, optional, must be keyword
The proper motion in Declination for this object (``pm_ra_cosdec`` must
also be given).
radial_velocity : `~astropy.units.Quantity`, optional, keyword-only
The radial velocity of this object.
"""
default_representation = coord.SphericalRepresentation
default_differential = coord.SphericalCosLatDifferential
frame_specific_representation_info = {
coord.SphericalRepresentation: [
coord.RepresentationMapping('lon', 'Lambda'),
coord.RepresentationMapping('lat', 'Beta'),
coord.RepresentationMapping('distance', 'distance')]
}
Breaking this down line-by-line, we define the class as a subclass of
BaseCoordinateFrame
. Then we include a descriptive
docstring. The final lines are class-level attributes that specify the
default representation for the data, default differential for the velocity
information, and mappings from the attribute names used by representation
objects to the names that are to be used by the Sagittarius
frame. In this
case we override the names in the spherical representations but don’t do
anything with other representations like cartesian or cylindrical.
Next we have to define the transformation from this coordinate system to some
other built-in coordinate system; we will use Galactic coordinates. We can do
this by defining functions that return transformation matrices, or by simply
defining a function that accepts a coordinate and returns a new coordinate in
the new system. Because the transformation to the Sagittarius coordinate
system is just a spherical rotation from Galactic coordinates, we’ll just
define a function that returns this matrix. We’ll start by constructing the
transformation matrix using pre-determined Euler angles and the
rotation_matrix
helper function:
SGR_PHI = (180 + 3.75) * u.degree # Euler angles (from Law & Majewski 2010)
SGR_THETA = (90 - 13.46) * u.degree
SGR_PSI = (180 + 14.111534) * u.degree
# Generate the rotation matrix using the x-convention (see Goldstein)
SGR_MATRIX = (
np.diag([1.,1.,-1.])
@ rotation_matrix(SGR_PSI, "z")
@ rotation_matrix(SGR_THETA, "x")
@ rotation_matrix(SGR_PHI, "z")
)
Since we already constructed the transformation (rotation) matrix above, and the inverse of a rotation matrix is just its transpose, the required transformation functions are very simple:
@frame_transform_graph.transform(coord.StaticMatrixTransform, coord.Galactic, Sagittarius)
def galactic_to_sgr():
""" Compute the transformation matrix from Galactic spherical to
heliocentric Sgr coordinates.
"""
return SGR_MATRIX
The decorator @frame_transform_graph.transform(coord.StaticMatrixTransform,
coord.Galactic, Sagittarius)
registers this function on the
frame_transform_graph
as a coordinate transformation. Inside the function,
we simply return the previously defined rotation matrix.
We then register the inverse transformation by using the transpose of the rotation matrix (which is faster to compute than the inverse):
@frame_transform_graph.transform(coord.StaticMatrixTransform, Sagittarius, coord.Galactic)
def sgr_to_galactic():
""" Compute the transformation matrix from heliocentric Sgr coordinates to
spherical Galactic.
"""
return matrix_transpose(SGR_MATRIX)
Now that we’ve registered these transformations between Sagittarius
and
Galactic
, we can transform between any coordinate
system and Sagittarius
(as long as the other system has a path to
transform to Galactic
). For example, to transform from
ICRS coordinates to Sagittarius
, we would do:
icrs = coord.SkyCoord(280.161732*u.degree, 11.91934*u.degree, frame='icrs')
sgr = icrs.transform_to(Sagittarius)
print(sgr)
Out:
<SkyCoord (Sagittarius): (Lambda, Beta) in deg
(346.81830652, -39.28360407)>
Or, to transform from the Sagittarius
frame to ICRS coordinates (in this
case, a line along the Sagittarius
x-y plane):
sgr = coord.SkyCoord(Lambda=np.linspace(0, 2*np.pi, 128)*u.radian,
Beta=np.zeros(128)*u.radian, frame='sagittarius')
icrs = sgr.transform_to(coord.ICRS)
print(icrs)
Out:
<SkyCoord (ICRS): (ra, dec) in deg
[(284.03876751, -29.00408353), (287.24685769, -29.44848352),
(290.48068369, -29.81535572), (293.7357366 , -30.1029631 ),
(297.00711066, -30.30991693), (300.28958688, -30.43520293),
(303.57772919, -30.47820084), (306.86598944, -30.43869669),
(310.14881715, -30.31688708), (313.42076929, -30.11337526),
(316.67661568, -29.82915917), (319.91143548, -29.46561215),
(323.12070147, -29.02445708), (326.30034928, -28.50773532),
(329.44683007, -27.9177717 ), (332.55714589, -27.257137 ),
(335.62886847, -26.52860943), (338.66014233, -25.73513624),
(341.64967439, -24.87979679), (344.59671212, -23.96576781),
(347.50101283, -22.99629167), (350.36280652, -21.97464811),
(353.18275454, -20.90412969), (355.96190618, -19.78802107),
(358.70165491, -18.62958199), ( 1.40369557, -17.43203397),
( 4.06998374, -16.19855028), ( 6.70269788, -14.93224899),
( 9.30420479, -13.63618882), ( 11.87702861, -12.31336727),
( 14.42382347, -10.96672102), ( 16.94734952, -9.59912794),
( 19.45045241, -8.21341071), ( 21.93604568, -6.81234162),
( 24.40709589, -5.39864845), ( 26.86661004, -3.97502106),
( 29.31762493, -2.54411871), ( 31.76319801, -1.10857781),
( 34.20639942, 0.32898001), ( 36.65030466, 1.76593955),
( 39.09798768, 3.19968374), ( 41.55251374, 4.6275852 ),
( 44.01693189, 6.04699804), ( 46.49426651, 7.45524993),
( 48.98750752, 8.84963453), ( 51.4995989 , 10.22740448),
( 54.03342512, 11.58576509), ( 56.59179508, 12.92186896),
( 59.17742314, 14.23281165), ( 61.79290712, 15.51562883),
( 64.44070278, 16.76729487), ( 67.12309478, 17.98472356),
( 69.84216409, 19.16477088), ( 72.59975183, 20.30424045),
( 75.39742013, 21.3998918 ), ( 78.23641033, 22.44845192),
( 81.11759966, 23.44663022), ( 84.04145735, 24.39113719),
( 87.00800203, 25.27870692), ( 90.01676196, 26.10612335),
( 93.06674057, 26.87025019), ( 96.15638947, 27.56806406),
( 99.28359159, 28.19669038), (102.44565666, 28.75344107),
(105.63933131, 29.23585315), (108.86082534, 29.64172698),
(112.105855 , 29.96916281), (115.36970341, 30.21659414),
(118.64729687, 30.38281659), (121.93329519, 30.46701088),
(125.22219273, 30.46875885), (128.50842634, 30.38805179),
(131.78648572, 30.22529063), (135.05102157, 29.98127794),
(138.29694697, 29.6572022 ), (141.51952827, 29.2546151 ),
(144.71446203, 28.77540295), (147.87793614, 28.22175338),
(151.00667382, 27.59611901), (154.09796066, 26.90117914),
(157.14965528, 26.13980125), (160.16018547, 25.31500315),
(163.12853176, 24.42991703), (166.05420084, 23.48775622),
(168.93719133, 22.49178507), (171.77795423, 21.44529257),
(174.57735037, 20.35156967), (177.33660656, 19.21389046),
(180.05727218, 18.03549704), (182.74117737, 16.81958784),
(185.39039367, 15.56930924), (188.00719783, 14.28774998),
(190.59403895, 12.97793826), (193.15350938, 11.64284103),
(195.68831902, 10.28536518), (198.20127316, 8.90836046),
(200.69525342, 7.51462369), (203.17320154, 6.10690412),
(205.63810576, 4.6879097 ), (208.09298919, 3.26031403),
(210.54090002, 1.82676397), (212.984903 , 0.38988751),
(215.42807182, -1.04769799), (217.87348209, -2.48337744),
(220.32420429, -3.91452965), (222.7832966 , -5.338519 ),
(225.25379684, -6.75268736), (227.73871349, -8.15434631),
(230.24101506, -9.54076983), (232.76361762, -10.90918763),
(235.30937003, -12.25677927), (237.88103647, -13.58066929),
(240.48127601, -14.87792359), (243.11261883, -16.14554723),
(245.777439 , -17.38048408), (248.47792364, -18.57961852),
(251.2160385 , -19.7397795 ), (253.9934903 , -20.85774736),
(256.81168612, -21.93026371), (259.67169071, -22.95404466),
(262.57418275, -23.92579758), (265.51941137, -24.84224172),
(268.50715471, -25.70013256), (271.53668252, -26.49628998),
(274.6067251 , -27.22762983), (277.71545113, -27.89119849),
(280.86045662, -28.48420985), (284.03876751, -29.00408353)]>
As an example, we’ll now plot the points in both coordinate systems:
fig, axes = plt.subplots(2, 1, figsize=(8, 10),
subplot_kw={'projection': 'aitoff'})
axes[0].set_title("Sagittarius")
axes[0].plot(sgr.Lambda.wrap_at(180*u.deg).radian, sgr.Beta.radian,
linestyle='none', marker='.')
axes[1].set_title("ICRS")
axes[1].plot(icrs.ra.wrap_at(180*u.deg).radian, icrs.dec.radian,
linestyle='none', marker='.')
plt.show()
This particular transformation is just a spherical rotation, which is a special case of an Affine transformation with no vector offset. The transformation of velocity components is therefore natively supported as well:
sgr = coord.SkyCoord(Lambda=np.linspace(0, 2*np.pi, 128)*u.radian,
Beta=np.zeros(128)*u.radian,
pm_Lambda_cosBeta=np.random.uniform(-5, 5, 128)*u.mas/u.yr,
pm_Beta=np.zeros(128)*u.mas/u.yr,
frame='sagittarius')
icrs = sgr.transform_to(coord.ICRS)
print(icrs)
fig, axes = plt.subplots(3, 1, figsize=(8, 10), sharex=True)
axes[0].set_title("Sagittarius")
axes[0].plot(sgr.Lambda.degree,
sgr.pm_Lambda_cosBeta.value,
linestyle='none', marker='.')
axes[0].set_xlabel(r"$\Lambda$ [deg]")
axes[0].set_ylabel(
fr"$\mu_\Lambda \, \cos B$ [{sgr.pm_Lambda_cosBeta.unit.to_string('latex_inline')}]")
axes[1].set_title("ICRS")
axes[1].plot(icrs.ra.degree, icrs.pm_ra_cosdec.value,
linestyle='none', marker='.')
axes[1].set_ylabel(
fr"$\mu_\alpha \, \cos\delta$ [{icrs.pm_ra_cosdec.unit.to_string('latex_inline')}]")
axes[2].set_title("ICRS")
axes[2].plot(icrs.ra.degree, icrs.pm_dec.value,
linestyle='none', marker='.')
axes[2].set_xlabel("RA [deg]")
axes[2].set_ylabel(
fr"$\mu_\delta$ [{icrs.pm_dec.unit.to_string('latex_inline')}]")
plt.show()
Out:
<SkyCoord (ICRS): (ra, dec) in deg
[(284.03876751, -29.00408353), (287.24685769, -29.44848352),
(290.48068369, -29.81535572), (293.7357366 , -30.1029631 ),
(297.00711066, -30.30991693), (300.28958688, -30.43520293),
(303.57772919, -30.47820084), (306.86598944, -30.43869669),
(310.14881715, -30.31688708), (313.42076929, -30.11337526),
(316.67661568, -29.82915917), (319.91143548, -29.46561215),
(323.12070147, -29.02445708), (326.30034928, -28.50773532),
(329.44683007, -27.9177717 ), (332.55714589, -27.257137 ),
(335.62886847, -26.52860943), (338.66014233, -25.73513624),
(341.64967439, -24.87979679), (344.59671212, -23.96576781),
(347.50101283, -22.99629167), (350.36280652, -21.97464811),
(353.18275454, -20.90412969), (355.96190618, -19.78802107),
(358.70165491, -18.62958199), ( 1.40369557, -17.43203397),
( 4.06998374, -16.19855028), ( 6.70269788, -14.93224899),
( 9.30420479, -13.63618882), ( 11.87702861, -12.31336727),
( 14.42382347, -10.96672102), ( 16.94734952, -9.59912794),
( 19.45045241, -8.21341071), ( 21.93604568, -6.81234162),
( 24.40709589, -5.39864845), ( 26.86661004, -3.97502106),
( 29.31762493, -2.54411871), ( 31.76319801, -1.10857781),
( 34.20639942, 0.32898001), ( 36.65030466, 1.76593955),
( 39.09798768, 3.19968374), ( 41.55251374, 4.6275852 ),
( 44.01693189, 6.04699804), ( 46.49426651, 7.45524993),
( 48.98750752, 8.84963453), ( 51.4995989 , 10.22740448),
( 54.03342512, 11.58576509), ( 56.59179508, 12.92186896),
( 59.17742314, 14.23281165), ( 61.79290712, 15.51562883),
( 64.44070278, 16.76729487), ( 67.12309478, 17.98472356),
( 69.84216409, 19.16477088), ( 72.59975183, 20.30424045),
( 75.39742013, 21.3998918 ), ( 78.23641033, 22.44845192),
( 81.11759966, 23.44663022), ( 84.04145735, 24.39113719),
( 87.00800203, 25.27870692), ( 90.01676196, 26.10612335),
( 93.06674057, 26.87025019), ( 96.15638947, 27.56806406),
( 99.28359159, 28.19669038), (102.44565666, 28.75344107),
(105.63933131, 29.23585315), (108.86082534, 29.64172698),
(112.105855 , 29.96916281), (115.36970341, 30.21659414),
(118.64729687, 30.38281659), (121.93329519, 30.46701088),
(125.22219273, 30.46875885), (128.50842634, 30.38805179),
(131.78648572, 30.22529063), (135.05102157, 29.98127794),
(138.29694697, 29.6572022 ), (141.51952827, 29.2546151 ),
(144.71446203, 28.77540295), (147.87793614, 28.22175338),
(151.00667382, 27.59611901), (154.09796066, 26.90117914),
(157.14965528, 26.13980125), (160.16018547, 25.31500315),
(163.12853176, 24.42991703), (166.05420084, 23.48775622),
(168.93719133, 22.49178507), (171.77795423, 21.44529257),
(174.57735037, 20.35156967), (177.33660656, 19.21389046),
(180.05727218, 18.03549704), (182.74117737, 16.81958784),
(185.39039367, 15.56930924), (188.00719783, 14.28774998),
(190.59403895, 12.97793826), (193.15350938, 11.64284103),
(195.68831902, 10.28536518), (198.20127316, 8.90836046),
(200.69525342, 7.51462369), (203.17320154, 6.10690412),
(205.63810576, 4.6879097 ), (208.09298919, 3.26031403),
(210.54090002, 1.82676397), (212.984903 , 0.38988751),
(215.42807182, -1.04769799), (217.87348209, -2.48337744),
(220.32420429, -3.91452965), (222.7832966 , -5.338519 ),
(225.25379684, -6.75268736), (227.73871349, -8.15434631),
(230.24101506, -9.54076983), (232.76361762, -10.90918763),
(235.30937003, -12.25677927), (237.88103647, -13.58066929),
(240.48127601, -14.87792359), (243.11261883, -16.14554723),
(245.777439 , -17.38048408), (248.47792364, -18.57961852),
(251.2160385 , -19.7397795 ), (253.9934903 , -20.85774736),
(256.81168612, -21.93026371), (259.67169071, -22.95404466),
(262.57418275, -23.92579758), (265.51941137, -24.84224172),
(268.50715471, -25.70013256), (271.53668252, -26.49628998),
(274.6067251 , -27.22762983), (277.71545113, -27.89119849),
(280.86045662, -28.48420985), (284.03876751, -29.00408353)]
(pm_ra_cosdec, pm_dec) in mas / yr
[(-0.38551314, 6.65925646e-02), (-1.16934017, 1.69206826e-01),
(-3.28764514, 3.82399377e-01), (-3.74713923, 3.28401789e-01),
( 2.97159267, -1.74589389e-01), (-2.59736063, 7.71965518e-02),
(-0.35400241, 2.18339082e-04), (-2.39376096, -6.81961345e-02),
(-1.60008466, -9.20453304e-02), ( 2.65113509, 2.29112911e-01),
(-0.7927522 , -9.12494220e-02), (-0.49060679, -7.04054246e-02),
( 0.68295777, 1.17166783e-01), ( 2.00555891, 3.99485453e-01),
( 2.63426559, 5.96220854e-01), (-4.82231071, -1.21967207e+00),
(-3.84942946, -1.07358160e+00), ( 2.04576608, 6.22285508e-01),
( 0.7302608 , 2.40055860e-01), ( 2.8276443 , 9.96647968e-01),
( 3.31649584, 1.24482314e+00), ( 0.88930036, 3.53320310e-01),
(-1.91632309, -8.01574417e-01), (-2.23921766, -9.81337992e-01),
( 0.07829119, 3.57900802e-02), (-0.9117652 , -4.33008289e-01),
( 2.58842827, 1.27226779e+00), ( 3.56404742, 1.80671481e+00),
( 1.18468662, 6.17332798e-01), (-3.77850766, -2.01766444e+00),
(-3.7994815 , -2.07287725e+00), (-1.15948522, -6.44461983e-01),
(-1.64494755, -9.28911803e-01), (-0.28592231, -1.63608406e-01),
( 0.93712981, 5.41960227e-01), ( 2.28036372, 1.32947598e+00),
( 3.82531264, 2.24268514e+00), (-2.4940563 , -1.46676643e+00),
(-3.99576513, -2.35148930e+00), (-3.28274284, -1.92843358e+00),
(-3.33252966, -1.94939220e+00), (-1.95927958, -1.13841995e+00),
( 3.11476025, 1.79316608e+00), (-0.21826068, -1.24179458e-01),
(-1.40384962, -7.87290828e-01), ( 3.30926272, 1.82437818e+00),
(-4.31781868, -2.33348427e+00), ( 2.32056225, 1.22582674e+00),
(-2.28077456, -1.17406786e+00), (-3.69434825, -1.84728141e+00),
( 2.43883048, 1.18055798e+00), ( 4.34451499, 2.02854476e+00),
( 2.70813081, 1.21495379e+00), ( 4.24426996, 1.82183158e+00),
( 1.51189179, 6.18072943e-01), (-0.95131627, -3.68511305e-01),
( 3.39872207, 1.24047990e+00), ( 2.39097685, 8.17009539e-01),
(-3.82902353, -1.21606051e+00), (-3.42646498, -1.00291929e+00),
(-4.42368602, -1.18151959e+00), (-2.7750979 , -6.68319077e-01),
(-0.41177562, -8.81099704e-02), (-1.55033935, -2.89294135e-01),
( 1.46481801, 2.32567283e-01), ( 2.67441128, 3.49140076e-01),
(-2.14328437, -2.18633388e-01), (-0.10320695, -7.55672332e-03),
(-2.96029729, -1.30994678e-01), ( 4.64608681, 7.04978861e-02),
( 0.73048056, -1.01832180e-02), (-3.17180752, 1.36452313e-01),
( 3.2663708 , -2.35162408e-01), (-3.807436 , 3.83764717e-01),
(-1.03587135, 1.33985075e-01), (-1.43876966, 2.26722132e-01),
( 1.35073329, -2.50466761e-01), (-1.11251726, 2.36772693e-01),
(-3.35487668, 8.04167752e-01), ( 3.5491918 , -9.44048101e-01),
( 1.74025872, -5.07506673e-01), ( 2.70999922, -8.57853019e-01),
( 3.94151721, -1.34287589e+00), (-2.57325396, 9.36704904e-01),
(-1.25501651, 4.84989074e-01), (-0.1164969 , 4.75213438e-02),
( 4.52822272, -1.93989114e+00), ( 1.25944546, -5.64020594e-01),
(-1.99333032, 9.29228678e-01), ( 4.11400302, -1.98856340e+00),
(-1.17605099, 5.87293258e-01), (-4.21821525, 2.16887171e+00),
(-2.06907901, 1.09185416e+00), ( 1.4574322 , -7.86928358e-01),
(-0.03925275, 2.16228268e-02), (-0.64674174, 3.62455508e-01),
( 3.77715937, -2.14781892e+00), (-0.53265295, 3.06510925e-01),
(-3.08520411, 1.79201794e+00), ( 1.40277775, -8.20374822e-01),
( 2.84978178, -1.67387515e+00), ( 3.70665161, -2.18129078e+00),
(-3.79845032, 2.23406061e+00), (-3.66884328, 2.15134241e+00),
(-3.99887122, 2.33205294e+00), ( 1.08137684, -6.25627128e-01),
(-1.84892015, 1.05850683e+00), (-3.73365294, 2.10970892e+00),
(-2.16192775, 1.20251104e+00), ( 4.03591898, -2.20373252e+00),
( 1.45044937, -7.75268026e-01), (-3.59311754, 1.87440716e+00),
(-0.60706961, 3.08120239e-01), (-2.42729342, 1.19471075e+00),
(-1.44928303, 6.89336382e-01), (-1.22191557, 5.59535266e-01),
(-3.04315601, 1.33616690e+00), ( 2.70739408, -1.13481657e+00),
( 3.68138371, -1.46596608e+00), (-3.70761188, 1.39514539e+00),
(-2.40216839, 8.49053406e-01), (-1.49847902, 4.94121467e-01),
(-3.59816497, 1.09829482e+00), (-0.4626614 , 1.29535628e-01),
(-2.71444557, 6.89566944e-01), (-2.75913106, 6.27622406e-01),
( 2.38369255, -4.77571368e-01), (-2.54583374, 4.39760880e-01)]>
Total running time of the script: ( 0 minutes 0.693 seconds)