from_matrix#
- static Rotation.from_matrix(matrix, *, assume_valid=False)[source]#
Initialize from rotation matrix.
Rotations in 3 dimensions can be represented with 3 x 3 orthogonal matrices [1]. If the input is not orthogonal, an approximation is created by orthogonalizing the input matrix using the method described in [2], and then converting the orthogonal rotation matrices to quaternions using the algorithm described in [3]. Matrices must be right-handed.
- Parameters:
- matrixarray_like, shape (…, 3, 3)
A single matrix or an ND array of matrices, where the last two dimensions contain the rotation matrices.
- assume_validbool, optional
Must be False unless users can guarantee the input is a valid rotation matrix, i.e. it is orthogonal, rows and columns have unit norm and the determinant is 1. Setting this to True without ensuring these properties is unsafe and will silently lead to incorrect results. If True, normalization steps are skipped, which can improve runtime performance. Default is False.
- Returns:
- rotation
Rotationinstance Object containing the rotations represented by the rotation matrices.
- rotation
Notes
This function was called from_dcm before.
Added in version 1.4.0.
Array API Standard Support
from_matrixhas experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variableSCIPY_ARRAY_API=1and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.Library
CPU
GPU
NumPy
✅
n/a
CuPy
n/a
✅
PyTorch
✅
✅
JAX
✅
✅
Dask
⛔
n/a
See Support for the array API standard for more information.
References
[3]F. Landis Markley, “Unit Quaternion from Rotation Matrix”, Journal of guidance, control, and dynamics vol. 31.2, pp. 440-442, 2008.
Examples
>>> from scipy.spatial.transform import Rotation as R >>> import numpy as np
Initialize a single rotation:
>>> r = R.from_matrix([ ... [0, -1, 0], ... [1, 0, 0], ... [0, 0, 1]]) >>> r.single True >>> r.as_matrix().shape (3, 3)
Initialize multiple rotations in a single object:
>>> r = R.from_matrix([ ... [ ... [0, -1, 0], ... [1, 0, 0], ... [0, 0, 1], ... ], ... [ ... [1, 0, 0], ... [0, 0, -1], ... [0, 1, 0], ... ]]) >>> r.as_matrix().shape (2, 3, 3) >>> r.single False >>> len(r) 2
If input matrices are not special orthogonal (orthogonal with determinant equal to +1), then a special orthogonal estimate is stored:
>>> a = np.array([ ... [0, -0.5, 0], ... [0.5, 0, 0], ... [0, 0, 0.5]]) >>> np.linalg.det(a) 0.125 >>> r = R.from_matrix(a) >>> matrix = r.as_matrix() >>> matrix array([[ 0., -1., 0.], [ 1., 0., 0.], [ 0., 0., 1.]]) >>> np.linalg.det(matrix) 1.0
It is also possible to have a stack containing a single rotation:
>>> r = R.from_matrix([[ ... [0, -1, 0], ... [1, 0, 0], ... [0, 0, 1]]]) >>> r.as_matrix() array([[[ 0., -1., 0.], [ 1., 0., 0.], [ 0., 0., 1.]]]) >>> r.as_matrix().shape (1, 3, 3)
We can also create an N-dimensional array of rotations:
>>> r = R.from_matrix(np.tile(np.eye(3), (2, 3, 1, 1))) >>> r.shape (2, 3)