BaseCoordinateFrame

class astropy.coordinates.BaseCoordinateFrame(*args, copy=True, representation_type=None, differential_type=None, **kwargs)[source]

Bases: ShapedLikeNDArray

The base class for coordinate frames.

This class is intended to be subclassed to create instances of specific systems. Subclasses can implement the following attributes:

  • default_representation

    A subclass of BaseRepresentation that will be treated as the default representation of this frame. This is the representation assumed by default when the frame is created.

  • default_differential

    A subclass of BaseDifferential that will be treated as the default differential class of this frame. This is the differential class assumed by default when the frame is created.

  • Attribute class attributes

    Frame attributes such as FK4.equinox or FK4.obstime are defined using a descriptor class. See the narrative documentation or built-in classes code for details.

  • frame_specific_representation_info

    A dictionary mapping the name or class of a representation to a list of RepresentationMapping objects that tell what names and default units should be used on this frame for the components of that representation.

Unless overridden via frame_specific_representation_info, velocity name defaults are:

where {lon} and {lat} are the frame names of the angular components.

Parameters:
dataBaseRepresentation subclass instance

A representation object or None to have no data (or use the coordinate component arguments, see below).

*args, **kwargs

Coordinate components, with names that depend on the subclass.

representation_typeBaseRepresentation subclass, python:str, optional

A representation class or string name of a representation class. This sets the expected input representation class, thereby changing the expected keyword arguments for the data passed in. For example, passing representation_type='cartesian' will make the classes expect position data with cartesian names, i.e. x, y, z in most cases unless overridden via frame_specific_representation_info. To see this frame’s names, check out <this frame>().representation_info.

differential_typeBaseDifferential subclass, python:str, python:dict, optional

A differential class or dictionary of differential classes (currently only a velocity differential with key ‘s’ is supported). This sets the expected input differential class, thereby changing the expected keyword arguments of the data passed in. For example, passing differential_type='cartesian' will make the classes expect velocity data with the argument names v_x, v_y, v_z unless overridden via frame_specific_representation_info. To see this frame’s names, check out <this frame>().representation_info.

copybool, optional

If True (default), make copies of the input coordinate arrays. Can only be passed in as a keyword argument.

Attributes Summary

cache

Cache for this frame, a dict. It stores anything that should be computed from the coordinate data (not from the frame attributes). This can be used in functions to store anything that might be expensive to compute but might be re-used by some other function. E.g.::.

cartesian

Shorthand for a cartesian representation of the coordinates in this object.

cylindrical

Shorthand for a cylindrical representation of the coordinates in this object.

data

The coordinate data for this object.

default_differential

default_representation

differential_type

The differential used for this frame's data.

frame_attributes

frame_specific_representation_info

has_data

True if this frame has data, False otherwise.

isscalar

proper_motion

Shorthand for the two-dimensional proper motion as a Quantity object with angular velocity units.

radial_velocity

Shorthand for the radial or line-of-sight velocity as a Quantity object.

representation_component_names

representation_component_units

representation_info

A dictionary with the information of what attribute names for this frame apply to particular representations.

representation_type

The representation class used for this frame's data.

shape

The shape of the underlying data.

size

The size of the object, as calculated from its shape.

spherical

Shorthand for a spherical representation of the coordinates in this object.

sphericalcoslat

Shorthand for a spherical representation of the positional data and a SphericalCosLatDifferential for the velocity data in this object.

velocity

Shorthand for retrieving the Cartesian space-motion as a CartesianDifferential object.

Methods Summary

get_frame_attr_defaults()

Return a dict with the defaults for each frame attribute

get_frame_attr_names()

Deprecated since version 5.2.

get_representation_cls([which])

The class used for part of this frame's data.

get_representation_component_names([which])

get_representation_component_units([which])

is_equivalent_frame(other)

Checks if this object is the same frame as the other object.

is_frame_attr_default(attrnm)

Determine whether or not a frame attribute has its value because it's the default value, or because this frame was created with that value explicitly requested.

