SkyCoord

class astropy.coordinates.SkyCoord(*args, copy=True, **kwargs)[source] [edit on github]

Bases: astropy.utils.misc.ShapedLikeNDArray

High-level object providing a flexible interface for celestial coordinate representation, manipulation, and transformation between systems.

The SkyCoord class accepts a wide variety of inputs for initialization. At a minimum these must provide one or more celestial coordinate values with unambiguous units. Inputs may be scalars or lists/tuples/arrays, yielding scalar or array coordinates (can be checked via SkyCoord.isscalar). Typically one also specifies the coordinate frame, though this is not required. The general pattern for spherical representations is:

SkyCoord(COORD, [FRAME], keyword_args ...)
SkyCoord(LON, LAT, [FRAME], keyword_args ...)
SkyCoord(LON, LAT, [DISTANCE], frame=FRAME, unit=UNIT, keyword_args ...)
SkyCoord([FRAME], <lon_attr>=LON, <lat_attr>=LAT, keyword_args ...)

It is also possible to input coordinate values in other representations such as cartesian or cylindrical. In this case one includes the keyword argument representation_type='cartesian' (for example) along with data in x, y, and z.

See also: http://docs.astropy.org/en/stable/coordinates/

Parameters:
frame : BaseCoordinateFrame class or string, optional

Type of coordinate frame this SkyCoord should represent. Defaults to to ICRS if not given or given as None.

unit : Unit, string, or tuple of Unit or str, optional

Units for supplied LON and LAT values, respectively. If only one unit is supplied then it applies to both LON and LAT.

obstime : valid Time initializer, optional

Time of observation

equinox : valid Time initializer, optional

Coordinate frame equinox

representation_type : str or Representation class

Specifies the representation, e.g. ‘spherical’, ‘cartesian’, or ‘cylindrical’. This affects the positional args and other keyword args which must correspond to the given representation.

copy : bool, optional

If True (default), a copy of any coordinate data is made. This argument can only be passed in as a keyword argument.

**keyword_args

Other keyword arguments as applicable for user-defined coordinate frames. Common options include:

ra, dec : valid Angle initializer, optional

RA and Dec for frames where ra and dec are keys in the frame’s representation_component_names, including ICRS, FK5, FK4, and FK4NoETerms.

pm_ra_cosdec, pm_dec : Quantity, optional

Proper motion components, in angle per time units.

l, b : valid Angle initializer, optional

Galactic l and b for for frames where l and b are keys in the frame’s representation_component_names, including the Galactic frame.

pm_l_cosb, pm_b : Quantity, optional

Proper motion components in the Galactic frame, in angle per time units.

x, y, z : float or Quantity, optional

Cartesian coordinates values

u, v, w : float or Quantity, optional

Cartesian coordinates values for the Galactic frame.

radial_velocity : Quantity, optional

The component of the velocity along the line-of-sight (i.e., the radial direction), in velocity units.

Examples

The examples below illustrate common ways of initializing a SkyCoord object. For a complete description of the allowed syntax see the full coordinates documentation. First some imports:

>>> from astropy.coordinates import SkyCoord  # High-level coordinates
>>> from astropy.coordinates import ICRS, Galactic, FK4, FK5  # Low-level frames
>>> from astropy.coordinates import Angle, Latitude, Longitude  # Angles
>>> import astropy.units as u

The coordinate values and frame specification can now be provided using positional and keyword arguments:

>>> c = SkyCoord(10, 20, unit="deg")  # defaults to ICRS frame
>>> c = SkyCoord([1, 2, 3], [-30, 45, 8], frame="icrs", unit="deg")  # 3 coords

>>> coords = ["1:12:43.2 +1:12:43", "1 12 43.2 +1 12 43"]
>>> c = SkyCoord(coords, frame=FK4, unit=(u.deg, u.hourangle), obstime="J1992.21")

