rowan

Overview

rowan.conjugate Conjugates an array of quaternions.
rowan.inverse Computes the inverse of an array of quaternions.
rowan.exp Computes the natural exponential function \(e^q\).
rowan.expb Computes the exponential function \(b^q\).
rowan.exp10 Computes the exponential function \(10^q\).
rowan.log Computes the quaternion natural logarithm.
rowan.logb Computes the quaternion logarithm to some base b.
rowan.log10 Computes the quaternion logarithm base 10.
rowan.multiply Multiplies two arrays of quaternions.
rowan.divide Divides two arrays of quaternions.
rowan.norm Compute the quaternion norm.
rowan.normalize Normalize quaternions.
rowan.rotate Rotate a list of vectors by a corresponding set of quaternions.
rowan.vector_vector_rotation Find the quaternion to rotate one vector onto another.
rowan.from_euler Convert Euler angles to quaternions.
rowan.to_euler Convert quaternions to Euler angles.
rowan.from_matrix Convert the rotation matrices mat to quaternions.
rowan.to_matrix Convert quaternions into rotation matrices.
rowan.from_axis_angle Find quaternions to rotate a specified angle about a specified axis.
rowan.to_axis_angle Convert the quaternions in q to axis angle representations.
rowan.from_mirror_plane Generate quaternions from mirror plane equations.
rowan.reflect Reflect a list of vectors by a corresponding set of quaternions.
rowan.equal Check whether two sets of quaternions are equal.
rowan.not_equal Check whether two sets of quaternions are not equal.
rowan.isfinite Test element-wise for finite quaternions.
rowan.isinf Test element-wise for infinite quaternions.
rowan.isnan Test element-wise for NaN quaternions.

Details

The core rowan package contains functions for operating on quaternions. The core package is focused on robust implementations of key functions like multiplication, exponentiation, norms, and others. Simple functionality such as addition is inherited directly from NumPy due to the representation of quaternions as NumPy arrays. Many core NumPy functions implemented for normal arrays are reimplemented to work on quaternions ( such as allclose() and isfinite()). Additionally, NumPy broadcasting is enabled throughout rowan unless otherwise specified. This means that any function of 2 (or more) quaternions can take arrays of shapes that do not match and return results according to NumPy’s broadcasting rules.

rowan.allclose(p, q, **kwargs)

Check whether two sets of quaternions are all close.

This is a direct wrapper of the corresponding NumPy function.

Parameters:
  • p ((..,4) np.array) – First array of quaternions.
  • q ((..,4) np.array) – Second array of quaternions.
  • **kwargs – Keyword arguments to pass to np.allclose.
Returns:

Boolean indicating whether or not all quaternions are close.

Example:

rowan.allclose([1, 0, 0, 0], [1, 0, 0, 0])
rowan.conjugate(q)

Conjugates an array of quaternions.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing conjugates of q.

Example:

q_star = rowan.conjugate([1, 0, 0, 0])
rowan.divide(qi, qj)

Divides two arrays of quaternions.

Division is non-commutative; this function returns \(q_i q_j^{-1}\).

Parameters:
  • qi ((..,4) np.array) – Dividend quaternions.
  • qj ((..,4) np.array) – Divisor quaternions.
Returns:

Array of shape (…) containing element-wise quotients of qi and qj.

Example:

quot = rowan.divide([1, 0, 0, 0], [2, 0, 0, 0])
rowan.exp(q)

Computes the natural exponential function \(e^q\).

The exponential of a quaternion in terms of its scalar and vector parts \(q = a + \boldsymbol{v}\) is defined by exponential power series: formula \(e^x = \sum_{k=0}^{\infty} \frac{x^k}{k!}\) as follows:

\[\begin{split}\begin{align} e^q &= e^{a+v} \\ &= e^a \left(\sum_{k=0}^{\infty} \frac{v^k}{k!} \right) \\ &= e^a \left(\cos \lvert \lvert \boldsymbol{v} \rvert \rvert + \frac{\boldsymbol{v}}{\lvert \lvert \boldsymbol{v} \rvert \rvert} \sin \lvert \lvert \boldsymbol{v} \rvert \rvert \right) \end{align}\end{split}\]
Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing exponentials of q.

Example:

q_exp = rowan.exp([1, 0, 0, 0])
rowan.expb(q, b)

Computes the exponential function \(b^q\).

We define the exponential of a quaternion to an arbitrary base relative to the exponential function \(e^q\) using the change of base formula as follows:

\[\begin{split}\begin{align} b^q &= y \\ q &= \log_b y = \frac{\ln y}{\ln b}\\ y &= e^{q\ln b} \end{align}\end{split}\]
Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing exponentials of q.

