1# -*- coding: utf-8 -*- 2# Licensed under a 3-clause BSD style license - see LICENSE.rst 3 4from astropy import units as u 5from astropy.utils.decorators import format_doc 6from astropy.coordinates import representation as r 7from astropy.coordinates.baseframe import BaseCoordinateFrame, RepresentationMapping, base_doc 8from .galactic import Galactic 9 10__all__ = ['Supergalactic'] 11 12 13doc_components = """ 14 sgl : `~astropy.coordinates.Angle`, optional, keyword-only 15 The supergalactic longitude for this object (``sgb`` must also be given and 16 ``representation`` must be None). 17 sgb : `~astropy.coordinates.Angle`, optional, keyword-only 18 The supergalactic latitude for this object (``sgl`` must also be given and 19 ``representation`` must be None). 20 distance : `~astropy.units.Quantity` ['speed'], optional, keyword-only 21 The Distance for this object along the line-of-sight. 22 23 pm_sgl_cossgb : `~astropy.units.Quantity` ['angular speed'], optional, keyword-only 24 The proper motion in Right Ascension for this object (``pm_sgb`` must 25 also be given). 26 pm_sgb : `~astropy.units.Quantity` ['angular speed'], optional, keyword-only 27 The proper motion in Declination for this object (``pm_sgl_cossgb`` must 28 also be given). 29 radial_velocity : `~astropy.units.Quantity` ['speed'], optional, keyword-only 30 The radial velocity of this object. 31""" 32 33 34@format_doc(base_doc, components=doc_components, footer="") 35class Supergalactic(BaseCoordinateFrame): 36 """ 37 Supergalactic Coordinates 38 (see Lahav et al. 2000, <https://ui.adsabs.harvard.edu/abs/2000MNRAS.312..166L>, 39 and references therein). 40 """ 41 42 frame_specific_representation_info = { 43 r.SphericalRepresentation: [ 44 RepresentationMapping('lon', 'sgl'), 45 RepresentationMapping('lat', 'sgb') 46 ], 47 r.CartesianRepresentation: [ 48 RepresentationMapping('x', 'sgx'), 49 RepresentationMapping('y', 'sgy'), 50 RepresentationMapping('z', 'sgz') 51 ], 52 r.CartesianDifferential: [ 53 RepresentationMapping('d_x', 'v_x', u.km/u.s), 54 RepresentationMapping('d_y', 'v_y', u.km/u.s), 55 RepresentationMapping('d_z', 'v_z', u.km/u.s) 56 ], 57 } 58 59 default_representation = r.SphericalRepresentation 60 default_differential = r.SphericalCosLatDifferential 61 62 # North supergalactic pole in Galactic coordinates. 63 # Needed for transformations to/from Galactic coordinates. 64 _nsgp_gal = Galactic(l=47.37*u.degree, b=+6.32*u.degree) 65