CartesianRepresentation¶
-
class
astropy.coordinates.CartesianRepresentation(x, y=None, z=None, unit=None, xyz_axis=None, copy=True)[source] [edit on github]¶ Bases:
astropy.coordinates.BaseRepresentationRepresentation of points in 3D cartesian coordinates.
Parameters: x, y, z :
Quantityor arrayThe x, y, and z coordinates of the point(s). If
x,y, andzhave different shapes, they should be broadcastable. If not quantity,unitshould be set. If onlyxis given, it is assumed that it contains an array with the 3 coordinates are stored alongxyz_axis.unit :
Unitor strIf given, the coordinates will be converted to this unit (or taken to be in this unit if not given.
xyz_axis : int, optional
The axis along which the coordinates are stored when a single array is provided rather than distinct
x,y, andz(default: 0).copy : bool, optional
If True (default), arrays will be copied rather than referenced.
Attributes Summary
attr_classesxThe x component of the point(s). xyzReturn a vector array of the x, y, and z coordinates. yThe y component of the point(s). zThe z component of the point(s). Methods Summary
cross(other)Cross product of two representations. dot(other)Dot product of two representations. from_cartesian(other)get_xyz([xyz_axis])Return a vector array of the x, y, and z coordinates. mean(\*args, \*\*kwargs)Vector mean. sum(\*args, \*\*kwargs)Vector sum. to_cartesian()transform(matrix)Transform the cartesian coordinates using a 3x3 matrix. Attributes Documentation
-
attr_classes= OrderedDict([(u'x', <class 'astropy.units.quantity.Quantity'>), (u'y', <class 'astropy.units.quantity.Quantity'>), (u'z', <class 'astropy.units.quantity.Quantity'>)])¶
-
x¶ The x component of the point(s).
-
xyz¶ Return a vector array of the x, y, and z coordinates.
Parameters: xyz_axis : int, optional
The axis in the final array along which the x, y, z components should be stored (default: 0).
Returns: xyz :
QuantityWith dimension 3 along
xyz_axis.
-
y¶ The y component of the point(s).
-
z¶ The z component of the point(s).
Methods Documentation
-
cross(other)[source] [edit on github]¶ Cross product of two representations.
Parameters: other : representation
If not already cartesian, it is converted.
Returns: cross_product :
CartesianRepresentationWith vectors perpendicular to both
selfandother.
-
dot(other)[source] [edit on github]¶ Dot product of two representations.
Parameters: other : representation
If not already cartesian, it is converted.
Returns: dot_product :
QuantityThe sum of the product of the x, y, and z components of
selfandother.
-
classmethod
from_cartesian(other)[source] [edit on github]¶
-
get_xyz(xyz_axis=0)[source] [edit on github]¶ Return a vector array of the x, y, and z coordinates.
Parameters: xyz_axis : int, optional
The axis in the final array along which the x, y, z components should be stored (default: 0).
Returns: xyz :
QuantityWith dimension 3 along
xyz_axis.
-
mean(*args, **kwargs)[source] [edit on github]¶ Vector mean.
Returns a new CartesianRepresentation instance with the means of the x, y, and z components.
Refer to
meanfor full documentation of the arguments, noting thataxisis the entry in theshapeof the representation, and that theoutargument cannot be used.
-
sum(*args, **kwargs)[source] [edit on github]¶ Vector sum.
Returns a new CartesianRepresentation instance with the sums of the x, y, and z components.
Refer to
sumfor full documentation of the arguments, noting thataxisis the entry in theshapeof the representation, and that theoutargument cannot be used.
-
to_cartesian()[source] [edit on github]¶
-
transform(matrix)[source] [edit on github]¶ Transform the cartesian coordinates using a 3x3 matrix.
This returns a new representation and does not modify the original one.
Parameters: matrix :
ndarrayA 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>
-