trimesh.transformations module

Homogeneous Transformation Matrices and Quaternions.

A library for calculating 4x4 matrices for translating, rotating, reflecting, scaling, shearing, projecting, orthogonalizing, and superimposing arrays of 3D homogeneous coordinates as well as for converting between rotation matrices, Euler angles, and quaternions. Also includes an Arcball control object and functions to decompose transformation matrices.

Author:

Christoph Gohlke

Organization:

Laboratory for Fluorescence Dynamics, University of California, Irvine

Version:

2017.02.17

Requirements

Notes

The API is not stable yet and is expected to change between revisions.

This Python code is not optimized for speed. Refer to the transformations.c module for a faster implementation of some functions.

Documentation in HTML format can be generated with epydoc.

Matrices (M) can be inverted using np.linalg.inv(M), be concatenated using np.dot(M0, M1), or transform homogeneous coordinate arrays (v) using np.dot(M, v) for shape (4, ) column vectors, respectively np.dot(v, M.T) for shape (, 4) row vectors (“array of points”).

This module follows the “column vectors on the right” and “row major storage” (C contiguous) conventions. The translation components are in the right column of the transformation matrix, i.e. M[:3, 3]. The transpose of the transformation matrices may have to be used to interface with other graphics systems, e.g. with OpenGL’s glMultMatrixd(). See also [16].

Calculations are carried out with np.float64 precision.

Vector, point, quaternion, and matrix function arguments are expected to be “array like”, i.e. tuple, list, or numpy arrays.

Return types are numpy arrays unless specified otherwise.

Angles are in radians unless specified otherwise.

Quaternions w+ix+jy+kz are represented as [w, x, y, z].

A triple of Euler angles can be applied/interpreted in 24 ways, which can be specified using a 4 character string or encoded 4-tuple:

Axes 4-string: e.g. ‘sxyz’ or ‘ryxy’

  • first character : rotations are applied to ‘s’tatic or ‘r’otating frame

  • remaining characters : successive rotation axis ‘x’, ‘y’, or ‘z’

Axes 4-tuple: e.g. (0, 0, 0, 0) or (1, 1, 1, 1)

  • inner axis: code of axis (‘x’:0, ‘y’:1, ‘z’:2) of rightmost matrix.

  • parity : even (0) if inner axis ‘x’ is followed by ‘y’, ‘y’ is followed by ‘z’, or ‘z’ is followed by ‘x’. Otherwise odd (1).

  • repetition : first and last axis are same (1) or different (0).

  • frame : rotations are applied to static (0) or rotating (1) frame.

Other Python packages and modules for 3D transformations and quaternions:

References

  1. Matrices and transformations. Ronald Goldman. In “Graphics Gems I”, pp 472-475. Morgan Kaufmann, 1990.

  2. More matrices and transformations: shear and pseudo-perspective. Ronald Goldman. In “Graphics Gems II”, pp 320-323. Morgan Kaufmann, 1991.

  3. Decomposing a matrix into simple transformations. Spencer Thomas. In “Graphics Gems II”, pp 320-323. Morgan Kaufmann, 1991.

  4. Recovering the data from the transformation matrix. Ronald Goldman. In “Graphics Gems II”, pp 324-331. Morgan Kaufmann, 1991.

  5. Euler angle conversion. Ken Shoemake. In “Graphics Gems IV”, pp 222-229. Morgan Kaufmann, 1994.

  6. Arcball rotation control. Ken Shoemake. In “Graphics Gems IV”, pp 175-192. Morgan Kaufmann, 1994.

  7. Representing attitude: Euler angles, unit quaternions, and rotation vectors. James Diebel. 2006.

  8. A discussion of the solution for the best rotation to relate two sets of vectors. W Kabsch. Acta Cryst. 1978. A34, 827-828.

  9. Closed-form solution of absolute orientation using unit quaternions. BKP Horn. J Opt Soc Am A. 1987. 4(4):629-642.

  10. Quaternions. Ken Shoemake. http://www.sfu.ca/~jwa3/cmpt461/files/quatut.pdf

  11. From quaternion to matrix and back. JMP van Waveren. 2005. http://www.intel.com/cd/ids/developer/asmo-na/eng/293748.htm

  12. Uniform random rotations. Ken Shoemake. In “Graphics Gems III”, pp 124-132. Morgan Kaufmann, 1992.

  13. Quaternion in molecular modeling. CFF Karney. J Mol Graph Mod, 25(5):595-604

  14. New method for extracting the quaternion from a rotation matrix. Itzhack Y Bar-Itzhack, J Guid Contr Dynam. 2000. 23(6): 1085-1087.

  15. Multiple View Geometry in Computer Vision. Hartley and Zissermann. Cambridge University Press; 2nd Ed. 2004. Chapter 4, Algorithm 4.7, p 130.

  16. Column Vectors vs. Row Vectors. http://steve.hollasch.net/cgindex/math/matrix/column-vec.html

Examples

