bundles / scipy latest / scipy / spatial / transform / _rotation / Rotation / from_matrix
staticmethod
scipy.spatial.transform._rotation:Rotation.from_matrix
Signature
def from_matrix ( matrix : ArrayLike , * , assume_valid : bool = False ) → Rotation Summary
Initialize from rotation matrix.
Extended Summary
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
matrix: array_like, shape (..., 3, 3)A single matrix or an ND array of matrices, where the last two dimensions contain the rotation matrices.
assume_valid: bool, optionalMust 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: `Rotation` instanceObject containing the rotations represented by the rotation matrices.
Notes
This function was called from_dcm before.
Array API Standard Support
from_matrix has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and 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
dev-arrayapifor more information.
Examples
from scipy.spatial.transform import Rotation as R import numpy as np✓
r = R.from_matrix([ [0, -1, 0], [1, 0, 0], [0, 0, 1]]) r.single r.as_matrix().shape✓
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 r.single len(r)✓
a = np.array([ [0, -0.5, 0], [0.5, 0, 0], [0, 0, 0.5]])✓
np.linalg.det(a)
✗r = R.from_matrix(a) matrix = r.as_matrix() matrix✓
np.linalg.det(matrix)
✗r = R.from_matrix([[ [0, -1, 0], [1, 0, 0], [0, 0, 1]]]) r.as_matrix() r.as_matrix().shape✓
r = R.from_matrix(np.tile(np.eye(3), (2, 3, 1, 1))) r.shape✓
Aliases
-
scipy.spatial.transform.Rotation.from_matrix