13.2.8. Homogeneous Transformation Matrices and Quaternions — MDAnalysis.lib.transformations
¶
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.
Authors: | Christoph Gohlke, Laboratory for Fluorescence Dynamics, University of California, Irvine |
---|---|
Version: | 2010.05.10 |
Licence: | BSD 3-clause |
13.2.8.1. Requirements¶
- Python 2.6 or 3.1
- Numpy 1.4
- transformations.c 2010.04.10 (optional implementation of some functions in C)
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 numpy.linalg.inv(M)
, concatenated using
numpy.dot(M0, M1)
, or used to transform homogeneous coordinates (v) using
numpy.dot(M, v)
for shape (4, *)
“point of arrays”, respectively
numpy.dot(v, M.T)
for shape (*, 4)
“array of points”.
Use the transpose of transformation matrices for OpenGL glMultMatrixd()
.
Calculations are carried out with numpy.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.
References
- Matrices and transformations. Ronald Goldman. In “Graphics Gems I”, pp 472-475. Morgan Kaufmann, 1990.
- More matrices and transformations: shear and pseudo-perspective. Ronald Goldman. In “Graphics Gems II”, pp 320-323. Morgan Kaufmann, 1991.
- Decomposing a matrix into simple transformations. Spencer Thomas. In “Graphics Gems II”, pp 320-323. Morgan Kaufmann, 1991.
- Recovering the data from the transformation matrix. Ronald Goldman. In “Graphics Gems II”, pp 324-331. Morgan Kaufmann, 1991.
- Euler angle conversion. Ken Shoemake. In “Graphics Gems IV”, pp 222-229. Morgan Kaufmann, 1994.
- Arcball rotation control. Ken Shoemake. In “Graphics Gems IV”, pp 175-192. Morgan Kaufmann, 1994.
- Representing attitude: Euler angles, unit quaternions, and rotation vectors. James Diebel. 2006.
- A discussion of the solution for the best rotation to relate two sets of vectors. W Kabsch. Acta Cryst. 1978. A34, 827-828.
- Closed-form solution of absolute orientation using unit quaternions. BKP Horn. J Opt Soc Am A. 1987. 4(4):629-642.
- Quaternions. Ken Shoemake. http://www.sfu.ca/~jwa3/cmpt461/files/quatut.pdf
- From quaternion to matrix and back. JMP van Waveren. 2005. http://www.intel.com/cd/ids/developer/asmo-na/eng/293748.htm
- Uniform random rotations. Ken Shoemake. In “Graphics Gems III”, pp 124-132. Morgan Kaufmann, 1992.
- Quaternion in molecular modeling. CFF Karney. J Mol Graph Mod, 25(5):595-604
- New method for extracting the quaternion from a rotation matrix. Itzhack Y Bar-Itzhack, J Guid Contr Dynam. 2000. 23(6): 1085-1087.
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')
>>> numpy.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(numpy.random.rand(3))
>>> M = concatenate_matrices(T, R, Z, S)
>>> scale, shear, angles, trans, persp = decompose_matrix(M)
>>> numpy.allclose(scale, 1.23)
True
>>> numpy.allclose(trans, (1, 2, 3))
True
>>> numpy.allclose(shear, (0, math.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
13.2.8.2. Functions¶
Changed in version 0.11.0: Transformations library moved from MDAnalysis.core.transformations to MDAnalysis.lib.transformations
-
class
MDAnalysis.lib.transformations.
Arcball
(initial=None)[source]¶ 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.setconstrain(True) >>> ball.down([400, 200]) >>> ball.drag([200, 400]) >>> R = ball.matrix() >>> np.allclose(np.sum(R), 0.2055924) True >>> ball.next()
Initialize virtual trackball control.
initial : quaternion or rotation matrix
-
MDAnalysis.lib.transformations.
arcball_nearest_axis
(point, axes)[source]¶ Return axis, which arc is nearest to point.
-
MDAnalysis.lib.transformations.
compose_matrix
(scale=None, shear=None, angles=None, translate=None, perspective=None)[source]¶ 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*math.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
-
MDAnalysis.lib.transformations.
concatenate_matrices
(*matrices)[source]¶ 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
-
MDAnalysis.lib.transformations.
decompose_matrix
(matrix)[source]¶ Return sequence of transformations from transformation matrix.
- matrix : array_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) >>> scale[0] 0.123 >>> R0 = euler_matrix(1, 2, 3) >>> scale, shear, angles, trans, persp = decompose_matrix(R0) >>> R1 = euler_matrix(*angles) >>> np.allclose(R0, R1) True
-
MDAnalysis.lib.transformations.
euler_from_quaternion
(quaternion, axes='sxyz')[source]¶ 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
-
MDAnalysis.lib.transformations.
projection_from_matrix
(matrix, pseudo=False)[source]¶ 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
-
MDAnalysis.lib.transformations.
quaternion_imag
(quaternion)[source]¶ Return imaginary part of quaternion.
>>> quaternion_imag([3.0, 0.0, 1.0, 2.0]) [0.0, 1.0, 2.0]
-
MDAnalysis.lib.transformations.
quaternion_real
(quaternion)[source]¶ Return real part of quaternion.
>>> quaternion_real([3.0, 0.0, 1.0, 2.0]) 3.0
-
MDAnalysis.lib.transformations.
reflection_from_matrix
(matrix)[source]¶ 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
-
MDAnalysis.lib.transformations.
rotation_from_matrix
(matrix)[source]¶ Return rotation angle and axis from rotation matrix.
>>> angle = (random.random() - 0.5) * (2*math.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
-
MDAnalysis.lib.transformations.
rotaxis
(a, b)[source]¶ Return the rotation axis to rotate vector a into b.
Parameters: b (a,) – two vectors Returns: c – vector to rotate a into b Return type: np.ndarray Note
If a == b this will always return [1, 0, 0]
-
MDAnalysis.lib.transformations.
scale_from_matrix
(matrix)[source]¶ 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
-
MDAnalysis.lib.transformations.
shear_from_matrix
(matrix)[source]¶ Return shear angle, direction and plane from shear matrix.
>>> angle = (random.random() - 0.5) * 4*math.pi >>> direct = np.random.random(3) - 0.5 >>> point = np.random.random(3) - 0.5 >>> normal = np.cross(direct, np.random.random(3)) >>> 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