>>> alpha, beta, gamma = 0.123, -1.234, 2.345
>>> origin, xaxis, yaxis, zaxis = [0, 0, 0], [1, 0, 0], [0, 1, 0], [0, 0, 1]
>>> I = identity_matrix()
>>> Rx = rotation_matrix(alpha, xaxis)
>>> Ry = rotation_matrix(beta, yaxis)
>>> Rz = rotation_matrix(gamma, zaxis)
>>> R = concatenate_matrices(Rx, Ry, Rz)
>>> euler = euler_from_matrix(R, 'rxyz')
>>> np.allclose([alpha, beta, gamma], euler)
True
>>> Re = euler_matrix(alpha, beta, gamma, 'rxyz')
>>> is_same_transform(R, Re)
True
>>> al, be, ga = euler_from_matrix(Re, 'rxyz')
>>> is_same_transform(Re, euler_matrix(al, be, ga, 'rxyz'))
True
>>> qx = quaternion_about_axis(alpha, xaxis)
>>> qy = quaternion_about_axis(beta, yaxis)
>>> qz = quaternion_about_axis(gamma, zaxis)
>>> q = quaternion_multiply(qx, qy)
>>> q = quaternion_multiply(q, qz)
>>> Rq = quaternion_matrix(q)
>>> is_same_transform(R, Rq)
True
>>> S = scale_matrix(1.23, origin)
>>> T = translation_matrix([1, 2, 3])
>>> Z = shear_matrix(beta, xaxis, origin, zaxis)
>>> R = random_rotation_matrix(np.random.rand(3))
>>> M = concatenate_matrices(T, R, Z, S)
>>> scale, shear, angles, trans, persp = decompose_matrix(M)
>>> np.allclose(scale, 1.23)
True
>>> np.allclose(trans, [1, 2, 3])
True
>>> np.allclose(shear, [0, np.tan(beta), 0])
True
>>> is_same_transform(R, euler_matrix(axes='sxyz', *angles))
True
>>> M1 = compose_matrix(scale, shear, angles, trans, persp)
>>> is_same_transform(M, M1)
True
>>> v0, v1 = random_vector(3), random_vector(3)
>>> M = rotation_matrix(angle_between_vectors(v0, v1), vector_product(v0, v1))
>>> v2 = np.dot(v0, M[:3,:3].T)
>>> np.allclose(unit_vector(v1), unit_vector(v2))
True
class trimesh.transformations.Arcball(initial=None)

Bases: object

Virtual Trackball Control.

>>> ball = Arcball()
>>> ball = Arcball(initial=np.identity(4))
>>> ball.place([320, 320], 320)
>>> ball.down([500, 250])
>>> ball.drag([475, 275])
>>> R = ball.matrix()
>>> np.allclose(np.sum(R), 3.90583455)
True
>>> ball = Arcball(initial=[1, 0, 0, 0])
>>> ball.place([320, 320], 320)
>>> ball.setaxes([1, 1, 0], [-1, 1, 0])
>>> ball.constrain = True
>>> ball.down([400, 200])
>>> ball.drag([200, 400])
>>> R = ball.matrix()
>>> np.allclose(np.sum(R), 0.2055924)
True
>>> ball.next()
__init__(initial=None)

Initialize virtual trackball control.

initial : quaternion or rotation matrix

property constrain

Return state of constrain to axis mode.

down(point)

Set initial cursor window coordinates and pick constrain-axis.

drag(point)

Update current cursor window coordinates.

matrix()

Return homogeneous rotation matrix.

next(acceleration=0.0)

Continue rotation in direction of last drag.

place(center, radius)

Place Arcball, e.g. when window size changes.

centersequence[2]

Window coordinates of trackball center.

radiusfloat

Radius of trackball in window coordinates.

setaxes(*axes)

Set axes to constrain rotations.

trimesh.transformations.affine_matrix_from_points(v0, v1, shear=True, scale=True, usesvd=True)

Return affine transform matrix to register two point sets.

v0 and v1 are shape (ndims, *) arrays of at least ndims non-homogeneous coordinates, where ndims is the dimensionality of the coordinate space.

If shear is False, a similarity transformation matrix is returned. If also scale is False, a rigid/Euclidean transformation matrix is returned.

By default the algorithm by Hartley and Zissermann [15] is used. If usesvd is True, similarity and Euclidean transformation matrices are calculated by minimizing the weighted sum of squared deviations (RMSD) according to the algorithm by Kabsch [8]. Otherwise, and if ndims is 3, the quaternion based algorithm by Horn [9] is used, which is slower when using this Python implementation.

The returned matrix performs rotation, translation and uniform scaling (if specified).

>>> v0 = [[0, 1031, 1031, 0], [0, 0, 1600, 1600]]
>>> v1 = [[675, 826, 826, 677], [55, 52, 281, 277]]
>>> mat = affine_matrix_from_points(v0, v1)
>>> T = translation_matrix(np.random.random(3)-0.5)
>>> R = random_rotation_matrix(np.random.random(3))
>>> S = scale_matrix(random.random())
>>> M = concatenate_matrices(T, R, S)
>>> v0 = (np.random.rand(4, 100) - 0.5) * 20
>>> v0[3] = 1
>>> v1 = np.dot(M, v0)
>>> v0[:3] += np.random.normal(0, 1e-8, 300).reshape(3, -1)
>>> M = affine_matrix_from_points(v0[:3], v1[:3])
>>> check = np.allclose(v1, np.dot(M, v0))

More examples in superimposition_matrix()

trimesh.transformations.angle_between_vectors(v0, v1, directed=True, axis=0)

Return angle between vectors.

If directed is False, the input vectors are interpreted as undirected axes, i.e. the maximum angle is pi/2.