Example:

q_exp = rowan.expb([1, 0, 0, 0], 2)
rowan.exp10(q)

Computes the exponential function \(10^q\).

Wrapper around expb().

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing exponentials of q.

Example:

q_exp = rowan.exp10([1, 0, 0, 0])
rowan.equal(p, q)

Check whether two sets of quaternions are equal.

This function is a simple wrapper that checks array equality and then aggregates along the quaternion axis.

Parameters:
  • p ((..,4) np.array) – First array of quaternions.
  • q ((..,4) np.array) – Second array of quaternions.
Returns:

A boolean array of shape (…) indicating equality.

Example:

rowan.equal([1, 0, 0, 0], [1, 0, 0, 0])
rowan.from_axis_angle(axes, angles)

Find quaternions to rotate a specified angle about a specified axis.

Parameters:
  • axes ((..,3) np.array) – An array of vectors (the axes).
  • angles (float or (..,1) np.array) – An array of angles in radians. Will be broadcast to match shape of v as needed.
Returns:

Array of shape (…, 4) containing the corresponding rotation quaternions.

Example:

quat = rowan.from_axis_angle([[1, 0, 0]], np.pi/3)
rowan.from_euler(alpha, beta, gamma, convention='zyx', axis_type='intrinsic')

Convert Euler angles to quaternions.

For generality, the rotations are computed by composing a sequence of quaternions corresponding to axis-angle rotations. While more efficient implementations are possible, this method was chosen to prioritize flexibility since it works for essentially arbitrary Euler angles as long as intrinsic and extrinsic rotations are not intermixed.

Parameters:
  • alpha ((..) np.array) – Array of \(\alpha\) values in radians.
  • beta ((..) np.array) – Array of \(\beta\) values in radians.
  • gamma ((..) np.array) – Array of \(\gamma\) values in radians.
  • convention (str) – One of the 12 valid conventions xzx, xyx, yxy, yzy, zyz, zxz, xzy, xyz, yxz, yzx, zyx, zxy.
  • axes (str) – Whether to use extrinsic or intrinsic rotations.
Returns:

Array of shape (…, 4) containing quaternions corresponding to the input angles.

Example:

ql = rowan.from_euler(0.3, 0.5, 0.7)
rowan.from_matrix(mat, require_orthogonal=True)

Convert the rotation matrices mat to quaternions.

This method uses the algorithm described by Bar-Itzhack in [Itzhack00]. The idea is to construct a matrix K whose largest eigenvalue corresponds to the desired quaternion. One of the strengths of the algorithm is that for nonorthogonal matrices it gives the closest quaternion representation rather than failing outright.

[Itzhack00]Itzhack Y. Bar-Itzhack. “New Method for Extracting the Quaternion from a Rotation Matrix”, Journal of Guidance, Control, and Dynamics, Vol. 23, No. 6 (2000), pp. 1085-1087 https://doi.org/10.2514/2.4654
Parameters:mat ((..,3,3) np.array) – An array of rotation matrices.
Returns:Array of shape (…, 4) containing the corresponding rotation quaternions.

Example:

ql = rowan.from_matrix([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
rowan.from_mirror_plane(x, y, z)

Generate quaternions from mirror plane equations.

Reflection quaternions can be constructed from the form \((0, x, y, z)\), i.e. with zero real component. The vector \((x, y, z)\) is the normal to the mirror plane.

Parameters:
  • x ((..) np.array) – First planar component.
  • y ((..) np.array) – Second planar component.
  • z ((..) np.array) – Third planar component.
Returns:

Array of shape (…) containing quaternions reflecting about the input plane \((x, y, z)\).

Example:

quat_ref = rowan.from_mirror_plane(*(1, 2, 3))
rowan.inverse(q)

Computes the inverse of an array of quaternions.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing inverses of q.

Example:

q_inv = rowan.inverse([1, 0, 0, 0])
rowan.isclose(p, q, **kwargs)

Element-wise check of whether two sets of quaternions are close.

This function is a simple wrapper that checks using the corresponding NumPy function and then aggregates along the quaternion axis.

Parameters:
  • p ((..,4) np.array) – First array of quaternions.
  • q ((..,4) np.array) – Second array of quaternions.
  • **kwargs – Keyword arguments to pass to np.isclose.
Returns:

A boolean array of shape (…) indicating which quaternions are close.

Example:

rowan.allclose([[1, 0, 0, 0]], [[1, 0, 0, 0]])
rowan.isinf(q)

Test element-wise for infinite quaternions.

A quaternion is defined as infinite if any elements are infinite.

Parameters:q ((..,4) np.array) – Array of quaternions
Returns:A boolean array of shape (…) indicating infinite quaternions.

Example:

import numpy as np
rowan.isinf([np.nan, 0, 0, 0])
rowan.isfinite(q)

Test element-wise for finite quaternions.

A quaternion is defined as finite if all elements are finite.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:A boolean array of shape (…) indicating finite quaternions.

Example:

rowan.isfinite([1, 0, 0, 0])
rowan.isnan(q)

Test element-wise for NaN quaternions.

A quaternion is defined as NaN if any elements are NaN.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:A boolean array of shape (…) indicating whether or not the input quaternions were NaN.

Example:

import numpy as np
rowan.isnan([np.nan, 0, 0, 0])
rowan.is_unit(q)

Check if all input quaternions have unit norm.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Whether or not all inputs are unit quaternions
Return type:bool

Example:

rowan.is_unit([10, 0, 0, 0])
rowan.log(q)

Computes the quaternion natural logarithm.

The natural of a quaternion in terms of its scalar and vector parts \(q = a + \boldsymbol{v}\) is defined by inverting the exponential formula (see exp()), and is defined by the formula \(\frac{x^k}{k!}\) as follows:

\[\begin{equation} \ln(q) = \ln\lvert\lvert q \rvert\rvert + \frac{\boldsymbol{v}}{\lvert\lvert \boldsymbol{v} \rvert\rvert} \arccos\left(\frac{a}{q}\right) \end{equation}\]
Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing logarithms of q.

Example:

ln_q  = rowan.log([1, 0, 0, 0])
rowan.logb(q, b)

Computes the quaternion logarithm to some base b.

The quaternion logarithm for arbitrary bases is defined using the standard change of basis formula relative to the natural logarithm.

\[\begin{split}\begin{align} \log_b q &= y \\ q &= b^y \\ \ln q &= y \ln b \\ y &= \log_b q = \frac{\ln q}{\ln b} \end{align}\end{split}\]
Parameters:
  • q ((..,4) np.array) – Array of quaternions.
  • n ((..) np.array) – Scalars to use as log bases.
Returns:

Array of shape (…) containing logarithms of q.

Example:

log2_q = rowan.logb([1, 0, 0, 0], 2)
rowan.log10(q)

Computes the quaternion logarithm base 10.

Wrapper around logb().

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing logarithms of q.

Example:

log10_q = rowan.log10([1, 0, 0, 0])
rowan.multiply(qi, qj)

Multiplies two arrays of quaternions.

Note that quaternion multiplication is generally non-commutative, so the first and second set of quaternions must be passed in the correct order.

Parameters:
  • qi ((..,4) np.array) – Array of left quaternions.
  • qj ((..,4) np.array) – Array of right quaternions.
Returns:

Array of shape (…) containing element-wise products of q.

Example:

prod = rowan.multiply([1, 0, 0, 0], [2, 0, 0, 0])
rowan.norm(q)

Compute the quaternion norm.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) containing norms of q.

Example:

norms = rowan.norm([10, 0, 0, 0])
rowan.normalize(q)

Normalize quaternions.

Parameters:q ((..,4) np.array) – Array of quaternions.
Returns:Array of shape (…) of normalized quaternions.

Example:

u = rowan.normalize([10, 0, 0, 0])
rowan.not_equal(p, q)

Check whether two sets of quaternions are not equal.

This function is a simple wrapper that checks array equality and then aggregates along the quaternion axis.

Parameters:
  • p ((..,4) np.array) – First array of quaternions.
  • q ((..,4) np.array) – Second array of quaternions.
Returns:

A boolean array of shape (…) indicating inequality.

Example:

rowan.not_equal([-1, 0, 0, 0], [1, 0, 0, 0])
rowan.power(q, n)

Computes the power of a quaternion \(q^n\).

Quaternions raised to a scalar power are defined according to the polar decomposition angle \(\theta\) and vector \(\hat{u}\): \(q^n = \lvert\lvert q \rvert\rvert^n \left( \cos(n\theta) + \hat{u} \sin(n\theta)\right)\). However, this can be computed more efficiently by noting that \(q^n = \exp(n \ln(q))\).

Parameters:
  • q ((..,4) np.array) – Array of quaternions.
  • n ((..) np.arrray) – Scalars to exponentiate quaternions with.
Returns:

Array of shape (…) containing powers of q.

Example:

q_5 = rowan.power([1, 0, 0, 0], 5)
rowan.reflect(q, v)

Reflect a list of vectors by a corresponding set of quaternions.

For help constructing a mirror plane, see from_mirror_plane().

Parameters:
  • q ((..,4) np.array) – Array of quaternions.
  • v ((..,3) np.array) – Array of vectors.
Returns:

Array of shape (…, 3) containing reflections of v.

Example:

v_reflected = rowan.reflect([1, 0, 0, 0], [1, 1, 1])
rowan.rotate(q, v)

Rotate a list of vectors by a corresponding set of quaternions.

Parameters:
  • q ((..,4) np.array) – Array of quaternions.
  • v ((..,3) np.array) – Array of vectors.
Returns:

Array of shape (…, 3) containing rotations of v.

Example:

v_rot = rowan.reflect([1, 0, 0, 0], [1, 1, 1])
rowan.to_axis_angle(q)

Convert the quaternions in q to axis angle representations.

Parameters:q ((..,4) np.array) – An array of quaternions.
Returns:A tuple of np.arrays (axes, angles) where axes has shape (…,3) and angles has shape (…,1). The angles are in radians.

Example:

quat = rowan.to_axis_angle([[1, 0, 0, 0]])
rowan.to_euler(q, convention='zyx', axis_type='intrinsic')

Convert quaternions to Euler angles.

Euler angles are returned in the sequence provided, so in, e.g., the default case (‘zyx’), the angles returned are for a rotation \(Z(\alpha) Y(\beta) X(\gamma)\).

Note

In all cases, the \(\alpha\) and \(\gamma\) angles are between \(\pm \pi\). For proper Euler angles, \(\beta\) is between \(0\) and \(pi\) degrees. For Tait-Bryan angles, \(\beta\) lies between \(\pm\pi/2\).

For simplicity, quaternions are converted to matrices, which are then converted to their Euler angle representations. All equations for rotations are derived by considering compositions of the three elemental rotations about the three Cartesian axes:

\begin{eqnarray*} R_x(\theta) =& \left(\begin{array}{ccc} 1 & 0 & 0 \\ 0 & \cos \theta & -\sin \theta \\ 0 & \sin \theta & \cos \theta \\ \end{array}\right)\\ R_y(\theta) =& \left(\begin{array}{ccc} \cos \theta & 0 & \sin \theta \\ 0 & 1 & 0 \\ -\sin \theta & 1 & \cos \theta \\ \end{array}\right)\\ R_z(\theta) =& \left(\begin{array}{ccc} \cos \theta & -\sin \theta & 0 \\ \sin \theta & \cos \theta & 0 \\ 0 & 0 & 1 \\ \end{array}\right)\\ \end{eqnarray*}

Extrinsic rotations are represented by matrix multiplications in the proper order, so \(z-y-x\) is represented by the multiplication \(XYZ\) so that the system is rotated first about \(Z\), then about \(Y\), then finally \(X\). For intrinsic rotations, the order of rotations is reversed, meaning that it matches the order in which the matrices actually appear i.e. the \(z-y'-x''\) convention (yaw, pitch, roll) corresponds to the multiplication of matrices \(ZYX\). For proof of the relationship between intrinsic and extrinsic rotations, see the Wikipedia page on Davenport chained rotations.

For more information, see the Wikipedia page for Euler angles (specifically the section on converting between representations).

Parameters:
  • q ((..,4) np.array) – Quaternions to transform.
  • convention (str) – One of the 6 valid conventions zxz, xyx, yzy, zyz, xzx, yxy.
  • axes (str) – Whether to use extrinsic or intrinsic.
Returns:

math:(alpha, beta, gamma) as the last dimension (in radians).

Return type:

Array of shape (.., 3) containing Euler angles

Example:

import numpy as np
rands = np.random.rand(100, 3)
alpha, beta, gamma = rands.T
ql = rowan.from_euler(alpha, beta, gamma)
alpha_return, beta_return, gamma_return = np.split(
    rowan.to_euler(ql), 3, axis = 1)
assert(np.allclose(alpha_return.flatten(), alpha))
assert(np.allclose(beta_return.flatten(), beta))
assert(np.allclose(gamma_return.flatten(), gamma))
rowan.to_matrix(q, require_unit=True)

Convert quaternions into rotation matrices.

Uses the conversion described on Wikipedia.

Parameters:q ((..,4) np.array) – An array of quaternions.
Returns:Array of shape (…, 3, 3) containing the corresponding rotation matrices.

Example:

ql = rowan.to_matrix([1, 0, 0, 0]])
rowan.vector_vector_rotation(v1, v2)

Find the quaternion to rotate one vector onto another.

Parameters:
  • v1 ((..,3) np.array) – Array of vectors to rotate.
  • v2 ((..,3) np.array) – Array of vector to rotate onto.
Returns:

Array of shape (…, 4) containing quaternions that rotate v1 onto v2.

Example:

q_rot = rowan.vector_vector_rotation([1, 0, 0], [0, 1, 0])