>>> c = SkyCoord("1h12m43.2s +1d12m43s", frame=Galactic)  # Units from string
>>> c = SkyCoord(frame="galactic", l="1h12m43.2s", b="+1d12m43s")

>>> ra = Longitude([1, 2, 3], unit=u.deg)  # Could also use Angle
>>> dec = np.array([4.5, 5.2, 6.3]) * u.deg  # Astropy Quantity
>>> c = SkyCoord(ra, dec, frame='icrs')
>>> c = SkyCoord(frame=ICRS, ra=ra, dec=dec, obstime='2001-01-02T12:34:56')

>>> c = FK4(1 * u.deg, 2 * u.deg)  # Uses defaults for obstime, equinox
>>> c = SkyCoord(c, obstime='J2010.11', equinox='B1965')  # Override defaults

>>> c = SkyCoord(w=0, u=1, v=2, unit='kpc', frame='galactic',
...              representation_type='cartesian')

>>> c = SkyCoord([ICRS(ra=1*u.deg, dec=2*u.deg), ICRS(ra=3*u.deg, dec=4*u.deg)])

Velocity components (proper motions or radial velocities) can also be provided in a similar manner:

>>> c = SkyCoord(ra=1*u.deg, dec=2*u.deg, radial_velocity=10*u.km/u.s)

>>> c = SkyCoord(ra=1*u.deg, dec=2*u.deg, pm_ra_cosdec=2*u.mas/u.yr, pm_dec=1*u.mas/u.yr)

As shown, the frame can be a BaseCoordinateFrame class or the corresponding string alias. The frame classes that are built in to astropy are ICRS, FK5, FK4, FK4NoETerms, and Galactic. The string aliases are simply lower-case versions of the class name, and allow for creating a SkyCoord object and transforming frames without explicitly importing the frame classes.

Attributes Summary

frame
info Container for meta information like name, description, format.
representation
representation_type
shape The shape of the instance and underlying arrays.

Methods Summary

apply_space_motion([new_obstime, dt]) Compute the position of the source represented by this coordinate object to a new time using the velocities stored in this object and assuming linear space motion (including relativistic corrections).
contained_by(wcs[, image]) Determines if the SkyCoord is contained in the given wcs footprint.
directional_offset_by(position_angle, separation) Computes coordinates at the given offset from this coordinate.
from_name(name[, frame, parse]) Given a name, query the CDS name resolver to attempt to retrieve coordinate information for that object.
from_pixel(xp, yp, wcs[, origin, mode]) Create a new SkyCoord from pixel coordinates using an WCS object.
get_constellation([short_name, …]) Determines the constellation(s) of the coordinates this SkyCoord contains.
guess_from_table(table, **coord_kwargs) A convenience method to create and return a new SkyCoord from the data in an astropy Table.
is_equivalent_frame(other) Checks if this object’s frame as the same as that of the other object.
match_to_catalog_3d(catalogcoord[, nthneighbor]) Finds the nearest 3-dimensional matches of this coordinate to a set of catalog coordinates.
match_to_catalog_sky(catalogcoord[, nthneighbor]) Finds the nearest on-sky matches of this coordinate in a set of catalog coordinates.
position_angle(other) Computes the on-sky position angle (East of North) between this SkyCoord and another.
radial_velocity_correction([kind, obstime, …]) Compute the correction required to convert a radial velocity at a given time and place on the Earth’s Surface to a barycentric or heliocentric velocity.
search_around_3d(searcharoundcoords, distlimit) Searches for all coordinates in this object around a supplied set of points within a given 3D radius.
search_around_sky(searcharoundcoords, seplimit) Searches for all coordinates in this object around a supplied set of points within a given on-sky separation.
separation(other) Computes on-sky separation between this coordinate and another.
separation_3d(other) Computes three dimensional separation between this coordinate and another.
skyoffset_frame([rotation]) Returns the sky offset frame with this SkyCoord at the origin.
spherical_offsets_to(tocoord) Computes angular offsets to go from this coordinate to another.
to_pixel(wcs[, origin, mode]) Convert this coordinate to pixel coordinates using a WCS object.
to_string([style]) A string representation of the coordinates.
transform_to(frame[, merge_attributes]) Transform this coordinate to a new frame.