>>> a = angle_between_vectors([1, -2, 3], [-1, 2, -3])
>>> np.allclose(a, np.pi)
True
>>> a = angle_between_vectors([1, -2, 3], [-1, 2, -3], directed=False)
>>> np.allclose(a, 0)
True
>>> v0 = [[2, 0, 0, 2], [0, 2, 0, 2], [0, 0, 2, 2]]
>>> v1 = [[3], [0], [0]]
>>> a = angle_between_vectors(v0, v1)
>>> np.allclose(a, [0, 1.5708, 1.5708, 0.95532])
True
>>> v0 = [[2, 0, 0], [2, 0, 0], [0, 2, 0], [2, 0, 0]]
>>> v1 = [[0, 3, 0], [0, 0, 3], [0, 0, 3], [3, 3, 3]]
>>> a = angle_between_vectors(v0, v1, axis=1)
>>> np.allclose(a, [1.5708, 1.5708, 1.5708, 0.95532])
True
trimesh.transformations.arcball_constrain_to_axis(point, axis)

Return sphere point perpendicular to axis.

trimesh.transformations.arcball_map_to_sphere(point, center, radius)

Return unit sphere coordinates from window coordinates.

trimesh.transformations.arcball_nearest_axis(point, axes)

Return axis, which arc is nearest to point.

trimesh.transformations.clip_matrix(left, right, bottom, top, near, far, perspective=False)

Return matrix to obtain normalized device coordinates from frustum.

The frustum bounds are axis-aligned along x (left, right), y (bottom, top) and z (near, far).

Normalized device coordinates are in range [-1, 1] if coordinates are inside the frustum.

If perspective is True the frustum is a truncated pyramid with the perspective point at origin and direction along z axis, otherwise an orthographic canonical view volume (a box).

Homogeneous coordinates transformed by the perspective clip matrix need to be dehomogenized (divided by w coordinate).

>>> frustum = np.random.rand(6)
>>> frustum[1] += frustum[0]
>>> frustum[3] += frustum[2]
>>> frustum[5] += frustum[4]
>>> M = clip_matrix(perspective=False, *frustum)
>>> a = np.dot(M, [frustum[0], frustum[2], frustum[4], 1])
>>> np.allclose(a, [-1., -1., -1.,  1.])
True
>>> b = np.dot(M, [frustum[1], frustum[3], frustum[5], 1])
>>> np.allclose(b, [ 1.,  1.,  1.,  1.])
True
>>> M = clip_matrix(perspective=True, *frustum)
>>> v = np.dot(M, [frustum[0], frustum[2], frustum[4], 1])
>>> c = v / v[3]
>>> np.allclose(c, [-1., -1., -1.,  1.])
True
>>> v = np.dot(M, [frustum[1], frustum[3], frustum[4], 1])
>>> d = v / v[3]
>>> np.allclose(d, [ 1.,  1., -1.,  1.])
True
trimesh.transformations.compose_matrix(scale=None, shear=None, angles=None, translate=None, perspective=None)

Return transformation matrix from sequence of transformations.

This is the inverse of the decompose_matrix function.

Sequence of transformations:

scale : vector of 3 scaling factors shear : list of shear factors for x-y, x-z, y-z axes angles : list of Euler angles about static x, y, z axes translate : translation vector along x, y, z axes perspective : perspective partition of matrix

>>> scale = np.random.random(3) - 0.5
>>> shear = np.random.random(3) - 0.5
>>> angles = (np.random.random(3) - 0.5) * (2*np.pi)
>>> trans = np.random.random(3) - 0.5
>>> persp = np.random.random(4) - 0.5
>>> M0 = compose_matrix(scale, shear, angles, trans, persp)
>>> result = decompose_matrix(M0)
>>> M1 = compose_matrix(*result)
>>> is_same_transform(M0, M1)
True
trimesh.transformations.concatenate_matrices(*matrices)

Return concatenation of series of transformation matrices.

>>> M = np.random.rand(16).reshape((4, 4)) - 0.5
>>> np.allclose(M, concatenate_matrices(M))
True
>>> np.allclose(np.dot(M, M.T), concatenate_matrices(M, M.T))
True
trimesh.transformations.decompose_matrix(matrix)

Return sequence of transformations from transformation matrix.

matrixarray_like

Non-degenerative homogeneous transformation matrix

Return tuple of:

scale : vector of 3 scaling factors shear : list of shear factors for x-y, x-z, y-z axes angles : list of Euler angles about static x, y, z axes translate : translation vector along x, y, z axes perspective : perspective partition of matrix

Raise ValueError if matrix is of wrong type or degenerative.

>>> T0 = translation_matrix([1, 2, 3])
>>> scale, shear, angles, trans, persp = decompose_matrix(T0)
>>> T1 = translation_matrix(trans)
>>> np.allclose(T0, T1)
True
>>> S = scale_matrix(0.123)
>>> scale, shear, angles, trans, persp = decompose_matrix(S)
>>> bool(np.isclose(scale[0], 0.123))
True
>>> R0 = euler_matrix(1, 2, 3)
>>> scale, shear, angles, trans, persp = decompose_matrix(R0)
>>> R1 = euler_matrix(*angles)
>>> np.allclose(R0, R1)
True
trimesh.transformations.euler_from_matrix(matrix, axes='sxyz')

Return Euler angles from rotation matrix for specified axis sequence.

axes : One of 24 axis sequences as string or encoded tuple

Note that many Euler angle triplets can describe one matrix.

