CartesianRepresentation#

class astropy.coordinates.CartesianRepresentation(x, y=None, z=None, unit=None, xyz_axis=None, differentials=None, copy=True)[source]#

Bases: BaseRepresentation

Representation of points in 3D cartesian coordinates.

Parameters:
x, y, zQuantity or array

The x, y, and z coordinates of the point(s). If x, y, and z have different shapes, they should be broadcastable. If not quantity, unit should be set. If only x is given, it is assumed that it contains an array with the 3 coordinates stored along xyz_axis.

unitastropy:unit-like

If given, the coordinates will be converted to this unit (or taken to be in this unit if not given.

xyz_axisint, optional

The axis along which the coordinates are stored when a single array is provided rather than distinct x, y, and z (default: 0).

differentialsdict, CartesianDifferential, optional

Any differential classes that should be associated with this representation. The input must either be a single CartesianDifferential instance, or a dictionary of CartesianDifferential s with keys set to a string representation of the SI unit with which the differential (derivative) is taken. For example, for a velocity differential on a positional representation, the key would be 's' for seconds, indicating that the derivative is a time derivative.

copybool, optional

If True (default), arrays will be copied. If False, arrays will be references, though possibly broadcast to ensure matching shapes.

Attributes Summary

attr_classes

x

The 'x' component of the points(s).

xyz

Return a vector array of the x, y, and z coordinates.

y

The 'y' component of the points(s).

z

The 'z' component of the points(s).

Methods Summary

cross(other)

Cross product of two representations.

dot(other)

Dot product of two representations.

from_cartesian(other)

Create a representation of this class from a supplied Cartesian one.

get_xyz([xyz_axis])

Return a vector array of the x, y, and z coordinates.

mean(*args, **kwargs)

Vector mean.

norm()

Vector norm.

scale_factors()

Scale factors for each component's direction.

sum(*args, **kwargs)

Vector sum.

to_cartesian()

Convert the representation to its Cartesian form.

transform(matrix)

Transform the cartesian coordinates using a 3x3 matrix.

unit_vectors()

Cartesian unit vectors in the direction of each component.

Attributes Documentation

attr_classes = {'x': <class 'astropy.units.quantity.Quantity'>, 'y': <class 'astropy.units.quantity.Quantity'>, 'z': <class 'astropy.units.quantity.Quantity'>}#
x#

The ‘x’ component of the points(s).

xyz#

Return a vector array of the x, y, and z coordinates.

Parameters:
xyz_axisint, optional

The axis in the final array along which the x, y, z components should be stored (default: 0).

Returns:
xyzQuantity

With dimension 3 along xyz_axis. Note that, if possible, this will be a view.

y#

The ‘y’ component of the points(s).

z#

The ‘z’ component of the points(s).

Methods Documentation

cross(other)[source]#

Cross product of two representations.

Parameters:
otherBaseRepresentation subclass instance

If not already cartesian, it is converted.

Returns:
cross_productCartesianRepresentation

With vectors perpendicular to both self and other.

dot(other)[source]#

Dot product of two representations.

Note that any associated differentials will be dropped during this operation.

Parameters:
otherBaseRepresentation subclass instance

If not already cartesian, it is converted.

Returns:
dot_productQuantity

The sum of the product of the x, y, and z components of self and other.

classmethod from_cartesian(other)[source]#

Create a representation of this class from a supplied Cartesian one.

Parameters:
otherCartesianRepresentation

The representation to turn into this class

Returns:
representationBaseRepresentation subclass instance

A new representation of this class’s type.

get_xyz(xyz_axis=0)[source]#

Return a vector array of the x, y, and z coordinates.

Parameters:
xyz_axisint, optional

The axis in the final array along which the x, y, z components should be stored (default: 0).

Returns:
xyzQuantity

With dimension 3 along xyz_axis. Note that, if possible, this will be a view.

mean(*args, **kwargs)[source]#

Vector mean.

Returns a new CartesianRepresentation instance with the means of the x, y, and z components.

Refer to mean for full documentation of the arguments, noting that axis is the entry in the shape of the representation, and that the out argument cannot be used.

norm()[source]#

Vector norm.

The norm is the standard Frobenius norm, i.e., the square root of the sum of the squares of all components with non-angular units.

Note that any associated differentials will be dropped during this operation.

Returns:
normastropy.units.Quantity

Vector norm, with the same shape as the representation.

scale_factors()[source]#

Scale factors for each component’s direction.

Given unit vectors \(\hat{e}_c\) and scale factors \(f_c\), a change in one component of \(\delta c\) corresponds to a change in representation of \(\delta c \times f_c \times \hat{e}_c\).

Returns:
scale_factorsdict of Quantity

The keys are the component names.

sum(*args, **kwargs)[source]#

Vector sum.

Returns a new CartesianRepresentation instance with the sums of the x, y, and z components.

Refer to sum for full documentation of the arguments, noting that axis is the entry in the shape of the representation, and that the out argument cannot be used.

to_cartesian()[source]#

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:
cartreprCartesianRepresentation

The representation in Cartesian form.

transform(matrix)[source]#

Transform the cartesian coordinates using a 3x3 matrix.

This returns a new representation and does not modify the original one. Any differentials attached to this representation will also be transformed.

Parameters:
matrixndarray

A 3x3 transformation matrix, such as a rotation matrix.

Examples

We can start off by creating a cartesian representation object:

>>> from astropy import units as u
>>> from astropy.coordinates import CartesianRepresentation
>>> rep = CartesianRepresentation([1, 2] * u.pc,
...                               [2, 3] * u.pc,
...                               [3, 4] * u.pc)

We now create a rotation matrix around the z axis:

>>> from astropy.coordinates.matrix_utilities import rotation_matrix
>>> rotation = rotation_matrix(30 * u.deg, axis='z')

Finally, we can apply this transformation:

>>> rep_new = rep.transform(rotation)
>>> rep_new.xyz  
<Quantity [[ 1.8660254 , 3.23205081],
           [ 1.23205081, 1.59807621],
           [ 3.        , 4.        ]] pc>
unit_vectors()[source]#

Cartesian unit vectors in the direction of each component.

Given unit vectors \(\hat{e}_c\) and scale factors \(f_c\), a change in one component of \(\delta c\) corresponds to a change in representation of \(\delta c \times f_c \times \hat{e}_c\).

Returns:
unit_vectorsdict of CartesianRepresentation

The keys are the component names.