InterpolatedRepresentation¶
- class interpolated_coordinates.InterpolatedRepresentation(representation: BaseRepresentation, *args: Any, **kwargs: Any)[source]¶
Bases:
InterpolatedBaseRepresentationOrDifferential
Wrapper for Representations, adding affine interpolations.
- Parameters:
- representation
BaseRepresentation
instance - affine
Quantity
array_like The affine interpolation parameter.
- interps
Mapping
orNone
(optional, keyword-only) Has same structure as a Representation
dict(component name: interpolation, ... "differentials": dict( "s" : dict(component name: interpolation, ...), ...))
- **interp_kwargs
Only used if
interps
is None. keyword arguments into interpolation class
- representation
- Other Parameters:
- interp_cls
Callable
(optional, keyword-only) option for ‘interp_kwargs’. If not specified, default is
IntpUnivarSplUnits
.
- interp_cls
Attributes Summary
The class used when taking a derivative.
Return the underlying Representation.
Methods Summary
__call__
([affine])Evaluate interpolated representation.
Return self, clearing cached derivatives.
copy
(*args, **kwargs)Return an instance containing copies of the internal data.
derivative
([n])Construct a new spline representing the derivative of this spline.
from_cartesian
(other)Create a representation of this class from a Cartesian one.
Headless tangent vector at each point in affine.
represent_as
(other_class[, differential_class])Convert coordinates to another representation.
Tangent vectors along the curve, from the origin.
Convert the representation to its Cartesian form.
with_differentials
(differentials)Realize Representation, with new differentials.
Return a copy of the representation without attached differentials.
Attributes Documentation
- affine¶
- derivative_type¶
The class used when taking a derivative.
- uninterpolated¶
Return the underlying Representation.
Methods Documentation
- __call__(affine: Quantity | None = None) BaseRepresentation [source]¶
Evaluate interpolated representation.
- Parameters:
- affine
Quantity
array_like The affine interpolation parameter. If None, returns representation points.
- affine
- Returns:
BaseRepresenation
Representation of type
self.data
evaluated withaffine
- clear_derivatives() IRoDType ¶
Return self, clearing cached derivatives.
- copy(*args: Any, **kwargs: Any) IRoDType ¶
Return an instance containing copies of the internal data.
Parameters are as for
copy()
.- Returns:
interpolated_coordinates.InterpolatedBaseRepresentationOrDifferential
Same type as this instance.
- derivative(n: int = 1) InterpolatedDifferential [source]¶
Construct a new spline representing the derivative of this spline.
- Parameters:
- n
int
, optional Order of derivative to evaluate. Default: 1
- n
- from_cartesian(other: CartesianRepresentation | CartesianDifferential) IRoDType ¶
Create a representation of this class from a Cartesian one.
- Parameters:
- other
CartesianRepresentation
orCartesianDifferential
The representation to turn into this class
Note: the affine parameter of this class is used. The representation must be the same length as the affine parameter.
- other
- Returns:
- representation
object
ofthis
class A new representation of this class’s type.
- representation
- Raises:
ValueError
If
other
is not same length as the this instance’s affine parameter.
- headless_tangent_vectors() IRType [source]¶
Headless tangent vector at each point in affine.
\(\vec{x} + \partial_{\affine} \vec{x}(\affine) \Delta\affine\)
- represent_as(other_class: BaseRepresentation, differential_class: BaseDifferential | None = None) InterpolatedRepresentation [source]¶
Convert coordinates to another representation.
If the instance is of the requested class, it is returned unmodified. By default, conversion is done via Cartesian coordinates. Also note that orientation information at the origin is not preserved by conversions through Cartesian coordinates. See the docstring for
represent_as()
for an example.- Parameters:
- other_class
BaseRepresentation
subclass The type of representation to turn the coordinates into.
- differential_class
dict
ofBaseDifferential
, optional Classes in which the differentials should be represented. Can be a single class if only a single differential is attached, otherwise it should be a
dict
keyed by the same keys as the differentials.
- other_class
- tangent_vectors() IRType [source]¶
Tangent vectors along the curve, from the origin.
\(\vec{x} + \partial_{\affine} \vec{x}(\affine) \Delta\affine\)
- to_cartesian() IRoDType ¶
Convert the representation to its Cartesian form.
Note that any differentials get dropped. Also note that orientation information at the origin is not preserved by conversions through Cartesian coordinates. For example, transforming an angular position defined at distance=0 through cartesian coordinates and back will lose the original angular coordinates:
>>> import astropy.units as u >>> import astropy.coordinates as coord >>> rep = coord.SphericalRepresentation( ... lon=15*u.deg, ... lat=-11*u.deg, ... distance=0*u.pc) >>> rep.to_cartesian().represent_as(coord.SphericalRepresentation) <SphericalRepresentation (lon, lat, distance) in (rad, rad, pc) (0., 0., 0.)>
- Returns:
- cartrepr
CartesianRepresentation
orCartesianDifferential
The representation in Cartesian form. If starting from a Cart
- cartrepr
- with_differentials(differentials: Sequence[BaseDifferential]) IRType [source]¶
Realize Representation, with new differentials.
Create a new representation with the same positions as this representation, but with these new differentials.
Differential keys that already exist in this object’s differential dict are overwritten.
- Parameters:
- differentials
Sequence
ofBaseDifferential
The differentials for the new representation to have.
- differentials
- Returns:
newrepr
A copy of this representation, but with the
differentials
as its differentials.