>>> R0 = euler_matrix(1, 2, 3, 'syxz')
>>> al, be, ga = euler_from_matrix(R0, 'syxz')
>>> R1 = euler_matrix(al, be, ga, 'syxz')
>>> np.allclose(R0, R1)
True
>>> angles = (4*np.pi) * (np.random.random(3) - 0.5)
>>> for axes in _AXES2TUPLE.keys():
...    R0 = euler_matrix(axes=axes, *angles)
...    R1 = euler_matrix(axes=axes, *euler_from_matrix(R0, axes))
...    if not np.allclose(R0, R1): print(axes, "failed")
trimesh.transformations.euler_from_quaternion(quaternion, axes='sxyz')

Return Euler angles from quaternion for specified axis sequence.

>>> angles = euler_from_quaternion([0.99810947, 0.06146124, 0, 0])
>>> np.allclose(angles, [0.123, 0, 0])
True
trimesh.transformations.euler_matrix(ai, aj, ak, axes='sxyz')

Return homogeneous rotation matrix from Euler angles and axis sequence.

ai, aj, ak : Euler’s roll, pitch and yaw angles axes : One of 24 axis sequences as string or encoded tuple

>>> R = euler_matrix(1, 2, 3, 'syxz')
>>> np.allclose(np.sum(R[0]), -1.34786452)
True
>>> R = euler_matrix(1, 2, 3, (0, 1, 0, 1))
>>> np.allclose(np.sum(R[0]), -0.383436184)
True
>>> ai, aj, ak = (4*np.pi) * (np.random.random(3) - 0.5)
>>> for axes in _AXES2TUPLE.keys():
...    R = euler_matrix(ai, aj, ak, axes)
>>> for axes in _TUPLE2AXES.keys():
...    R = euler_matrix(ai, aj, ak, axes)
trimesh.transformations.fix_rigid(matrix, max_deviance=1e-05)

If a homogeneous transformation matrix is almost a rigid transform but many matrix-multiplies have accumulated some floating point error try to restore the matrix using SVD.

Parameters:
  • matrix ((4, 4) or (3, 3) float) – Homogeneous transformation matrix.

  • max_deviance (float) – Do not alter the matrix if it is not rigid by more than this amount.

Returns:

repaired – Repaired homogeneous transformation matrix

Return type:

(4, 4) or (3, 3) float

trimesh.transformations.flips_winding(matrix)

Check to see if a matrix will invert triangles.

Parameters:

matrix ((4, 4) float) – Homogeneous transformation matrix

Returns:

flip – True if matrix will flip winding of triangles.

Return type:

bool

trimesh.transformations.identity_matrix()

Return 4x4 identity/unit matrix.

>>> I = identity_matrix()
>>> np.allclose(I, np.dot(I, I))
True
>>> float(np.sum(I)), float(np.trace(I))
(4.0, 4.0)
>>> np.allclose(I, np.identity(4))
True
trimesh.transformations.inverse_matrix(matrix)

Return inverse of square transformation matrix.

>>> M0 = random_rotation_matrix()
>>> M1 = inverse_matrix(M0.T)
>>> np.allclose(M1, np.linalg.inv(M0.T))
True
>>> for size in range(1, 7):
...     M0 = np.random.rand(size, size)
...     M1 = inverse_matrix(M0)
...     if not np.allclose(M1, np.linalg.inv(M0)): print(size)
trimesh.transformations.is_rigid(matrix, epsilon=1e-08)

Check to make sure a homogeonous transformation matrix is a rigid transform.

Parameters:

matrix ((4, 4) float) – A transformation matrix

Returns:

check – True if matrix is a a transform with only translation, scale, and rotation

Return type:

bool

trimesh.transformations.is_same_quaternion(q0, q1)

Return True if two quaternions are equal.

trimesh.transformations.is_same_transform(matrix0, matrix1)

Return True if two matrices perform same transformation.

>>> is_same_transform(np.identity(4), np.identity(4))
True
>>> is_same_transform(np.identity(4), random_rotation_matrix())
False
trimesh.transformations.orthogonalization_matrix(lengths, angles)

Return orthogonalization matrix for crystallographic cell coordinates.

Angles are expected in degrees.

The de-orthogonalization matrix is the inverse.

>>> O = orthogonalization_matrix([10, 10, 10], [90, 90, 90])
>>> np.allclose(O[:3, :3], np.identity(3, float) * 10)
True
>>> O = orthogonalization_matrix([9.8, 12.0, 15.5], [87.2, 80.7, 69.7])
>>> np.allclose(np.sum(O), 43.063229)
True
trimesh.transformations.planar_matrix(offset=None, theta=None, point=None, scale=None)

2D homogeonous transformation matrix.

Parameters:
  • offset ((2,) float) – XY offset

  • theta (float) – Rotation around Z in radians

  • point ((2, ) float) – Point to rotate around

  • scale ((2,) float or None) – Scale to apply

Returns:

matrix – Homogeneous 2D transformation matrix

Return type:

(3, 3) flat

trimesh.transformations.planar_matrix_to_3D(matrix_2D)

Given a 2D homogeneous rotation matrix convert it to a 3D rotation matrix that is rotating around the Z axis

Parameters:

matrix_2D ((3,3) float, homogeneous 2D rotation matrix)

Returns:

matrix_3D

Return type:

(4,4) float, homogeneous 3D rotation matrix

trimesh.transformations.projection_from_matrix(matrix, pseudo=False)

Return projection plane and perspective point from projection matrix.

Return values are same as arguments for projection_matrix function: point, normal, direction, perspective, and pseudo.