Attributes Documentation

frame
info

Container for meta information like name, description, format. This is required when the object is used as a mixin column within a table, but can be used as a general way to store meta information.

representation
representation_type
shape

The shape of the instance and underlying arrays.

Methods Documentation

apply_space_motion(new_obstime=None, dt=None)[source] [edit on github]

Compute the position of the source represented by this coordinate object to a new time using the velocities stored in this object and assuming linear space motion (including relativistic corrections). This is sometimes referred to as an “epoch transformation.”

The initial time before the evolution is taken from the obstime attribute of this coordinate. Note that this method currently does not support evolving coordinates where the frame has an obstime frame attribute, so the obstime is only used for storing the before and after times, not actually as an attribute of the frame. Alternatively, if dt is given, an obstime need not be provided at all.

Parameters:
new_obstime : Time, optional

The time at which to evolve the position to. Requires that the obstime attribute be present on this frame.

dt : Quantity, TimeDelta, optional

An amount of time to evolve the position of the source. Cannot be given at the same time as new_obstime.

Returns:
new_coord : SkyCoord

A new coordinate object with the evolved location of this coordinate at the new time. obstime will be set on this object to the new time only if self also has obstime.

contained_by(wcs, image=None, **kwargs)[source] [edit on github]

Determines if the SkyCoord is contained in the given wcs footprint.

Parameters:
wcs : WCS

The coordinate to check if it is within the wcs coordinate.

image : array

Optional. The image associated with the wcs object that the cooordinate is being checked against. If not given the naxis keywords will be used to determine if the coordinate falls within the wcs footprint.

**kwargs :

Additional arguments to pass to to_pixel

Returns:
response : bool

True means the WCS footprint contains the coordinate, False means it does not.

directional_offset_by(position_angle, separation)[source] [edit on github]

Computes coordinates at the given offset from this coordinate.

Parameters:
position_angle : Angle

position_angle of offset

separation : Angle

offset angular separation

Returns:
newpoints : SkyCoord

The coordinates for the location that corresponds to offsetting by the given position_angle and separation.

See also

position_angle
inverse operation for the position_angle component
separation
inverse operation for the separation component

Notes

Returned SkyCoord frame retains only the frame attributes that are for the resulting frame type. (e.g. if the input frame is ICRS, an equinox value will be retained, but an obstime will not.)

For a more complete set of transform offsets, use WCS. skyoffset_frame() can also be used to create a spherical frame with (lat=0, lon=0) at a reference point, approximating an xy cartesian system for small offsets. This method is distinct in that it is accurate on the sphere.

classmethod from_name(name, frame='icrs', parse=False)[source] [edit on github]

Given a name, query the CDS name resolver to attempt to retrieve coordinate information for that object. The search database, sesame url, and query timeout can be set through configuration items in astropy.coordinates.name_resolve – see docstring for get_icrs_coordinates for more information.

Parameters:
name : str

The name of the object to get coordinates for, e.g. 'M42'.

frame : str or BaseCoordinateFrame class or instance

The frame to transform the object to.

parse: bool

Whether to attempt extracting the coordinates from the name by parsing with a regex. For objects catalog names that have J-coordinates embedded in their names eg: ‘CRTS SSS100805 J194428-420209’, this may be much faster than a sesame query for the same object name. The coordinates extracted in this way may differ from the database coordinates by a few deci-arcseconds, so only use this option if you do not need sub-arcsecond accuracy for coordinates.

Returns:
coord : SkyCoord

Instance of the SkyCoord class.

classmethod from_pixel(xp, yp, wcs, origin=0, mode='all')[source] [edit on github]

Create a new SkyCoord from pixel coordinates using an WCS object.

