__mul__#
- RigidTransform.__mul__()#
Compose this transform with the other.
If p and q are two transforms, then the composition of ‘q followed by p’ is equivalent to p * q. In terms of transformation matrices, the composition can be expressed as
p.as_matrix() @ q.as_matrix()
.In terms of translations and rotations, the composition when applied to a vector
v
is equivalent top.translation + p.rotation.apply(q.translation) + (p.rotation * q.rotation).apply(v)
.This function supports composition of multiple transforms at a time. The following cases are possible:
Either
p
orq
contains a single or length 1 transform. In this case the result contains the result of composing each transform in the other object with the one transform. If both are single transforms, the result is a single transform.Both
p
andq
containN
transforms. In this case each transformp[i]
is composed with the corresponding transformq[i]
and the result containsN
transforms.
- Parameters:
- other
RigidTransform
instance Object containing the transforms to be composed with this one.
- other
- Returns:
RigidTransform
instanceThe composed transform.
Examples
>>> from scipy.spatial.transform import RigidTransform as Tf >>> from scipy.spatial.transform import Rotation as R >>> import numpy as np
Compose two transforms:
>>> tf1 = Tf.from_translation([1, 0, 0]) >>> tf2 = Tf.from_translation([0, 1, 0]) >>> tf = tf1 * tf2 >>> tf.translation array([1., 1., 0.]) >>> tf.single True
When applied to a vector, the composition of two transforms is applied in right-to-left order.
>>> t1, r1 = [1, 2, 3], R.from_euler('z', 60, degrees=True) >>> t2, r2 = [0, 1, 0], R.from_euler('x', 30, degrees=True) >>> tf1 = Tf.from_components(t1, r1) >>> tf2 = Tf.from_components(t2, r2) >>> tf = tf1 * tf2 >>> tf.apply([1, 0, 0]) array([0.6339746, 3.3660254, 3. ]) >>> tf1.apply(tf2.apply([1, 0, 0])) array([0.6339746, 3.3660254, 3. ])
When at least one of the transforms is not single, the result is a stack of transforms.
>>> tf1 = Tf.from_translation([1, 0, 0]) >>> tf2 = Tf.from_translation([[0, 2, 0], [0, 0, 3]]) >>> tf = tf1 * tf2 >>> tf.translation array([[1., 2., 0.], [1., 0., 3.]]) >>> tf.single False >>> len(tf) 2