>>> point = np.random.random(3) - 0.5
>>> normal = np.random.random(3) - 0.5
>>> direct = np.random.random(3) - 0.5
>>> persp = np.random.random(3) - 0.5
>>> P0 = projection_matrix(point, normal)
>>> result = projection_from_matrix(P0)
>>> P1 = projection_matrix(*result)
>>> is_same_transform(P0, P1)
True
>>> P0 = projection_matrix(point, normal, direct)
>>> result = projection_from_matrix(P0)
>>> P1 = projection_matrix(*result)
>>> is_same_transform(P0, P1)
True
>>> P0 = projection_matrix(point, normal, perspective=persp, pseudo=False)
>>> result = projection_from_matrix(P0, pseudo=False)
>>> P1 = projection_matrix(*result)
>>> is_same_transform(P0, P1)
True
>>> P0 = projection_matrix(point, normal, perspective=persp, pseudo=True)
>>> result = projection_from_matrix(P0, pseudo=True)
>>> P1 = projection_matrix(*result)
>>> is_same_transform(P0, P1)
True
trimesh.transformations.projection_matrix(point, normal, direction=None, perspective=None, pseudo=False)

Return matrix to project onto plane defined by point and normal.

Using either perspective point, projection direction, or none of both.

If pseudo is True, perspective projections will preserve relative depth such that Perspective = dot(Orthogonal, PseudoPerspective).

>>> P = projection_matrix([0, 0, 0], [1, 0, 0])
>>> np.allclose(P[1:, 1:], np.identity(4)[1:, 1:])
True
>>> point = np.random.random(3) - 0.5
>>> normal = np.random.random(3) - 0.5
>>> direct = np.random.random(3) - 0.5
>>> persp = np.random.random(3) - 0.5
>>> P0 = projection_matrix(point, normal)
>>> P1 = projection_matrix(point, normal, direction=direct)
>>> P2 = projection_matrix(point, normal, perspective=persp)
>>> P3 = projection_matrix(point, normal, perspective=persp, pseudo=True)
>>> is_same_transform(P2, np.dot(P0, P3))
True
>>> P = projection_matrix([3, 0, 0], [1, 1, 0], [1, 0, 0])
>>> v0 = (np.random.rand(4, 5) - 0.5) * 20
>>> v0[3] = 1
>>> v1 = np.dot(P, v0)
>>> np.allclose(v1[1], v0[1])
True
>>> np.allclose(v1[0], 3-v1[1])
True
trimesh.transformations.quaternion_about_axis(angle, axis)

Return quaternion for rotation about axis.

>>> q = quaternion_about_axis(0.123, [1, 0, 0])
>>> np.allclose(q, [0.99810947, 0.06146124, 0, 0])
True
trimesh.transformations.quaternion_conjugate(quaternion)

Return conjugate of quaternion.

>>> q0 = random_quaternion()
>>> q1 = quaternion_conjugate(q0)
>>> q1[0] == q0[0] and all(q1[1:] == -q0[1:])
True
trimesh.transformations.quaternion_from_euler(ai, aj, ak, axes='sxyz')

Return quaternion from Euler angles and axis sequence.

ai, aj, ak : Euler’s roll, pitch and yaw angles axes : One of 24 axis sequences as string or encoded tuple

>>> q = quaternion_from_euler(1, 2, 3, 'ryxz')
>>> np.allclose(q, [0.435953, 0.310622, -0.718287, 0.444435])
True
trimesh.transformations.quaternion_from_matrix(matrix, isprecise=False)

Return quaternion from rotation matrix.

If isprecise is True, the input matrix is assumed to be a precise rotation matrix and a faster algorithm is used.

>>> q = quaternion_from_matrix(np.identity(4), True)
>>> np.allclose(q, [1, 0, 0, 0])
True
>>> q = quaternion_from_matrix(np.diag([1, -1, -1, 1]))
>>> np.allclose(q, [0, 1, 0, 0]) or np.allclose(q, [0, -1, 0, 0])
True
>>> R = rotation_matrix(0.123, (1, 2, 3))
>>> q = quaternion_from_matrix(R, True)
>>> np.allclose(q, [0.9981095, 0.0164262, 0.0328524, 0.0492786])
True
>>> R = [[-0.545, 0.797, 0.260, 0], [0.733, 0.603, -0.313, 0],
...      [-0.407, 0.021, -0.913, 0], [0, 0, 0, 1]]
>>> q = quaternion_from_matrix(R)
>>> np.allclose(q, [0.19069, 0.43736, 0.87485, -0.083611])
True
>>> R = [[0.395, 0.362, 0.843, 0], [-0.626, 0.796, -0.056, 0],
...      [-0.677, -0.498, 0.529, 0], [0, 0, 0, 1]]
>>> q = quaternion_from_matrix(R)
>>> np.allclose(q, [0.82336615, -0.13610694, 0.46344705, -0.29792603])
True
>>> R = random_rotation_matrix()
>>> q = quaternion_from_matrix(R)
>>> is_same_transform(R, quaternion_matrix(q))
True
>>> is_same_quaternion(quaternion_from_matrix(R, isprecise=False),
...                    quaternion_from_matrix(R, isprecise=True))
True
>>> R = euler_matrix(0.0, 0.0, np.pi/2.0)
>>> is_same_quaternion(quaternion_from_matrix(R, isprecise=False),
...                    quaternion_from_matrix(R, isprecise=True))
True
trimesh.transformations.quaternion_imag(quaternion)

Return imaginary part of quaternion.