Parameters:
xp, yp : float or numpy.ndarray

The coordinates to convert.

wcs : WCS

The WCS to use for convert

origin : int

Whether to return 0 or 1-based pixel coordinates.

mode : ‘all’ or ‘wcs’

Whether to do the transformation including distortions ('all') or only including only the core WCS transformation ('wcs').

Returns:
coord : an instance of this class

A new object with sky coordinates corresponding to the input xp and yp.

See also

to_pixel
to do the inverse operation
astropy.wcs.utils.pixel_to_skycoord
the implementation of this method
get_constellation(short_name=False, constellation_list='iau')[source] [edit on github]

Determines the constellation(s) of the coordinates this SkyCoord contains.

Parameters:
short_name : bool

If True, the returned names are the IAU-sanctioned abbreviated names. Otherwise, full names for the constellations are used.

constellation_list : str

The set of constellations to use. Currently only 'iau' is supported, meaning the 88 “modern” constellations endorsed by the IAU.

Returns:
constellation : str or string array

If this is a scalar coordinate, returns the name of the constellation. If it is an array SkyCoord, it returns an array of names.

Notes

To determine which constellation a point on the sky is in, this first precesses to B1875, and then uses the Delporte boundaries of the 88 modern constellations, as tabulated by Roman 1987.

classmethod guess_from_table(table, **coord_kwargs)[source] [edit on github]

A convenience method to create and return a new SkyCoord from the data in an astropy Table.

This method matches table columns that start with the case-insensitive names of the the components of the requested frames, if they are also followed by a non-alphanumeric character. It will also match columns that end with the component name if a non-alphanumeric character is before it.

For example, the first rule means columns with names like 'RA[J2000]' or 'ra' will be interpreted as ra attributes for ICRS frames, but 'RAJ2000' or 'radius' are not. Similarly, the second rule applied to the Galactic frame means that a column named 'gal_l' will be used as the the l component, but gall or 'fill' will not.

The definition of alphanumeric here is based on Unicode’s definition of alphanumeric, except without _ (which is normally considered alphanumeric). So for ASCII, this means the non-alphanumeric characters are <space>_!"#$%&'()*+,-./\:;<=>?@[]^`{|}~).

Parameters:
table : astropy.Table

The table to load data from.

coord_kwargs

Any additional keyword arguments are passed directly to this class’s constructor.

Returns:
newsc : same as this class

The new SkyCoord (or subclass) object.

is_equivalent_frame(other)[source] [edit on github]

Checks if this object’s frame as the same as that of the other object.

To be the same frame, two objects must be the same frame class and have the same frame attributes. For two SkyCoord objects, all of the frame attributes have to match, not just those relevant for the object’s frame.

Parameters:
other : SkyCoord or BaseCoordinateFrame

The other object to check.

Returns:
isequiv : bool

True if the frames are the same, False if not.

Raises:
TypeError

If other isn’t a SkyCoord or a BaseCoordinateFrame or subclass.

match_to_catalog_3d(catalogcoord, nthneighbor=1)[source] [edit on github]

Finds the nearest 3-dimensional matches of this coordinate to a set of catalog coordinates.

This finds the 3-dimensional closest neighbor, which is only different from the on-sky distance if distance is set in this object or the catalogcoord object.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
catalogcoord : SkyCoord or BaseCoordinateFrame

The base catalog in which to search for matches. Typically this will be a coordinate object that is an array (i.e., catalogcoord.isscalar == False)

nthneighbor : int, optional

Which closest neighbor to search for. Typically 1 is desired here, as that is correct for matching one set of coordinates to another. The next likely use case is 2, for matching a coordinate catalog against itself (1 is inappropriate because each point will find itself as the closest match).

Returns:
idx : integer array

Indices into catalogcoord to get the matched points for each of this object’s coordinates. Shape matches this object.

sep2d : Angle

The on-sky separation between the closest match for each element in this object in catalogcoord. Shape matches this object.

dist3d : Quantity