is_transformable_to(new_frame)

Determines if this coordinate frame can be transformed to another given frame.

realize_frame(data, **kwargs)

Generates a new frame with new data from another frame (which may or may not have data).

replicate([copy])

Return a replica of the frame, optionally with new frame attributes.

replicate_without_data([copy])

Return a replica without data, optionally with new frame attributes.

represent_as(base[, s, in_frame_units])

Generate and return a new representation of this frame's data as a Representation object.

separation(other)

Computes on-sky separation between this coordinate and another.

separation_3d(other)

Computes three dimensional separation between this coordinate and another.

set_representation_cls([base, s])

Set representation and/or differential class for this frame's data.

transform_to(new_frame)

Transform this object's coordinate data to a new frame.

Attributes Documentation

cache

Cache for this frame, a dict. It stores anything that should be computed from the coordinate data (not from the frame attributes). This can be used in functions to store anything that might be expensive to compute but might be re-used by some other function. E.g.:

if 'user_data' in myframe.cache:
    data = myframe.cache['user_data']
else:
    myframe.cache['user_data'] = data = expensive_func(myframe.lat)

If in-place modifications are made to the frame data, the cache should be cleared:

myframe.cache.clear()
cartesian

Shorthand for a cartesian representation of the coordinates in this object.

cylindrical

Shorthand for a cylindrical representation of the coordinates in this object.

data

The coordinate data for this object. If this frame has no data, an ValueError will be raised. Use has_data to check if data is present on this frame object.

default_differential = None
default_representation = None
differential_type

The differential used for this frame’s data.

This will be a subclass from BaseDifferential. For simultaneous setting of representation and differentials, see the set_representation_cls method.

frame_attributes = {}
frame_specific_representation_info = {<class 'astropy.coordinates.representation.SphericalRepresentation'>: [RepresentationMapping(reprname='lon', framename='lon', defaultunit='recommended'), RepresentationMapping(reprname='lat', framename='lat', defaultunit='recommended')], <class 'astropy.coordinates.representation.SphericalCosLatDifferential'>: [RepresentationMapping(reprname='d_lon_coslat', framename='pm_lon_coslat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_lat', framename='pm_lat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_distance', framename='radial_velocity', defaultunit=Unit("km / s"))], <class 'astropy.coordinates.representation.SphericalDifferential'>: [RepresentationMapping(reprname='d_lon', framename='pm_lon', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_lat', framename='pm_lat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_distance', framename='radial_velocity', defaultunit=Unit("km / s"))], <class 'astropy.coordinates.representation.CartesianDifferential'>: [RepresentationMapping(reprname='d_x', framename='v_x', defaultunit=Unit("km / s")), RepresentationMapping(reprname='d_y', framename='v_y', defaultunit=Unit("km / s")), RepresentationMapping(reprname='d_z', framename='v_z', defaultunit=Unit("km / s"))], <class 'astropy.coordinates.representation.UnitSphericalRepresentation'>: [RepresentationMapping(reprname='lon', framename='lon', defaultunit='recommended'), RepresentationMapping(reprname='lat', framename='lat', defaultunit='recommended')], <class 'astropy.coordinates.representation.UnitSphericalCosLatDifferential'>: [RepresentationMapping(reprname='d_lon_coslat', framename='pm_lon_coslat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_lat', framename='pm_lat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_distance', framename='radial_velocity', defaultunit=Unit("km / s"))], <class 'astropy.coordinates.representation.UnitSphericalDifferential'>: [RepresentationMapping(reprname='d_lon', framename='pm_lon', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_lat', framename='pm_lat', defaultunit=Unit("mas / yr")), RepresentationMapping(reprname='d_distance', framename='radial_velocity', defaultunit=Unit("km / s"))]}
has_data

True if this frame has data, False otherwise.

isscalar
proper_motion

Shorthand for the two-dimensional proper motion as a Quantity object with angular velocity units. In the returned Quantity, axis=0 is the longitude/latitude dimension so that .proper_motion[0] is the longitudinal proper motion and .proper_motion[1] is latitudinal. The longitudinal proper motion already includes the cos(latitude) term.