>>> quaternion_imag([3, 0, 1, 2])
array([0., 1., 2.])
trimesh.transformations.quaternion_inverse(quaternion)

Return inverse of quaternion.

>>> q0 = random_quaternion()
>>> q1 = quaternion_inverse(q0)
>>> np.allclose(quaternion_multiply(q0, q1), [1, 0, 0, 0])
True
trimesh.transformations.quaternion_matrix(quaternion)

Return a homogeneous rotation matrix from quaternion.

>>> M = quaternion_matrix([0.99810947, 0.06146124, 0, 0])
>>> np.allclose(M, rotation_matrix(0.123, [1, 0, 0]))
True
>>> M = quaternion_matrix([1, 0, 0, 0])
>>> np.allclose(M, np.identity(4))
True
>>> M = quaternion_matrix([0, 1, 0, 0])
>>> np.allclose(M, np.diag([1, -1, -1, 1]))
True
>>> M = quaternion_matrix([[1, 0, 0, 0],[0, 1, 0, 0]])
>>> np.allclose(M, np.array([np.identity(4), np.diag([1, -1, -1, 1])]))
True
trimesh.transformations.quaternion_multiply(quaternion1, quaternion0)

Return multiplication of two quaternions.

>>> q = quaternion_multiply([4, 1, -2, 3], [8, -5, 6, 7])
>>> np.allclose(q, [28, -44, -14, 48])
True
trimesh.transformations.quaternion_real(quaternion)

Return real part of quaternion.

>>> quaternion_real([3, 0, 1, 2])
3.0
trimesh.transformations.quaternion_slerp(quat0, quat1, fraction, spin=0, shortestpath=True)

Return spherical linear interpolation between two quaternions.

>>> q0 = random_quaternion()
>>> q1 = random_quaternion()
>>> q = quaternion_slerp(q0, q1, 0)
>>> np.allclose(q, q0)
True
>>> q = quaternion_slerp(q0, q1, 1, 1)
>>> np.allclose(q, q1)
True
>>> q = quaternion_slerp(q0, q1, 0.5)
>>> angle = np.arccos(np.dot(q0, q))
>>> np.allclose(2, np.arccos(np.dot(q0, q1)) / angle) or         np.allclose(2, np.arccos(-np.dot(q0, q1)) / angle)
True
trimesh.transformations.random_quaternion(rand=None, num=1)

Return uniform random unit quaternion.

rand: array like or None

Three independent random variables that are uniformly distributed between 0 and 1.

>>> q = random_quaternion()
>>> np.allclose(1, vector_norm(q))
True
>>> q = random_quaternion(num=10)
>>> np.allclose(1, vector_norm(q, axis=1))
True
>>> q = random_quaternion(np.random.random(3))
>>> len(q.shape), q.shape[0]==4
(1, True)
trimesh.transformations.random_rotation_matrix(rand=None, num=1, translate=False)

Return uniform random rotation matrix.

rand: array like

Three independent random variables that are uniformly distributed between 0 and 1 for each returned quaternion.

>>> R = random_rotation_matrix()
>>> np.allclose(np.dot(R.T, R), np.identity(4))
True
>>> R = random_rotation_matrix(num=10)
>>> np.allclose(np.einsum('...ji,...jk->...ik', R, R), np.identity(4))
True
trimesh.transformations.random_vector(size)

Return array of random doubles in the half-open interval [0.0, 1.0).

>>> v = random_vector(10000)
>>> bool(np.all(v >= 0) and np.all(v < 1))
True
>>> v0 = random_vector(10)
>>> v1 = random_vector(10)
>>> bool(np.any(v0 == v1))
False
trimesh.transformations.reflection_from_matrix(matrix)

Return mirror plane point and normal vector from reflection matrix.

>>> v0 = np.random.random(3) - 0.5
>>> v1 = np.random.random(3) - 0.5
>>> M0 = reflection_matrix(v0, v1)
>>> point, normal = reflection_from_matrix(M0)
>>> M1 = reflection_matrix(point, normal)
>>> is_same_transform(M0, M1)
True
trimesh.transformations.reflection_matrix(point, normal)

Return matrix to mirror at plane defined by point and normal vector.

>>> v0 = np.random.random(4) - 0.5
>>> v0[3] = 1.
>>> v1 = np.random.random(3) - 0.5
>>> R = reflection_matrix(v0, v1)
>>> np.allclose(2, np.trace(R))
True
>>> np.allclose(v0, np.dot(R, v0))
True
>>> v2 = v0.copy()
>>> v2[:3] += v1
>>> v3 = v0.copy()
>>> v2[:3] -= v1
>>> np.allclose(v2, np.dot(R, v3))
True
trimesh.transformations.rotation_from_matrix(matrix)

Return rotation angle and axis from rotation matrix.

>>> angle = (random.random() - 0.5) * (2*np.pi)
>>> direc = np.random.random(3) - 0.5
>>> point = np.random.random(3) - 0.5
>>> R0 = rotation_matrix(angle, direc, point)
>>> angle, direc, point = rotation_from_matrix(R0)
>>> R1 = rotation_matrix(angle, direc, point)
>>> is_same_transform(R0, R1)
True
trimesh.transformations.rotation_matrix(angle, direction, point=None)

Return matrix to rotate about axis defined by point and direction.

Parameters:
  • angle (float, or sympy.Symbol) – Angle, in radians or symbolic angle

  • direction ((3,) float) – Any vector along rotation axis

  • point ((3, ) float, or None) – Origin point of rotation axis