The 3D distance between the closest match for each element in this object in catalogcoord. Shape matches this object.

Notes

This method requires SciPy to be installed or it will fail.

match_to_catalog_sky(catalogcoord, nthneighbor=1)[source] [edit on github]

Finds the nearest on-sky matches of this coordinate in a set of catalog coordinates.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
catalogcoord : SkyCoord or BaseCoordinateFrame

The base catalog in which to search for matches. Typically this will be a coordinate object that is an array (i.e., catalogcoord.isscalar == False)

nthneighbor : int, optional

Which closest neighbor to search for. Typically 1 is desired here, as that is correct for matching one set of coordinates to another. The next likely use case is 2, for matching a coordinate catalog against itself (1 is inappropriate because each point will find itself as the closest match).

Returns:
idx : integer array

Indices into catalogcoord to get the matched points for each of this object’s coordinates. Shape matches this object.

sep2d : Angle

The on-sky separation between the closest match for each element in this object in catalogcoord. Shape matches this object.

dist3d : Quantity

The 3D distance between the closest match for each element in this object in catalogcoord. Shape matches this object. Unless both this and catalogcoord have associated distances, this quantity assumes that all sources are at a distance of 1 (dimensionless).

Notes

This method requires SciPy to be installed or it will fail.

position_angle(other)[source] [edit on github]

Computes the on-sky position angle (East of North) between this SkyCoord and another.

Parameters:
other : SkyCoord

The other coordinate to compute the position angle to. It is treated as the “head” of the vector of the position angle.

Returns:
pa : Angle

The (positive) position angle of the vector pointing from self to other. If either self or other contain arrays, this will be an array following the appropriate numpy broadcasting rules.

Examples

>>> c1 = SkyCoord(0*u.deg, 0*u.deg)
>>> c2 = SkyCoord(1*u.deg, 0*u.deg)
>>> c1.position_angle(c2).degree
90.0
>>> c3 = SkyCoord(1*u.deg, 1*u.deg)
>>> c1.position_angle(c3).degree  # doctest: +FLOAT_CMP
44.995636455344844
radial_velocity_correction(kind='barycentric', obstime=None, location=None)[source] [edit on github]

Compute the correction required to convert a radial velocity at a given time and place on the Earth’s Surface to a barycentric or heliocentric velocity.

Parameters:
kind : str

The kind of velocity correction. Must be ‘barycentric’ or ‘heliocentric’.

obstime : Time or None, optional

The time at which to compute the correction. If None, the obstime frame attribute on the SkyCoord will be used.

location : EarthLocation or None, optional

The observer location at which to compute the correction. If None, the location frame attribute on the passed-in obstime will be used, and if that is None, the location frame attribute on the SkyCoord will be used.

Returns:
vcorr : Quantity with velocity units

The correction with a positive sign. I.e., add this to an observed radial velocity to get the barycentric (or heliocentric) velocity. If m/s precision or better is needed, see the notes below.

Raises:
ValueError

If either obstime or location are passed in (not None) when the frame attribute is already set on this SkyCoord.

TypeError

If obstime or location aren’t provided, either as arguments or as frame attributes.

Notes

The barycentric correction is calculated to higher precision than the heliocentric correction and includes additional physics (e.g time dilation). Use barycentric corrections if m/s precision is required.

The algorithm here is sufficient to perform corrections at the mm/s level, but care is needed in application. Strictly speaking, the barycentric correction is multiplicative and should be applied as:

sc = SkyCoord(1*u.deg, 2*u.deg)
vcorr = sc.rv_correction(kind='barycentric', obstime=t, location=loc)
rv = rv + vcorr + rv * vcorr / consts.c

If your target is nearby and/or has finite proper motion you may need to account for terms arising from this. See Wright & Eastmann (2014) for details.

The default is for this method to use the builtin ephemeris for computing the sun and earth location. Other ephemerides can be chosen by setting the solar_system_ephemeris variable, either directly or via with statement. For example, to use the JPL ephemeris, do:

sc = SkyCoord(1*u.deg, 2*u.deg)
with coord.solar_system_ephemeris.set('jpl'):
    rv += sc.rv_correction(obstime=t, location=loc)
search_around_3d(searcharoundcoords, distlimit)[source] [edit on github]

Searches for all coordinates in this object around a supplied set of points within a given 3D radius.

This is intended for use on SkyCoord objects with coordinate arrays, rather than a scalar coordinate. For a scalar coordinate, it is better to use separation_3d.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
searcharoundcoords : SkyCoord or BaseCoordinateFrame

The coordinates to search around to try to find matching points in this SkyCoord. This should be an object with array coordinates, not a scalar coordinate object.

distlimit : Quantity with distance units

The physical radius to search within.

Returns:
idxsearcharound : integer array

Indices into self that matches to the corresponding element of idxself. Shape matches idxself.

idxself : integer array

Indices into searcharoundcoords that matches to the corresponding element of idxsearcharound. Shape matches idxsearcharound.

sep2d : Angle

The on-sky separation between the coordinates. Shape matches idxsearcharound and idxself.

dist3d : Quantity

The 3D distance between the coordinates. Shape matches idxsearcharound and idxself.

Notes

This method requires SciPy (>=0.12.0) to be installed or it will fail.

In the current implementation, the return values are always sorted in the same order as the searcharoundcoords (so idxsearcharound is in ascending order). This is considered an implementation detail, though, so it could change in a future release.

search_around_sky(searcharoundcoords, seplimit)[source] [edit on github]

Searches for all coordinates in this object around a supplied set of points within a given on-sky separation.

This is intended for use on SkyCoord objects with coordinate arrays, rather than a scalar coordinate. For a scalar coordinate, it is better to use separation.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
searcharoundcoords : SkyCoord or BaseCoordinateFrame

The coordinates to search around to try to find matching points in this SkyCoord. This should be an object with array coordinates, not a scalar coordinate object.

seplimit : Quantity with angle units

The on-sky separation to search within.

Returns:
idxsearcharound : integer array

Indices into self that matches to the corresponding element of idxself. Shape matches idxself.

idxself : integer array

Indices into searcharoundcoords that matches to the corresponding element of idxsearcharound. Shape matches idxsearcharound.

sep2d : Angle

The on-sky separation between the coordinates. Shape matches idxsearcharound and idxself.

dist3d : Quantity

The 3D distance between the coordinates. Shape matches idxsearcharound and idxself.

Notes

This method requires SciPy (>=0.12.0) to be installed or it will fail.

In the current implementation, the return values are always sorted in the same order as the searcharoundcoords (so idxsearcharound is in ascending order). This is considered an implementation detail, though, so it could change in a future release.

separation(other)[source] [edit on github]

Computes on-sky separation between this coordinate and another.

Note

If the other coordinate object is in a different frame, it is first transformed to the frame of this object. This can lead to unintuitive behavior if not accounted for. Particularly of note is that self.separation(other) and other.separation(self) may not give the same answer in this case.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
other : SkyCoord or BaseCoordinateFrame

The coordinate to get the separation to.

Returns:
sep : Angle

The on-sky separation between this and the other coordinate.

Notes

The separation is calculated using the Vincenty formula, which is stable at all locations, including poles and antipodes [1].

[1]https://en.wikipedia.org/wiki/Great-circle_distance
separation_3d(other)[source] [edit on github]

Computes three dimensional separation between this coordinate and another.

For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.

Parameters:
other : SkyCoord or BaseCoordinateFrame

The coordinate to get the separation to.

Returns:
sep : Distance

The real-space distance between these two coordinates.

Raises:
ValueError

If this or the other coordinate do not have distances.

skyoffset_frame(rotation=None)[source] [edit on github]

Returns the sky offset frame with this SkyCoord at the origin.

Returns:
astrframe : SkyOffsetFrame

