InterpolatedCoordinateFrame¶
- class interpolated_coordinates.InterpolatedCoordinateFrame(data: CoordinateType, affine: Quantity | None = None, *, interps: dict[str, Any] | None = None, **interp_kwargs: Any)[source]¶
Bases:
object
Wrapper for Coordinate Frame, adding affine interpolations.
- Parameters:
- data
InterpolatedRepresentation
orBaseRepresentation
orBaseCoordinateFrame
For either an InterpolatedRepresentation or Representation the kwarg ‘frame’ must also be specified. If CoordinateFrame, then ‘frame’ is ignored.
- affine
Quantity
array_like (optional) if not a Quantity, one is assigned. Only used if data is not already interpolated. If data is NOT interpolated, this is required.
- data
- Other Parameters:
- frame
str
orBaseCoordinateFrame
only used if
data
is an InterpolatedRepresentation or Representation
- frame
- Raises:
Exception
if
frame
has no errorValueError
if
data
is not an interpolated type andaffine
is NoneTypeError
if
data
is not one of types specified in Parameters.
Attributes Summary
Return uninterpolated CoordinateFrame.
Methods Summary
__call__
([affine])Evaluate interpolated coordinate frame.
copy
()derivative
([n])Construct a new spline representing the derivative of this spline.
Headless tangent vector at each point in affine.
realize_frame
(data[, affine])Generates a new frame with new data from another frame.
represent_as
(base[, s, in_frame_units])Generate and return a new representation of this frame's
data
as a Representation object.separation
(point, *[, interpolate, affine])Computes on-sky separation between this coordinate and another.
separation_3d
(point, *[, interpolate, affine])Compute three dimensional separation between coordinates.
Tangent vectors along the curve, from the origin.
transform_to
(new_frame)Transform this object's coordinate data to a new frame.
Attributes Documentation
- affine¶
- representation_type¶
- uninterpolated¶
Return uninterpolated CoordinateFrame.
Methods Documentation
- __call__(affine: Quantity | None = None) BaseRepresentation [source]¶
Evaluate interpolated coordinate frame.
- Parameters:
- affine
Quantity
array_like The affine interpolation parameter. If None, returns representation points.
- affine
- Returns:
BaseRepresenation
Representation of type
self.data
evaluated withaffine
- derivative(n: int = 1) BaseRepresentationOrDifferential [source]¶
Construct a new spline representing the derivative of this spline.
- Parameters:
- n
int
, optional Order of derivative to evaluate. Default: 1
- n
- headless_tangent_vectors() InterpolatedCoordinateFrame [source]¶
Headless tangent vector at each point in affine.
\(\vec{x} + \partial_{\affine} \vec{x}(\affine) \Delta\affine\)
- realize_frame(data: BaseRepresentation, affine: Quantity | None = None, **kwargs: Any) InterpolatedCoordinateFrame [source]¶
Generates a new frame with new data from another frame.
The new frame may or may not have data. Roughly speaking, the converse of
replicate_without_data
.- Parameters:
- data
BaseRepresentation
The representation to use as the data for the new frame.
- affine
Quantity
|None
, optional The affine interpolation parameter.
- **kwargs: Any, optional
Any additional keywords are treated as frame attributes to be set on the new frame object. In particular,
representation_type
can be specified.
- data
- Returns:
- frameobj
same
as
this
frame
A new object with the same frame attributes as this one, but with the
data
as the coordinate data.
- frameobj
- represent_as(base: RepLikeType, s: str | BaseDifferential = 'base', in_frame_units: bool = False) InterpolatedRepresentation [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 theset_representation_cls
method to also set the differential.- Parameters:
- basesubclass of
BaseRepresentation
orstr
The type of representation to generate. Must be a class (not an instance), or the string name of the representation class.
- ssubclass of
BaseDifferential
,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
- basesubclass of
- Returns:
- newrepBaseRepresentation-derived
object
A new representation object of this frame’s
data
.
- newrepBaseRepresentation-derived
- 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(point: CoordinateType, *, interpolate: bool = True, affine: Quantity | None = None) Angle | IUSUnits [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 thatself.separation(other)
andother.separation(self)
may not give the same answer in this case.- Parameters:
- point
BaseCoordinateFrame
The coordinate to get the separation to.
- interpolatebool, optional keyword-only
Whether to return the separation as an interpolation, or as points.
- affine
Quantity
array_like orNone
, optional The affine interpolation parameter. If None (default), return angular width evaluated at all “tick” interpolation points.
- point
- Returns:
- sep
Angle
orInterpolatedUnivariateSplinewithUnits
The on-sky separation between this and the
other
coordinate.
- sep
Notes
The separation is calculated using the Vincenty formula, which is stable at all locations, including poles and antipodes [1].
- separation_3d(point: CoordinateType, *, interpolate: bool = True, affine: Quantity | None = None) Distance | IUSUnits [source]¶
Compute three dimensional separation between coordinates.
For more on how to use this (and related) functionality, see the examples in Separations, Offsets, Catalog Matching, and Related Functionality.
- Parameters:
- point
BaseCoordinateFrame
The coordinate to get the separation to.
- interpolatebool, optional keyword-only
Whether to return the separation as an interpolation, or as points.
- affine
Quantity
array_like orNone
, optional The affine interpolation parameter. If None (default), return angular width evaluated at all “tick” interpolation points.
- point
- Returns:
- sep
Distance
|InterpolatedUnivariateSplinewithUnits
The real-space distance between these two coordinates.
- sep
- Raises:
ValueError
If this or the other coordinate do not have distances.
- tangent_vectors() InterpolatedCoordinateFrame [source]¶
Tangent vectors along the curve, from the origin.
\(\vec{x} + \partial_{\affine} \vec{x}(\affine) \Delta\affine\)
- transform_to(new_frame: BaseCoordinateFrame | SkyCoord) InterpolatedCoordinateFrame [source]¶
Transform this object’s coordinate data to a new frame.
- Parameters:
- Returns:
transframe
A new object with the coordinate data represented in the
newframe
system.
- Raises:
ValueError
If there is no possible transformation route.