Returns:

matrix – Homogeneous transformation matrix

Return type:

(4, 4) float, or (4, 4) sympy.Matrix

Examples

>>> R = rotation_matrix(np.pi/2, [0, 0, 1], [1, 0, 0])
>>> np.allclose(np.dot(R, [0, 0, 0, 1]), [1, -1, 0, 1])
True
>>> angle = (random.random() - 0.5) * (2*np.pi)
>>> direc = np.random.random(3) - 0.5
>>> point = np.random.random(3) - 0.5
>>> R0 = rotation_matrix(angle, direc, point)
>>> R1 = rotation_matrix(angle-2*np.pi, direc, point)
>>> is_same_transform(R0, R1)
True
>>> R0 = rotation_matrix(angle, direc, point)
>>> R1 = rotation_matrix(-angle, -direc, point)
>>> is_same_transform(R0, R1)
True
>>> I = np.identity(4, np.float64)
>>> np.allclose(I, rotation_matrix(np.pi*2, direc))
True
>>> np.allclose(2, np.trace(rotation_matrix(np.pi/2,direc,point)))
True
trimesh.transformations.scale_and_translate(scale=None, translate=None)

Optimized version of compose_matrix for just scaling then translating.

Scalar args are broadcast to arrays of shape (3,)

Parameters:
  • scale (float or (3,) float) – Scale factor

  • translate (float or (3,) float) – Translation

trimesh.transformations.scale_from_matrix(matrix)

Return scaling factor, origin and direction from scaling matrix.

>>> factor = random.random() * 10 - 5
>>> origin = np.random.random(3) - 0.5
>>> direct = np.random.random(3) - 0.5
>>> S0 = scale_matrix(factor, origin)
>>> factor, origin, direction = scale_from_matrix(S0)
>>> S1 = scale_matrix(factor, origin, direction)
>>> is_same_transform(S0, S1)
True
>>> S0 = scale_matrix(factor, origin, direct)
>>> factor, origin, direction = scale_from_matrix(S0)
>>> S1 = scale_matrix(factor, origin, direction)
>>> is_same_transform(S0, S1)
True
trimesh.transformations.scale_matrix(factor, origin=None, direction=None)

Return matrix to scale by factor around origin in direction.

Use factor -1 for point symmetry.

>>> v = (np.random.rand(4, 5) - 0.5) * 20
>>> v[3] = 1
>>> S = scale_matrix(-1.234)
>>> np.allclose(np.dot(S, v)[:3], -1.234*v[:3])
True
>>> factor = random.random() * 10 - 5
>>> origin = np.random.random(3) - 0.5
>>> direct = np.random.random(3) - 0.5
>>> S = scale_matrix(factor, origin)
>>> S = scale_matrix(factor, origin, direct)
trimesh.transformations.shear_from_matrix(matrix)

Return shear angle, direction and plane from shear matrix.

>>> angle  = np.pi / 2.0
>>> direct = [0.0, 1.0, 0.0]
>>> point  = [0.0, 0.0, 0.0]
>>> normal = np.cross(direct, np.roll(direct,1))
>>> S0 = shear_matrix(angle, direct, point, normal)
>>> angle, direct, point, normal = shear_from_matrix(S0)
>>> S1 = shear_matrix(angle, direct, point, normal)
>>> is_same_transform(S0, S1)
True
trimesh.transformations.shear_matrix(angle, direction, point, normal)

Return matrix to shear by angle along direction vector on shear plane.

The shear plane is defined by a point and normal vector. The direction vector must be orthogonal to the plane’s normal vector.

A point P is transformed by the shear matrix into P” such that the vector P-P” is parallel to the direction vector and its extent is given by the angle of P-P’-P”, where P’ is the orthogonal projection of P onto the shear plane.

>>> angle = (random.random() - 0.5) * 4*np.pi
>>> direct = np.random.random(3) - 0.5
>>> point = np.random.random(3) - 0.5
>>> normal = np.cross(direct, np.random.random(3))
>>> S = shear_matrix(angle, direct, point, normal)
>>> np.allclose(1, np.linalg.det(S))
True
trimesh.transformations.spherical_matrix(theta, phi, axes='sxyz')

Give a spherical coordinate vector, find the rotation that will transform a [0,0,1] vector to those coordinates

Parameters:
  • theta (float, rotation angle in radians)

  • phi (float, rotation angle in radians)

Returns:

matrix – be a cartesian vector in the direction of the input spherical coordinates:

np.dot(matrix, [0,0,1,0])

Return type:

(4,4) rotation matrix where the following will

trimesh.transformations.superimposition_matrix(v0, v1, scale=False, usesvd=True)

Return matrix to transform given 3D point set into second point set.

v0 and v1 are shape (3, *) or (4, *) arrays of at least 3 points.

The parameters scale and usesvd are explained in the more general affine_matrix_from_points function.

The returned matrix is a similarity or Euclidean transformation matrix. This function has a fast C implementation in transformations.c.