A sky offset frame of the same type as this SkyCoord (e.g., if this object has an ICRS coordinate, the resulting frame is SkyOffsetICRS, with the origin set to this object)

rotation : Angle or Quantity with angle units

The final rotation of the frame about the origin. The sign of the rotation is the left-hand rule. That is, an object at a particular position angle in the un-rotated system will be sent to the positive latitude (z) direction in the final frame.

spherical_offsets_to(tocoord)[source] [edit on github]

Computes angular offsets to go from this coordinate to another.

Parameters:
tocoord : BaseCoordinateFrame

The coordinate to find the offset to.

Returns:
lon_offset : Angle

The angular offset in the longitude direction (i.e., RA for equatorial coordinates).

lat_offset : Angle

The angular offset in the latitude direction (i.e., Dec for equatorial coordinates).

Raises:
ValueError

If the tocoord is not in the same frame as this one. This is different from the behavior of the separation/separation_3d methods because the offset components depend critically on the specific choice of frame.

See also

separation
for the total angular offset (not broken out into components).
position_angle
for the direction of the offset.

Notes

This uses the sky offset frame machinery, and hence will produce a new sky offset frame if one does not already exist for this object’s frame class.

to_pixel(wcs, origin=0, mode='all')[source] [edit on github]

Convert this coordinate to pixel coordinates using a WCS object.

Parameters:
wcs : WCS

The WCS to use for convert

origin : int

Whether to return 0 or 1-based pixel coordinates.

mode : ‘all’ or ‘wcs’

Whether to do the transformation including distortions ('all') or only including only the core WCS transformation ('wcs').

Returns:
xp, yp : numpy.ndarray

The pixel coordinates

See also

astropy.wcs.utils.skycoord_to_pixel
the implementation of this method
to_string(style='decimal', **kwargs)[source] [edit on github]

A string representation of the coordinates.

The default styles definitions are:

'decimal': 'lat': {'decimal': True, 'unit': "deg"}
           'lon': {'decimal': True, 'unit': "deg"}
'dms': 'lat': {'unit': "deg"}
       'lon': {'unit': "deg"}
'hmsdms': 'lat': {'alwayssign': True, 'pad': True, 'unit': "deg"}
          'lon': {'pad': True, 'unit': "hour"}

See to_string() for details and keyword arguments (the two angles forming the coordinates are are both Angle instances). Keyword arguments have precedence over the style defaults and are passed to to_string().

Parameters:
style : {‘hmsdms’, ‘dms’, ‘decimal’}

The formatting specification to use. These encode the three most common ways to represent coordinates. The default is decimal.

kwargs

Keyword args passed to to_string().

transform_to(frame, merge_attributes=True)[source] [edit on github]

Transform this coordinate to a new frame.

The precise frame transformed to depends on merge_attributes. If False, the destination frame is used exactly as passed in. But this is often not quite what one wants. E.g., suppose one wants to transform an ICRS coordinate that has an obstime attribute to FK4; in this case, one likely would want to use this information. Thus, the default for merge_attributes is True, in which the precedence is as follows: (1) explicitly set (i.e., non-default) values in the destination frame; (2) explicitly set values in the source; (3) default value in the destination frame.

Note that in either case, any explicitly set attributes on the source SkyCoord that are not part of the destination frame’s definition are kept (stored on the resulting SkyCoord), and thus one can round-trip (e.g., from FK4 to ICRS to FK4 without loosing obstime).

Parameters:
frame : str, BaseCoordinateFrame class or instance, or SkyCoord instance

The frame to transform this coordinate into. If a SkyCoord, the underlying frame is extracted, and all other information ignored.

merge_attributes : bool, optional

Whether the default attributes in the destination frame are allowed to be overridden by explicitly set attributes in the source (see note above; default: True).

Returns:
coord : SkyCoord

A new object with this coordinate represented in the frame frame.

Raises:
ValueError

If there is no possible transformation route.