spacepy.coordinates

Implementation of Coords class functions for coordinate transformations

Authors: Steven Morley and Josef Koller Institution: Los ALamos National Laboratory Contact: smorley@lanl.gov

For additional documentation Coordinates - Implementation of Coords class functions for coordinate transformations

Copyright 2010-2016 Los Alamos National Security, LLC.

Functions

car2sph(car_in)

Coordinate transformation from Cartesian to spherical

quaternionConjugate(Qin[, scalarPos])

Given an input quaternion (or array of quaternions), return the conjugate

quaternionFromMatrix(matrix[, scalarPos])

Given an input rotation matrix, return the equivalent quaternion

quaternionMultiply(Qin1, Qin2[, scalarPos])

Given quaternions, return the product, i.e. Qin1*Qin2.

quaternionNormalize(Qin[, scalarPos])

Given an input quaternion (or array of quaternions), return the unit quaternion

quaternionRotateVector(Qin, Vin[, ...])

Given quaternions and vectors, return the vectors rotated by the quaternions

quaternionToMatrix(Qin[, scalarPos, normalize])

Given an input quaternion, return the equivalent rotation matrix.

sph2car(sph_in)

Coordinate transformation from spherical to Cartesian

Classes

Coords(data, dtype, carsph, [units, ticks, ...)

A class holding spatial coordinates and enabling transformation between coordinate systems.

Attributes

spacepy.coordinates.car2sph(car_in)[source]

Coordinate transformation from Cartesian to spherical

Parameters:
- car_in (list or ndarray)coordinate points in (n,3) shape with n coordinate points in

units of [Re, Re, Re] = [x,y,z]

Returns:
- results (ndarray)values after conversion to spherical coordinates in

radius, latitude, longitude in units of [Re, deg, deg]

See also

sph2car

Examples

>>> sph2car([1,45,0])
array([ 0.70710678,  0.        ,  0.70710678])
spacepy.coordinates.quaternionConjugate(Qin, scalarPos='last')[source]

Given an input quaternion (or array of quaternions), return the conjugate

Parameters:
Qinarray_like

input quaternion to conjugate

Returns:
outarray_like

conjugate quaternion

Examples

>>> import spacepy.coordinates
>>> spacepy.coordinates.quaternionConjugate(
...     [0.707, 0, 0.707, 0.2], scalarPos='last')
array([-0.707, -0.   , -0.707,  0.2  ])
spacepy.coordinates.quaternionFromMatrix(matrix, scalarPos='last')[source]

Given an input rotation matrix, return the equivalent quaternion

The output has one fewer axis than the input (the last axis) and the shape is otherwise unchanged, allowing multi-dimensional matrix input.

Parameters:
matrixarray_like

input rotation matrix or array of matrices

Returns:
outarray_like

Quaternions representing the same rotation as the input rotation matrices.

Other Parameters:
scalarPosstr

Location of the scalar component of the output quaternion, either ‘last’ (default) or ‘first’.

Raises:
NotImplementedError

for invalid values of scalarPos

ValueError

for inputs which are obviously not valid 3D rotation matrices or arrays thereof: if the size doesn’t end in (3, 3), if the matrix is not orthogonal, or not a proper rotation.

Notes

Added in version 0.2.2.

No attempt is made to resolve the sign ambiguity; in particular, conversions of very similar matrices may result in equivalent quaternions with the opposite sign. This may have implications for interpolating a sequence of quaternions.

The conversion of a rotation matrix to a quaternion suffers from some of the same disadvantages inherent to rotation matrices, such as potential numerical instabilities. Working in quaternion space as much as possible is recommended.

There are several algorithms; the most well-known algorithm for this conversion is Shepperd’s [1], although the many “rediscoveries” indicate it is not sufficiently well-known. This function uses the method of Bar-Itzhack [2] (version 3), which should be resistant to small errors in the rotation matrix. As a result, the input checking is quite coarse and will likely accept many matrices that do not represent valid rotations.

Also potentially of interest, although not implemented here, is Sarabandi and Thomas [3].

References

Examples

>>> import spacepy.coordinates
>>> spacepy.coordinates.quaternionFromMatrix(
...     [[ 0.,  0.,  1.],
...      [ 1.,  0.,  0.],
...      [ 0.,  1.,  0.]])
array([0.5, 0.5, 0.5, 0.5])
spacepy.coordinates.quaternionMultiply(Qin1, Qin2, scalarPos='last')[source]

Given quaternions, return the product, i.e. Qin1*Qin2

Parameters:
Qin1array_like

input quaternion, first position

Qin2array-like

input quaternion, second position

Returns:
outarray_like

quaternion product

Examples

>>> import spacepy.coordinates
>>> import numpy as np
>>> vecX = [1, 0, 0] #shared X-axis
>>> vecZ = [0, 0, 1] #unshared, but similar, Z-axis
>>> quat_eci_to_gsm = [-0.05395384,  0.07589845, -0.15172533,  0.98402634]
>>> quat_eci_to_gse = [ 0.20016056,  0.03445775, -0.16611386,  0.96496352]
>>> quat_gsm_to_eci = spacepy.coordinates.quaternionConjugate(
...     quat_eci_to_gsm)
>>> quat_gse_to_gsm = spacepy.coordinates.quaternionMultiply(
...     quat_gsm_to_eci, quat_eci_to_gse)
>>> spacepy.coordinates.quaternionRotateVector(quat_gse_to_gsm, vecX)
array([  1.00000000e+00,   1.06536725e-09,  -1.16892107e-08])
>>> spacepy.coordinates.quaternionRotateVector(quat_gse_to_gsm, vecZ)
array([  1.06802834e-08,  -4.95669027e-01,   8.68511494e-01])
spacepy.coordinates.quaternionNormalize(Qin, scalarPos='last')[source]

Given an input quaternion (or array of quaternions), return the unit quaternion

Parameters:
vecarray_like

input quaternion to normalize

Returns:
outarray_like

normalized quaternion

Examples

>>> import spacepy.coordinates
>>> spacepy.coordinates.quaternionNormalize([0.707, 0, 0.707, 0.2])
array([ 0.69337122,  0.        ,  0.69337122,  0.19614462])
spacepy.coordinates.quaternionRotateVector(Qin, Vin, scalarPos='last', normalize=True)[source]

Given quaternions and vectors, return the vectors rotated by the quaternions

Parameters:
Qinarray_like

input quaternion to rotate by

Vinarray-like

input vector to rotate

Returns:
outarray_like

rotated vector

Examples

>>> import spacepy.coordinates
>>> import numpy as np
>>> vec = [1, 0, 0]
>>> quat_wijk = [np.sin(np.pi/4), 0, np.sin(np.pi/4), 0.0]
>>> quat_ijkw = [0.0, np.sin(np.pi/4), 0, np.sin(np.pi/4)]
>>> spacepy.coordinates.quaternionRotateVector(quat_ijkw, vec)
array([ 0.,  0., -1.])
>>> spacepy.coordinates.quaternionRotateVector(
...     quat_wijk, vec, scalarPos='first')
array([ 0.,  0., -1.])
spacepy.coordinates.quaternionToMatrix(Qin, scalarPos='last', normalize=True)[source]

Given an input quaternion, return the equivalent rotation matrix.

The output has one more axis than the input (the last axis) and the shape is otherwise unchanged, allowing multi-dimensional quaternion input.

Parameters:
Qinarray_like

input quaternion or array of quaternions, must be normalized.

Returns:
outarray_like

Rotation matrix

Other Parameters:
scalarPosstr

Location of the scalar component of the input quaternion, either ‘last’ (default) or ‘first’.

normalizeTrue

Normalize input quaternions before conversion (default). If False, raises error for non-normalized.

Raises:
NotImplementedError

for invalid values of scalarPos.

ValueError

for inputs which are not valid normalized quaternions or arrays thereof: if the size doesn’t end in (4), if the quaternion is not normalized and normalize is False.

Notes

Added in version 0.2.2.

Implementation of the Euler–Rodrigues formula.

Examples

>>> import spacepy.coordinates
>>> spacepy.coordinates.quaternionToMatrix([0.5, 0.5, 0.5, 0.5])
array([[ 0.,  0.,  1.],
       [ 1.,  0.,  0.],
       [ 0.,  1.,  0.]])
spacepy.coordinates.sph2car(sph_in)[source]

Coordinate transformation from spherical to Cartesian

Parameters:
- sph_in (list or ndarray)coordinate points in (n,3) shape with n coordinate points in

units of [Re, deg, deg] = [r, latitude, longitude]

Returns:
- results (ndarray)values after conversion to cartesian coordinates x,y,z

See also

car2sph

Examples

>>> sph2car([1,45,45])
array([ 0.5       ,  0.5       ,  0.70710678])
spacepy.coordinates.DEFAULTS = DEFAULTS(use_irbem=False, ellipsoid={'A': 6378.137, 'iFlat': 298.257223563, 'A2': 40680631.59076899, 'Flat': 0.0033528106647474805, 'E2': 0.0066943799901413165, 'E4': 4.481472345240445e-05, '1mE2': 0.9933056200098587, 'B': 6356.752314245179, 'B2': 40408299.98466144, 'A2mB2': 272331.60610755533, 'EP2': 0.0067394967422764514}, itol=30)

DEFAULTS attribute is a named tuple with default values for coordinates module