>>> v0 = np.random.rand(3, 10)
>>> M = superimposition_matrix(v0, v0)
>>> np.allclose(M, np.identity(4))
True
>>> R = random_rotation_matrix(np.random.random(3))
>>> v0 = [[1,0,0], [0,1,0], [0,0,1], [1,1,1]]
>>> v1 = np.dot(R, v0)
>>> M = superimposition_matrix(v0, v1)
>>> np.allclose(v1, np.dot(M, v0))
True
>>> v0 = (np.random.rand(4, 100) - 0.5) * 20
>>> v0[3] = 1
>>> v1 = np.dot(R, v0)
>>> M = superimposition_matrix(v0, v1)
>>> np.allclose(v1, np.dot(M, v0))
True
>>> S = scale_matrix(random.random())
>>> T = translation_matrix(np.random.random(3)-0.5)
>>> M = concatenate_matrices(T, R, S)
>>> v1 = np.dot(M, v0)
>>> v0[:3] += np.random.normal(0, 1e-9, 300).reshape(3, -1)
>>> M = superimposition_matrix(v0, v1, scale=True)
>>> np.allclose(v1, np.dot(M, v0))
True
>>> M = superimposition_matrix(v0, v1, scale=True, usesvd=False)
>>> np.allclose(v1, np.dot(M, v0))
True
>>> v = np.zeros((4, 100, 3))
>>> v[:, :, 0] = v0
>>> M = superimposition_matrix(v0, v1, scale=True, usesvd=False)
>>> np.allclose(v1, np.dot(M, v[:, :, 0]))
True
trimesh.transformations.transform_around(matrix, point)

Given a transformation matrix, apply its rotation around a point in space.

Parameters:
  • matrix ((4,4) or (3, 3) float, transformation matrix)

  • point ((3,) or (2,) float, point in space)

Returns:

result

Return type:

(4,4) transformation matrix

trimesh.transformations.transform_points(points: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], matrix: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], translate: bool = True) ndarray[tuple[int, ...], dtype[float64]]

Returns points rotated by a homogeneous transformation matrix.

If points are (n, 2) matrix must be (3, 3) If points are (n, 3) matrix must be (4, 4)

Parameters:
  • points ((n, dim) float) – Points where dim is 2 or 3.

  • matrix ((3, 3) or (4, 4) float) – Homogeneous rotation matrix.

  • translate (bool) – Apply translation from matrix or not.

Returns:

transformed – Transformed points.

Return type:

(n, dim) float

trimesh.transformations.translation_from_matrix(matrix)

Return translation vector from translation matrix.

>>> v0 = np.random.random(3) - 0.5
>>> v1 = translation_from_matrix(translation_matrix(v0))
>>> np.allclose(v0, v1)
True
trimesh.transformations.translation_matrix(direction)

Return matrix to translate by direction vector.

>>> v = np.random.random(3) - 0.5
>>> np.allclose(v, translation_matrix(v)[:3, 3])
True
trimesh.transformations.unit_vector(data, axis=None, out=None)

Return ndarray normalized by length, i.e. Euclidean norm, along axis.

>>> v0 = np.random.random(3)
>>> v1 = unit_vector(v0)
>>> np.allclose(v1, v0 / np.linalg.norm(v0))
True
>>> v0 = np.random.rand(5, 4, 3)
>>> v1 = unit_vector(v0, axis=-1)
>>> v2 = v0 / np.expand_dims(np.sqrt(np.sum(v0*v0, axis=2)), 2)
>>> np.allclose(v1, v2)
True
>>> v1 = unit_vector(v0, axis=1)
>>> v2 = v0 / np.expand_dims(np.sqrt(np.sum(v0*v0, axis=1)), 1)
>>> np.allclose(v1, v2)
True
>>> v1 = np.zeros((5, 4, 3))
>>> unit_vector(v0, axis=1, out=v1)
>>> np.allclose(v1, v2)
True
>>> list(unit_vector([]))
[]
>>> [float(i) for i in unit_vector([1])]
[1.0]
trimesh.transformations.vector_norm(data, axis=None, out=None)

Return length, i.e. Euclidean norm, of ndarray along axis.

>>> v = np.random.random(3)
>>> n = vector_norm(v)
>>> np.allclose(n, np.linalg.norm(v))
True
>>> v = np.random.rand(6, 5, 3)
>>> n = vector_norm(v, axis=-1)
>>> np.allclose(n, np.sqrt(np.sum(v*v, axis=2)))
True
>>> n = vector_norm(v, axis=1)
>>> np.allclose(n, np.sqrt(np.sum(v*v, axis=1)))
True
>>> v = np.random.rand(5, 4, 3)
>>> n = np.zeros((5, 3))
>>> vector_norm(v, axis=1, out=n)
>>> np.allclose(n, np.sqrt(np.sum(v*v, axis=1)))
True
>>> float(vector_norm([]))
0.0
>>> float(vector_norm([1]))
1.0
trimesh.transformations.vector_product(v0, v1, axis=0)

Return vector perpendicular to vectors.

>>> v = vector_product([2, 0, 0], [0, 3, 0])
>>> np.allclose(v, [0, 0, 6])
True
>>> v0 = [[2, 0, 0, 2], [0, 2, 0, 2], [0, 0, 2, 2]]
>>> v1 = [[3], [0], [0]]
>>> v = vector_product(v0, v1)
>>> np.allclose(v, [[0, 0, 0, 0], [0, 0, 6, 6], [0, -6, 0, -6]])
True
>>> v0 = [[2, 0, 0], [2, 0, 0], [0, 2, 0], [2, 0, 0]]
>>> v1 = [[0, 3, 0], [0, 0, 3], [0, 0, 3], [3, 3, 3]]
>>> v = vector_product(v0, v1, axis=1)
>>> np.allclose(v, [[0, 0, 6], [0, -6, 0], [6, 0, 0], [0, -6, 6]])
True