radial_velocity

Shorthand for the radial or line-of-sight velocity as a Quantity object.

representation_component_names
representation_component_units
representation_info

A dictionary with the information of what attribute names for this frame apply to particular representations.

representation_type

The representation class used for this frame’s data.

This will be a subclass from BaseRepresentation. Can also be set using the string name of the representation. If you wish to set an explicit differential class (rather than have it be inferred), use the set_representation_cls method.

shape
size
spherical

Shorthand for a spherical representation of the coordinates in this object.

sphericalcoslat

Shorthand for a spherical representation of the positional data and a SphericalCosLatDifferential for the velocity data in this object.

velocity

Shorthand for retrieving the Cartesian space-motion as a CartesianDifferential object.

This is equivalent to calling self.cartesian.differentials['s'].

Methods Documentation

classmethod get_frame_attr_defaults()[source]

Return a dict with the defaults for each frame attribute

classmethod get_frame_attr_names()[source]

Deprecated since version 5.2: The get_frame_attr_names() method is deprecated and may be removed in a future version. Use get_frame_attr_defaults() to obtain a dict of frame attribute names and default values. The fastest way to obtain the names is frame_attributes.keys()

Return a dict with the defaults for each frame attribute

get_representation_cls(which='base')[source]

The class used for part of this frame’s data.

Parameters:
which(‘base’, ‘s’, None)

The class of which part to return. ‘base’ means the class used to represent the coordinates; ‘s’ the first derivative to time, i.e., the class representing the proper motion and/or radial velocity. If None, return a dict with both.

Returns:
representationBaseRepresentation or BaseDifferential.
get_representation_component_names(which='base')[source]
get_representation_component_units(which='base')[source]
is_equivalent_frame(other)[source]

Checks if this object is the same frame as the other object.

To be the same frame, two objects must be the same frame class and have the same frame attributes. Note that it does not matter what, if any, data either object has.

Parameters:
otherBaseCoordinateFrame

the other frame to check

Returns:
isequivbool

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

Raises:
TypeError

If other isn’t a BaseCoordinateFrame or subclass.

is_frame_attr_default(attrnm)[source]

Determine whether or not a frame attribute has its value because it’s the default value, or because this frame was created with that value explicitly requested.

Parameters:
attrnmpython:str

The name of the attribute to check.

Returns:
isdefaultbool

True if the attribute attrnm has its value by default, False if it was specified at creation of this frame.

is_transformable_to(new_frame)[source]

Determines if this coordinate frame can be transformed to another given frame.

Parameters:
new_frameBaseCoordinateFrame subclass or instance

The proposed frame to transform into.

Returns:
transformablebool or python:str

True if this can be transformed to new_frame, False if not, or the string ‘same’ if new_frame is the same system as this object but no transformation is defined.

Notes

A return value of ‘same’ means the transformation will work, but it will just give back a copy of this object. The intended usage is:

if coord.is_transformable_to(some_unknown_frame):
    coord2 = coord.transform_to(some_unknown_frame)

This will work even if some_unknown_frame turns out to be the same frame class as coord. This is intended for cases where the frame is the same regardless of the frame attributes (e.g. ICRS), but be aware that it might also indicate that someone forgot to define the transformation between two objects of the same frame class but with different attributes.

realize_frame(data, **kwargs)[source]

Generates a new frame with new data from another frame (which may or may not have data). Roughly speaking, the converse of replicate_without_data.

Parameters:
dataBaseRepresentation

The representation to use as the data for the new frame.

**kwargs

Any additional keywords are treated as frame attributes to be set on the new frame object. In particular, representation_type can be specified.

Returns:
frameobjBaseCoordinateFrame subclass instance

A new object in this frame, with the same frame attributes as this one, but with the data as the coordinate data.

replicate(copy=False, **kwargs)[source]

Return a replica of the frame, optionally with new frame attributes.

