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#

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 Galactic spherical to heliocentric Sgr transformation matrix."""
    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 heliocentric Sgr to spherical Galactic transformation matrix."""
    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)
<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)
<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()
Sagittarius, ICRS

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()
Sagittarius, ICRS, ICRS
<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.17063273,  2.94746670e-02), (-3.71129155,  5.37034371e-01),
     ( 3.10594613, -3.61265225e-01), ( 0.52629461, -4.61248114e-02),
     ( 3.50567214, -2.05968053e-01), ( 2.70676517, -8.04481807e-02),
     (-3.3101294 ,  2.04159798e-03), (-0.3605259 , -1.02710644e-02),
     ( 2.93358773,  1.68755478e-01), (-3.67329528, -3.17448695e-01),
     ( 0.32158476,  3.70158837e-02), ( 3.61565611,  5.18871340e-01),
     (-1.52909831, -2.62328855e-01), ( 2.02144578,  4.02649945e-01),
     ( 1.06961413,  2.42088821e-01), (-1.2361072 , -3.12639629e-01),
     ( 0.1738598 ,  4.84884022e-02), (-2.58916718, -7.87578421e-01),
     (-1.29501067, -4.25703943e-01), ( 3.54347244,  1.24895292e+00),
     ( 2.71158709,  1.01777494e+00), ( 4.31227909,  1.71327468e+00),
     ( 0.24814692,  1.03796809e-01), (-1.40473037, -6.15623621e-01),
     (-2.31746077, -1.05940534e+00), (-3.80097197, -1.80512743e+00),
     (-0.56867973, -2.79518234e-01), ( 2.1127857 ,  1.07102986e+00),
     (-2.1440554 , -1.11725388e+00), (-3.55049238, -1.89590783e+00),
     (-1.31881351, -7.19503047e-01), (-3.32858557, -1.85008555e+00),
     ( 3.69964023,  2.08920915e+00), ( 0.73327895,  4.19591601e-01),
     ( 0.25972949,  1.50206569e-01), ( 1.59111921,  9.27639195e-01),
     ( 3.50828782,  2.05682142e+00), ( 2.30345642,  1.35467373e+00),
     ( 2.67209445,  1.57251523e+00), ( 1.8138608 ,  1.06554495e+00),
     (-3.05746149, -1.78848868e+00), (-2.84091056, -1.65068288e+00),
     ( 0.01468062,  8.45162734e-03), (-2.65756969, -1.51202479e+00),
     ( 1.56937003,  8.80116082e-01), ( 2.4927322 ,  1.37422944e+00),
     (-1.84307157, -9.96053523e-01), (-0.27658694, -1.46105830e-01),
     ( 0.00969674,  4.99156167e-03), ( 2.26009207,  1.13011167e+00),
     ( 4.06184892,  1.96620807e+00), ( 2.46547758,  1.15118296e+00),
     (-3.7491872 , -1.68200486e+00), ( 0.84150466,  3.61211654e-01),
     (-0.03909946, -1.59841606e-02), ( 4.23464549,  1.64037427e+00),
     ( 3.5391611 ,  1.29173793e+00), (-4.72891063, -1.61589398e+00),
     ( 2.19830678,  6.98160783e-01), ( 0.39572929,  1.15829155e-01),
     ( 3.51519225,  9.38870545e-01), (-3.65133356, -8.79340464e-01),
     (-4.02533094, -8.61322937e-01), (-3.05065807, -5.69254394e-01),
     ( 1.04489239,  1.65896228e-01), ( 0.15876567,  2.07266020e-02),
     (-1.09866212, -1.12072959e-01), (-0.26419866, -1.93443971e-02),
     ( 3.95276976,  1.74912095e-01), (-1.4733558 , -2.23561190e-02),
     (-0.12742795,  1.77640126e-03), ( 1.96813425, -8.46698515e-02),
     (-1.64552293,  1.18469444e-01), ( 2.29567545, -2.31389114e-01),
     (-1.85829319,  2.40361461e-01), ( 4.11586   , -6.48579531e-01),
     ( 0.15177933, -2.81444734e-02), ( 0.81051893, -1.72499572e-01),
     (-0.87601458,  2.09981689e-01), ( 2.169456  , -5.77052730e-01),
     ( 3.45581325, -1.00780893e+00), ( 0.80954826, -2.56263326e-01),
     (-2.32106679,  7.90788032e-01), (-3.81610187,  1.38912109e+00),
     (-0.86258771,  3.33338736e-01), (-4.51940783,  1.84355403e+00),
     ( 1.96531715, -8.41942097e-01), (-3.54615272,  1.58808239e+00),
     (-2.21327173,  1.03175853e+00), ( 0.16369656, -7.91251188e-02),
     (-2.02079127,  1.00913744e+00), (-3.70638273,  1.90570376e+00),
     (-0.39505523,  2.08470870e-01), ( 0.02686432, -1.45051651e-02),
     (-1.22827424,  6.76608897e-01), (-3.80168827,  2.13059210e+00),
     ( 1.59862005, -9.09028733e-01), ( 3.76879055, -2.16872070e+00),
     ( 2.77859752, -1.61392777e+00), (-0.25415393,  1.48634726e-01),
     ( 1.03325071, -6.06900044e-01), ( 2.71006723, -1.59482069e+00),
     ( 1.07751593, -6.33741579e-01), (-4.14452391,  2.43027281e+00),
     (-1.71700443,  1.00131888e+00), (-1.3896258 ,  8.03963582e-01),
     ( 1.91625479, -1.09705592e+00), (-2.65944222,  1.50272376e+00),
     ( 1.19416935, -6.64222859e-01), (-4.13627286,  2.25852874e+00),
     ( 2.80821258, -1.50099511e+00), ( 0.77797062, -4.05840798e-01),
     (-2.48751689,  1.26254763e+00), ( 4.09459944, -2.01535666e+00),
     (-0.48542139,  2.30885628e-01), ( 3.53896432, -1.62055005e+00),
     ( 3.9631996 , -1.74013297e+00), ( 1.32293104, -5.54512574e-01),
     ( 2.65443472, -1.05702408e+00), (-3.28526111,  1.23621809e+00),
     (-1.93249715,  6.83046738e-01), (-2.35252751,  7.75742821e-01),
     ( 2.09950528, -6.40847708e-01), ( 1.27403639, -3.56703853e-01),
     ( 0.07250097, -1.84178565e-02), ( 1.4801909 , -3.36700560e-01),
     (-3.19175598,  6.39466391e-01), ( 2.87176381, -4.96061216e-01)]>

Total running time of the script: (0 minutes 0.453 seconds)

Gallery generated by Sphinx-Gallery