The replica is a new frame object that has the same data as this frame object and with frame attributes overridden if they are provided as extra keyword arguments to this method. If copy is set to True then a copy of the internal arrays will be made. Otherwise the replica will use a reference to the original arrays when possible to save memory. The internal arrays are normally not changeable by the user so in most cases it should not be necessary to set copy to True.

Parameters:
copybool, optional

If True, the resulting object is a copy of the data. When False, references are used where possible. This rule also applies to the frame attributes.

**kwargs

Any additional keywords are treated as frame attributes to be set on the new frame object.

Returns:
frameobjBaseCoordinateFrame subclass instance

Replica of this object, but possibly with new frame attributes.

replicate_without_data(copy=False, **kwargs)[source]

Return a replica without data, optionally with new frame attributes.

The replica is a new frame object without data but with the same frame attributes as this object, except where overridden by extra keyword arguments to this method. The copy keyword determines if the frame attributes are truly copied vs being references (which saves memory for cases where frame attributes are large).

This method is essentially the converse of realize_frame.

Parameters:
copybool, optional

If True, the resulting object has copies of the frame attributes. When False, references are used where possible.

**kwargs

Any additional keywords are treated as frame attributes to be set on the new frame object.

Returns:
frameobjBaseCoordinateFrame subclass instance

Replica of this object, but without data and possibly with new frame attributes.

represent_as(base, s='base', in_frame_units=False)[source]

Generate and return a new representation of this frame’s data as a Representation object.

Note: In order to make an in-place change of the representation of a Frame or SkyCoord object, set the representation attribute of that object to the desired new representation, or use the set_representation_cls method to also set the differential.

Parameters:
basesubclass of BaseRepresentation or python:str

The type of representation to generate. Must be a class (not an instance), or the string name of the representation class.

ssubclass of BaseDifferential, python:str, optional

Class in which any velocities should be represented. Must be a class (not an instance), or the string name of the differential class. If equal to ‘base’ (default), inferred from the base class. If None, all velocity information is dropped.

in_frame_unitsbool, keyword-only

Force the representation units to match the specified units particular to this frame

Returns:
newrepBaseRepresentation-derived object

A new representation object of this frame’s data.

Raises:
AttributeError

If this object had no data

Examples

>>> from astropy import units as u
>>> from astropy.coordinates import SkyCoord, CartesianRepresentation
>>> coord = SkyCoord(0*u.deg, 0*u.deg)
>>> coord.represent_as(CartesianRepresentation)  
<CartesianRepresentation (x, y, z) [dimensionless]
        (1., 0., 0.)>
>>> coord.representation_type = CartesianRepresentation
>>> coord  
<SkyCoord (ICRS): (x, y, z) [dimensionless]
    (1., 0., 0.)>
separation(other)[source]

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.

Parameters:
otherBaseCoordinateFrame

The coordinate to get the separation to.

Returns:
sepAngle

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].

separation_3d(other)[source]

Computes three dimensional separation between this coordinate and another.

Parameters:
otherBaseCoordinateFrame

The coordinate system to get the distance to.

Returns:
sepDistance

The real-space distance between these two coordinates.

Raises:
ValueError

If this or the other coordinate do not have distances.

set_representation_cls(base=None, s='base')[source]

Set representation and/or differential class for this frame’s data.

Parameters:
basepython:str, BaseRepresentation subclass, optional

The name or subclass to use to represent the coordinate data.

sBaseDifferential subclass, optional

The differential subclass to use to represent any velocities, such as proper motion and radial velocity. If equal to ‘base’, which is the default, it will be inferred from the representation. If None, the representation will drop any differentials.

transform_to(new_frame)[source]

Transform this object’s coordinate data to a new frame.

Parameters:
new_frameastropy:coordinate-like or BaseCoordinateFrame subclass instance

The frame to transform this coordinate frame into. The frame class option is deprecated.

Returns:
transframeastropy:coordinate-like

A new object with the coordinate data represented in the newframe system.

Raises:
ValueError

If there is no possible transformation route.