Type Alias nalgebra::geometry::UnitQuaternion

source ·
pub type UnitQuaternion<T> = Unit<Quaternion<T>>;
Expand description

A unit quaternions. May be used to represent a rotation.

Aliased Type§

struct UnitQuaternion<T> { /* private fields */ }

Implementations§

source§

impl<T: SimdRealField> UnitQuaternion<T>

source

pub fn angle(&self) -> T

The rotation angle in [0; pi] of this unit quaternion.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let rot = UnitQuaternion::from_axis_angle(&axis, 1.78);
assert_eq!(rot.angle(), 1.78);
source

pub fn quaternion(&self) -> &Quaternion<T>

The underlying quaternion.

Same as self.as_ref().

§Example
let axis = UnitQuaternion::identity();
assert_eq!(*axis.quaternion(), Quaternion::new(1.0, 0.0, 0.0, 0.0));
source

pub fn conjugate(&self) -> Self

Compute the conjugate of this unit quaternion.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let rot = UnitQuaternion::from_axis_angle(&axis, 1.78);
let conj = rot.conjugate();
assert_eq!(conj, UnitQuaternion::from_axis_angle(&-axis, 1.78));
source

pub fn inverse(&self) -> Self

Inverts this quaternion if it is not zero.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let rot = UnitQuaternion::from_axis_angle(&axis, 1.78);
let inv = rot.inverse();
assert_eq!(rot * inv, UnitQuaternion::identity());
assert_eq!(inv * rot, UnitQuaternion::identity());
source

pub fn angle_to(&self, other: &Self) -> T

The rotation angle needed to make self and other coincide.

§Example
let rot1 = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), 1.0);
let rot2 = UnitQuaternion::from_axis_angle(&Vector3::x_axis(), 0.1);
assert_relative_eq!(rot1.angle_to(&rot2), 1.0045657, epsilon = 1.0e-6);
source

pub fn rotation_to(&self, other: &Self) -> Self

The unit quaternion needed to make self and other coincide.

The result is such that: self.rotation_to(other) * self == other.

§Example
let rot1 = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), 1.0);
let rot2 = UnitQuaternion::from_axis_angle(&Vector3::x_axis(), 0.1);
let rot_to = rot1.rotation_to(&rot2);
assert_relative_eq!(rot_to * rot1, rot2, epsilon = 1.0e-6);
source

pub fn lerp(&self, other: &Self, t: T) -> Quaternion<T>

Linear interpolation between two unit quaternions.

The result is not normalized.

§Example
let q1 = UnitQuaternion::new_normalize(Quaternion::new(1.0, 0.0, 0.0, 0.0));
let q2 = UnitQuaternion::new_normalize(Quaternion::new(0.0, 1.0, 0.0, 0.0));
assert_eq!(q1.lerp(&q2, 0.1), Quaternion::new(0.9, 0.1, 0.0, 0.0));
source

pub fn nlerp(&self, other: &Self, t: T) -> Self

Normalized linear interpolation between two unit quaternions.

This is the same as self.lerp except that the result is normalized.

§Example
let q1 = UnitQuaternion::new_normalize(Quaternion::new(1.0, 0.0, 0.0, 0.0));
let q2 = UnitQuaternion::new_normalize(Quaternion::new(0.0, 1.0, 0.0, 0.0));
assert_eq!(q1.nlerp(&q2, 0.1), UnitQuaternion::new_normalize(Quaternion::new(0.9, 0.1, 0.0, 0.0)));
source

pub fn slerp(&self, other: &Self, t: T) -> Self
where T: RealField,

Spherical linear interpolation between two unit quaternions.

Panics if the angle between both quaternion is 180 degrees (in which case the interpolation is not well-defined). Use .try_slerp instead to avoid the panic.

§Examples:

let q1 = UnitQuaternion::from_euler_angles(std::f32::consts::FRAC_PI_4, 0.0, 0.0);
let q2 = UnitQuaternion::from_euler_angles(-std::f32::consts::PI, 0.0, 0.0);

let q = q1.slerp(&q2, 1.0 / 3.0);

assert_eq!(q.euler_angles(), (std::f32::consts::FRAC_PI_2, 0.0, 0.0));
source

pub fn try_slerp(&self, other: &Self, t: T, epsilon: T) -> Option<Self>
where T: RealField,

Computes the spherical linear interpolation between two unit quaternions or returns None if both quaternions are approximately 180 degrees apart (in which case the interpolation is not well-defined).

§Arguments
  • self: the first quaternion to interpolate from.
  • other: the second quaternion to interpolate toward.
  • t: the interpolation parameter. Should be between 0 and 1.
  • epsilon: the value below which the sinus of the angle separating both quaternion must be to return None.
source

pub fn conjugate_mut(&mut self)

Compute the conjugate of this unit quaternion in-place.

source

pub fn inverse_mut(&mut self)

Inverts this quaternion if it is not zero.

§Example
let axisangle = Vector3::new(0.1, 0.2, 0.3);
let mut rot = UnitQuaternion::new(axisangle);
rot.inverse_mut();
assert_relative_eq!(rot * UnitQuaternion::new(axisangle), UnitQuaternion::identity());
assert_relative_eq!(UnitQuaternion::new(axisangle) * rot, UnitQuaternion::identity());
source

pub fn axis(&self) -> Option<Unit<Vector3<T>>>
where T: RealField,

The rotation axis of this unit quaternion or None if the rotation is zero.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let angle = 1.2;
let rot = UnitQuaternion::from_axis_angle(&axis, angle);
assert_eq!(rot.axis(), Some(axis));

// Case with a zero angle.
let rot = UnitQuaternion::from_axis_angle(&axis, 0.0);
assert!(rot.axis().is_none());
source

pub fn scaled_axis(&self) -> Vector3<T>
where T: RealField,

The rotation axis of this unit quaternion multiplied by the rotation angle.

§Example
let axisangle = Vector3::new(0.1, 0.2, 0.3);
let rot = UnitQuaternion::new(axisangle);
assert_relative_eq!(rot.scaled_axis(), axisangle, epsilon = 1.0e-6);
source

pub fn axis_angle(&self) -> Option<(Unit<Vector3<T>>, T)>
where T: RealField,

The rotation axis and angle in ]0, pi] of this unit quaternion.

Returns None if the angle is zero.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let angle = 1.2;
let rot = UnitQuaternion::from_axis_angle(&axis, angle);
assert_eq!(rot.axis_angle(), Some((axis, angle)));

// Case with a zero angle.
let rot = UnitQuaternion::from_axis_angle(&axis, 0.0);
assert!(rot.axis_angle().is_none());
source

pub fn exp(&self) -> Quaternion<T>

Compute the exponential of a quaternion.

Note that this function yields a Quaternion<T> because it loses the unit property.

source

pub fn ln(&self) -> Quaternion<T>
where T: RealField,

Compute the natural logarithm of a quaternion.

Note that this function yields a Quaternion<T> because it loses the unit property. The vector part of the return value corresponds to the axis-angle representation (divided by 2.0) of this unit quaternion.

§Example
let axisangle = Vector3::new(0.1, 0.2, 0.3);
let q = UnitQuaternion::new(axisangle);
assert_relative_eq!(q.ln().vector().into_owned(), axisangle, epsilon = 1.0e-6);
source

pub fn powf(&self, n: T) -> Self
where T: RealField,

Raise the quaternion to a given floating power.

This returns the unit quaternion that identifies a rotation with axis self.axis() and angle self.angle() × n.

§Example
let axis = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let angle = 1.2;
let rot = UnitQuaternion::from_axis_angle(&axis, angle);
let pow = rot.powf(2.0);
assert_relative_eq!(pow.axis().unwrap(), axis, epsilon = 1.0e-6);
assert_eq!(pow.angle(), 2.4);
source

pub fn to_rotation_matrix(self) -> Rotation<T, 3>

Builds a rotation matrix from this unit quaternion.

§Example
let q = UnitQuaternion::from_axis_angle(&Vector3::z_axis(), f32::consts::FRAC_PI_6);
let rot = q.to_rotation_matrix();
let expected = Matrix3::new(0.8660254, -0.5,      0.0,
                            0.5,       0.8660254, 0.0,
                            0.0,       0.0,       1.0);

assert_relative_eq!(*rot.matrix(), expected, epsilon = 1.0e-6);
source

pub fn to_euler_angles(self) -> (T, T, T)
where T: RealField,

👎Deprecated: This is renamed to use .euler_angles().

Converts this unit quaternion into its equivalent Euler angles.

The angles are produced in the form (roll, pitch, yaw).

source

pub fn euler_angles(&self) -> (T, T, T)
where T: RealField,

Retrieves the euler angles corresponding to this unit quaternion.

The angles are produced in the form (roll, pitch, yaw).

§Example
let rot = UnitQuaternion::from_euler_angles(0.1, 0.2, 0.3);
let euler = rot.euler_angles();
assert_relative_eq!(euler.0, 0.1, epsilon = 1.0e-6);
assert_relative_eq!(euler.1, 0.2, epsilon = 1.0e-6);
assert_relative_eq!(euler.2, 0.3, epsilon = 1.0e-6);
source

pub fn to_homogeneous(self) -> Matrix4<T>

Converts this unit quaternion into its equivalent homogeneous transformation matrix.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::z_axis(), f32::consts::FRAC_PI_6);
let expected = Matrix4::new(0.8660254, -0.5,      0.0, 0.0,
                            0.5,       0.8660254, 0.0, 0.0,
                            0.0,       0.0,       1.0, 0.0,
                            0.0,       0.0,       0.0, 1.0);

assert_relative_eq!(rot.to_homogeneous(), expected, epsilon = 1.0e-6);
source

pub fn transform_point(&self, pt: &Point3<T>) -> Point3<T>

Rotate a point by this unit quaternion.

This is the same as the multiplication self * pt.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), f32::consts::FRAC_PI_2);
let transformed_point = rot.transform_point(&Point3::new(1.0, 2.0, 3.0));

assert_relative_eq!(transformed_point, Point3::new(3.0, 2.0, -1.0), epsilon = 1.0e-6);
source

pub fn transform_vector(&self, v: &Vector3<T>) -> Vector3<T>

Rotate a vector by this unit quaternion.

This is the same as the multiplication self * v.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), f32::consts::FRAC_PI_2);
let transformed_vector = rot.transform_vector(&Vector3::new(1.0, 2.0, 3.0));

assert_relative_eq!(transformed_vector, Vector3::new(3.0, 2.0, -1.0), epsilon = 1.0e-6);
source

pub fn inverse_transform_point(&self, pt: &Point3<T>) -> Point3<T>

Rotate a point by the inverse of this unit quaternion. This may be cheaper than inverting the unit quaternion and transforming the point.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), f32::consts::FRAC_PI_2);
let transformed_point = rot.inverse_transform_point(&Point3::new(1.0, 2.0, 3.0));

assert_relative_eq!(transformed_point, Point3::new(-3.0, 2.0, 1.0), epsilon = 1.0e-6);
source

pub fn inverse_transform_vector(&self, v: &Vector3<T>) -> Vector3<T>

Rotate a vector by the inverse of this unit quaternion. This may be cheaper than inverting the unit quaternion and transforming the vector.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::y_axis(), f32::consts::FRAC_PI_2);
let transformed_vector = rot.inverse_transform_vector(&Vector3::new(1.0, 2.0, 3.0));

assert_relative_eq!(transformed_vector, Vector3::new(-3.0, 2.0, 1.0), epsilon = 1.0e-6);
source

pub fn inverse_transform_unit_vector( &self, v: &Unit<Vector3<T>>, ) -> Unit<Vector3<T>>

Rotate a vector by the inverse of this unit quaternion. This may be cheaper than inverting the unit quaternion and transforming the vector.

§Example
let rot = UnitQuaternion::from_axis_angle(&Vector3::z_axis(), f32::consts::FRAC_PI_2);
let transformed_vector = rot.inverse_transform_unit_vector(&Vector3::x_axis());

assert_relative_eq!(transformed_vector, -Vector3::y_axis(), epsilon = 1.0e-6);
source

pub fn append_axisangle_linearized(&self, axisangle: &Vector3<T>) -> Self

Appends to self a rotation given in the axis-angle form, using a linearized formulation.

This is faster, but approximate, way to compute UnitQuaternion::new(axisangle) * self.

source§

impl<T: SimdRealField> UnitQuaternion<T>

source

pub fn identity() -> Self

The rotation identity.

§Example
let q = UnitQuaternion::identity();
let q2 = UnitQuaternion::new(Vector3::new(1.0, 2.0, 3.0));
let v = Vector3::new_random();
let p = Point3::from(v);

assert_eq!(q * q2, q2);
assert_eq!(q2 * q, q2);
assert_eq!(q * v, v);
assert_eq!(q * p, p);
source

pub fn cast<To>(self) -> UnitQuaternion<To>
where To: SupersetOf<T> + Scalar,

Cast the components of self to another type.

§Example
let q = UnitQuaternion::from_euler_angles(1.0f64, 2.0, 3.0);
let q2 = q.cast::<f32>();
assert_relative_eq!(q2, UnitQuaternion::from_euler_angles(1.0f32, 2.0, 3.0), epsilon = 1.0e-6);
source

pub fn from_axis_angle<SB>(axis: &Unit<Vector<T, U3, SB>>, angle: T) -> Self
where SB: Storage<T, U3>,

Creates a new quaternion from a unit vector (the rotation axis) and an angle (the rotation angle).

§Example
let axis = Vector3::y_axis();
let angle = f32::consts::FRAC_PI_2;
// Point and vector being transformed in the tests.
let pt = Point3::new(4.0, 5.0, 6.0);
let vec = Vector3::new(4.0, 5.0, 6.0);
let q = UnitQuaternion::from_axis_angle(&axis, angle);

assert_eq!(q.axis().unwrap(), axis);
assert_eq!(q.angle(), angle);
assert_relative_eq!(q * pt, Point3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);
assert_relative_eq!(q * vec, Vector3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);

// A zero vector yields an identity.
assert_eq!(UnitQuaternion::from_scaled_axis(Vector3::<f32>::zeros()), UnitQuaternion::identity());
source

pub fn from_quaternion(q: Quaternion<T>) -> Self

Creates a new unit quaternion from a quaternion.

The input quaternion will be normalized.

source

pub fn from_euler_angles(roll: T, pitch: T, yaw: T) -> Self

Creates a new unit quaternion from Euler angles.

The primitive rotations are applied in order: 1 roll − 2 pitch − 3 yaw.

§Example
let rot = UnitQuaternion::from_euler_angles(0.1, 0.2, 0.3);
let euler = rot.euler_angles();
assert_relative_eq!(euler.0, 0.1, epsilon = 1.0e-6);
assert_relative_eq!(euler.1, 0.2, epsilon = 1.0e-6);
assert_relative_eq!(euler.2, 0.3, epsilon = 1.0e-6);
source

pub fn from_basis_unchecked(basis: &[Vector3<T>; 3]) -> Self

Builds an unit quaternion from a basis assumed to be orthonormal.

In order to get a valid unit-quaternion, the input must be an orthonormal basis, i.e., all vectors are normalized, and the are all orthogonal to each other. These invariants are not checked by this method.

source

pub fn from_rotation_matrix(rotmat: &Rotation3<T>) -> Self

Builds an unit quaternion from a rotation matrix.

§Example
let axis = Vector3::y_axis();
let angle = 0.1;
let rot = Rotation3::from_axis_angle(&axis, angle);
let q = UnitQuaternion::from_rotation_matrix(&rot);
assert_relative_eq!(q.to_rotation_matrix(), rot, epsilon = 1.0e-6);
assert_relative_eq!(q.axis().unwrap(), rot.axis().unwrap(), epsilon = 1.0e-6);
assert_relative_eq!(q.angle(), rot.angle(), epsilon = 1.0e-6);
source

pub fn from_matrix(m: &Matrix3<T>) -> Self
where T: RealField,

Builds an unit quaternion by extracting the rotation part of the given transformation m.

This is an iterative method. See .from_matrix_eps to provide mover convergence parameters and starting solution. This implements “A Robust Method to Extract the Rotational Part of Deformations” by Müller et al.

source

pub fn from_matrix_eps( m: &Matrix3<T>, eps: T, max_iter: usize, guess: Self, ) -> Self
where T: RealField,

Builds an unit quaternion by extracting the rotation part of the given transformation m.

This implements “A Robust Method to Extract the Rotational Part of Deformations” by Müller et al.

§Parameters
  • m: the matrix from which the rotational part is to be extracted.
  • eps: the angular errors tolerated between the current rotation and the optimal one.
  • max_iter: the maximum number of iterations. Loops indefinitely until convergence if set to 0.
  • guess: an estimate of the solution. Convergence will be significantly faster if an initial solution close to the actual solution is provided. Can be set to UnitQuaternion::identity() if no other guesses come to mind.
source

pub fn rotation_between<SB, SC>( a: &Vector<T, U3, SB>, b: &Vector<T, U3, SC>, ) -> Option<Self>
where T: RealField, SB: Storage<T, U3>, SC: Storage<T, U3>,

The unit quaternion needed to make a and b be collinear and point toward the same direction. Returns None if both a and b are collinear and point to opposite directions, as then the rotation desired is not unique.

§Example
let a = Vector3::new(1.0, 2.0, 3.0);
let b = Vector3::new(3.0, 1.0, 2.0);
let q = UnitQuaternion::rotation_between(&a, &b).unwrap();
assert_relative_eq!(q * a, b);
assert_relative_eq!(q.inverse() * b, a);
source

pub fn scaled_rotation_between<SB, SC>( a: &Vector<T, U3, SB>, b: &Vector<T, U3, SC>, s: T, ) -> Option<Self>
where T: RealField, SB: Storage<T, U3>, SC: Storage<T, U3>,

The smallest rotation needed to make a and b collinear and point toward the same direction, raised to the power s.

§Example
let a = Vector3::new(1.0, 2.0, 3.0);
let b = Vector3::new(3.0, 1.0, 2.0);
let q2 = UnitQuaternion::scaled_rotation_between(&a, &b, 0.2).unwrap();
let q5 = UnitQuaternion::scaled_rotation_between(&a, &b, 0.5).unwrap();
assert_relative_eq!(q2 * q2 * q2 * q2 * q2 * a, b, epsilon = 1.0e-6);
assert_relative_eq!(q5 * q5 * a, b, epsilon = 1.0e-6);
source

pub fn rotation_between_axis<SB, SC>( a: &Unit<Vector<T, U3, SB>>, b: &Unit<Vector<T, U3, SC>>, ) -> Option<Self>
where T: RealField, SB: Storage<T, U3>, SC: Storage<T, U3>,

The unit quaternion needed to make a and b be collinear and point toward the same direction.

§Example
let a = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let b = Unit::new_normalize(Vector3::new(3.0, 1.0, 2.0));
let q = UnitQuaternion::rotation_between(&a, &b).unwrap();
assert_relative_eq!(q * a, b);
assert_relative_eq!(q.inverse() * b, a);
source

pub fn scaled_rotation_between_axis<SB, SC>( na: &Unit<Vector<T, U3, SB>>, nb: &Unit<Vector<T, U3, SC>>, s: T, ) -> Option<Self>
where T: RealField, SB: Storage<T, U3>, SC: Storage<T, U3>,

The smallest rotation needed to make a and b collinear and point toward the same direction, raised to the power s.

§Example
let a = Unit::new_normalize(Vector3::new(1.0, 2.0, 3.0));
let b = Unit::new_normalize(Vector3::new(3.0, 1.0, 2.0));
let q2 = UnitQuaternion::scaled_rotation_between(&a, &b, 0.2).unwrap();
let q5 = UnitQuaternion::scaled_rotation_between(&a, &b, 0.5).unwrap();
assert_relative_eq!(q2 * q2 * q2 * q2 * q2 * a, b, epsilon = 1.0e-6);
assert_relative_eq!(q5 * q5 * a, b, epsilon = 1.0e-6);
source

pub fn face_towards<SB, SC>( dir: &Vector<T, U3, SB>, up: &Vector<T, U3, SC>, ) -> Self
where SB: Storage<T, U3>, SC: Storage<T, U3>,

Creates an unit quaternion that corresponds to the local frame of an observer standing at the origin and looking toward dir.

It maps the z axis to the direction dir.

§Arguments
  • dir - The look direction. It does not need to be normalized.
  • up - The vertical direction. It does not need to be normalized. The only requirement of this parameter is to not be collinear to dir. Non-collinearity is not checked.
§Example
let dir = Vector3::new(1.0, 2.0, 3.0);
let up = Vector3::y();

let q = UnitQuaternion::face_towards(&dir, &up);
assert_relative_eq!(q * Vector3::z(), dir.normalize());
source

pub fn new_observer_frames<SB, SC>( dir: &Vector<T, U3, SB>, up: &Vector<T, U3, SC>, ) -> Self
where SB: Storage<T, U3>, SC: Storage<T, U3>,

👎Deprecated: renamed to face_towards

Deprecated: Use UnitQuaternion::face_towards instead.

source

pub fn look_at_rh<SB, SC>( dir: &Vector<T, U3, SB>, up: &Vector<T, U3, SC>, ) -> Self
where SB: Storage<T, U3>, SC: Storage<T, U3>,

Builds a right-handed look-at view matrix without translation.

It maps the view direction dir to the negative z axis. This conforms to the common notion of right handed look-at matrix from the computer graphics community.

§Arguments
  • dir − The view direction. It does not need to be normalized.
  • up - A vector approximately aligned with required the vertical axis. It does not need to be normalized. The only requirement of this parameter is to not be collinear to dir.
§Example
let dir = Vector3::new(1.0, 2.0, 3.0);
let up = Vector3::y();

let q = UnitQuaternion::look_at_rh(&dir, &up);
assert_relative_eq!(q * dir.normalize(), -Vector3::z());
source

pub fn look_at_lh<SB, SC>( dir: &Vector<T, U3, SB>, up: &Vector<T, U3, SC>, ) -> Self
where SB: Storage<T, U3>, SC: Storage<T, U3>,

Builds a left-handed look-at view matrix without translation.

It maps the view direction dir to the positive z axis. This conforms to the common notion of left handed look-at matrix from the computer graphics community.

§Arguments
  • dir − The view direction. It does not need to be normalized.
  • up - A vector approximately aligned with required the vertical axis. The only requirement of this parameter is to not be collinear to dir.
§Example
let dir = Vector3::new(1.0, 2.0, 3.0);
let up = Vector3::y();

let q = UnitQuaternion::look_at_lh(&dir, &up);
assert_relative_eq!(q * dir.normalize(), Vector3::z());
source

pub fn new<SB>(axisangle: Vector<T, U3, SB>) -> Self
where SB: Storage<T, U3>,

Creates a new unit quaternion rotation from a rotation axis scaled by the rotation angle.

If axisangle has a magnitude smaller than T::default_epsilon(), this returns the identity rotation.

§Example
let axisangle = Vector3::y() * f32::consts::FRAC_PI_2;
// Point and vector being transformed in the tests.
let pt = Point3::new(4.0, 5.0, 6.0);
let vec = Vector3::new(4.0, 5.0, 6.0);
let q = UnitQuaternion::new(axisangle);

assert_relative_eq!(q * pt, Point3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);
assert_relative_eq!(q * vec, Vector3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);

// A zero vector yields an identity.
assert_eq!(UnitQuaternion::new(Vector3::<f32>::zeros()), UnitQuaternion::identity());
source

pub fn new_eps<SB>(axisangle: Vector<T, U3, SB>, eps: T) -> Self
where SB: Storage<T, U3>,

Creates a new unit quaternion rotation from a rotation axis scaled by the rotation angle.

If axisangle has a magnitude smaller than eps, this returns the identity rotation.

§Example
let axisangle = Vector3::y() * f32::consts::FRAC_PI_2;
// Point and vector being transformed in the tests.
let pt = Point3::new(4.0, 5.0, 6.0);
let vec = Vector3::new(4.0, 5.0, 6.0);
let q = UnitQuaternion::new_eps(axisangle, 1.0e-6);

assert_relative_eq!(q * pt, Point3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);
assert_relative_eq!(q * vec, Vector3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);

// An almost zero vector yields an identity.
assert_eq!(UnitQuaternion::new_eps(Vector3::new(1.0e-8, 1.0e-9, 1.0e-7), 1.0e-6), UnitQuaternion::identity());
source

pub fn from_scaled_axis<SB>(axisangle: Vector<T, U3, SB>) -> Self
where SB: Storage<T, U3>,

Creates a new unit quaternion rotation from a rotation axis scaled by the rotation angle.

If axisangle has a magnitude smaller than T::default_epsilon(), this returns the identity rotation. Same as Self::new(axisangle).

§Example
let axisangle = Vector3::y() * f32::consts::FRAC_PI_2;
// Point and vector being transformed in the tests.
let pt = Point3::new(4.0, 5.0, 6.0);
let vec = Vector3::new(4.0, 5.0, 6.0);
let q = UnitQuaternion::from_scaled_axis(axisangle);

assert_relative_eq!(q * pt, Point3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);
assert_relative_eq!(q * vec, Vector3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);

// A zero vector yields an identity.
assert_eq!(UnitQuaternion::from_scaled_axis(Vector3::<f32>::zeros()), UnitQuaternion::identity());
source

pub fn from_scaled_axis_eps<SB>(axisangle: Vector<T, U3, SB>, eps: T) -> Self
where SB: Storage<T, U3>,

Creates a new unit quaternion rotation from a rotation axis scaled by the rotation angle.

If axisangle has a magnitude smaller than eps, this returns the identity rotation. Same as Self::new_eps(axisangle, eps).

§Example
let axisangle = Vector3::y() * f32::consts::FRAC_PI_2;
// Point and vector being transformed in the tests.
let pt = Point3::new(4.0, 5.0, 6.0);
let vec = Vector3::new(4.0, 5.0, 6.0);
let q = UnitQuaternion::from_scaled_axis_eps(axisangle, 1.0e-6);

assert_relative_eq!(q * pt, Point3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);
assert_relative_eq!(q * vec, Vector3::new(6.0, 5.0, -4.0), epsilon = 1.0e-6);

// An almost zero vector yields an identity.
assert_eq!(UnitQuaternion::from_scaled_axis_eps(Vector3::new(1.0e-8, 1.0e-9, 1.0e-7), 1.0e-6), UnitQuaternion::identity());
source

pub fn mean_of(unit_quaternions: impl IntoIterator<Item = Self>) -> Self
where T: RealField,

Create the mean unit quaternion from a data structure implementing IntoIterator returning unit quaternions.

The method will panic if the iterator does not return any quaternions.

Algorithm from: Oshman, Yaakov, and Avishy Carmi. “Attitude estimation from vector observations using a genetic-algorithm-embedded quaternion particle filter.” Journal of Guidance, Control, and Dynamics 29.4 (2006): 879-891.

§Example
let q1 = UnitQuaternion::from_euler_angles(0.0, 0.0, 0.0);
let q2 = UnitQuaternion::from_euler_angles(-0.1, 0.0, 0.0);
let q3 = UnitQuaternion::from_euler_angles(0.1, 0.0, 0.0);

let quat_vec = vec![q1, q2, q3];
let q_mean = UnitQuaternion::mean_of(quat_vec);

let euler_angles_mean = q_mean.euler_angles();
assert_relative_eq!(euler_angles_mean.0, 0.0, epsilon = 1.0e-7)

Trait Implementations§

source§

impl<T: RealField + AbsDiffEq<Epsilon = T>> AbsDiffEq for UnitQuaternion<T>

source§

type Epsilon = T

Used for specifying relative comparisons.
source§

fn default_epsilon() -> Self::Epsilon

The default tolerance to use when testing values that are close together. Read more
source§

fn abs_diff_eq(&self, other: &Self, epsilon: Self::Epsilon) -> bool

A test for equality that uses the absolute difference to compute the approximate equality of two numbers.
source§

fn abs_diff_ne(&self, other: &Rhs, epsilon: Self::Epsilon) -> bool

The inverse of AbsDiffEq::abs_diff_eq.
source§

impl<T: SimdRealField> AbstractRotation<T, 3> for UnitQuaternion<T>

source§

fn identity() -> Self

The rotation identity.
source§

fn inverse(&self) -> Self

The rotation inverse.
source§

fn inverse_mut(&mut self)

Change self to its inverse.
source§

fn transform_vector(&self, v: &SVector<T, 3>) -> SVector<T, 3>

Apply the rotation to the given vector.
source§

fn transform_point(&self, p: &Point<T, 3>) -> Point<T, 3>

Apply the rotation to the given point.
source§

fn inverse_transform_vector(&self, v: &SVector<T, 3>) -> SVector<T, 3>

Apply the inverse rotation to the given vector.
source§

fn inverse_transform_point(&self, p: &Point<T, 3>) -> Point<T, 3>

Apply the inverse rotation to the given point.
source§

fn inverse_transform_unit_vector( &self, v: &Unit<SVector<T, D>>, ) -> Unit<SVector<T, D>>

Apply the inverse rotation to the given unit vector.
source§

impl<T: RealField> Default for UnitQuaternion<T>

source§

fn default() -> Self

Returns the “default value” for a type. Read more
source§

impl<T: RealField + Display> Display for UnitQuaternion<T>

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
source§

impl<'a, 'b, T: SimdRealField> Div<&'b Isometry<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: &'b Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> Div<&'b Isometry<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: &'b Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Div<&'b Rotation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b Rotation<T, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> Div<&'b Rotation<T, 3>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b Rotation<T, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Div<&'b Similarity<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: &'b Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> Div<&'b Similarity<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: &'b Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, 'b, T, C> Div<&'b Transform<T, C, 3>> for &'a UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b Transform<T, C, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T, C> Div<&'b Transform<T, C, 3>> for UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b Transform<T, C, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Div<&'b Unit<DualQuaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b UnitDualQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> Div<&'b Unit<DualQuaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b UnitDualQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Div<&'b Unit<Quaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b UnitQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> Div<&'b Unit<Quaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: &'b UnitQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T: SimdRealField> Div<Isometry<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<T: SimdRealField> Div<Isometry<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T: SimdRealField> Div<Rotation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: Rotation<T, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<T: SimdRealField> Div<Rotation<T, 3>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: Rotation<T, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T: SimdRealField> Div<Similarity<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<T: SimdRealField> Div<Similarity<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the / operator.
source§

fn div(self, right: Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T, C> Div<Transform<T, C, 3>> for &'a UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the / operator.
source§

fn div(self, rhs: Transform<T, C, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<T, C> Div<Transform<T, C, 3>> for UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the / operator.
source§

fn div(self, rhs: Transform<T, C, 3>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T: SimdRealField> Div<Unit<DualQuaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: UnitDualQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<T: SimdRealField> Div<Unit<DualQuaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: UnitDualQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'a, T: SimdRealField> Div<Unit<Quaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: UnitQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<T: SimdRealField> Div for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the / operator.
source§

fn div(self, rhs: UnitQuaternion<T>) -> Self::Output

Performs the / operation. Read more
source§

impl<'b, T: SimdRealField> DivAssign<&'b Rotation<T, 3>> for UnitQuaternion<T>

source§

fn div_assign(&mut self, rhs: &'b Rotation<T, 3>)

Performs the /= operation. Read more
source§

impl<'b, T: SimdRealField> DivAssign<&'b Unit<Quaternion<T>>> for UnitQuaternion<T>

source§

fn div_assign(&mut self, rhs: &'b UnitQuaternion<T>)

Performs the /= operation. Read more
source§

impl<T: SimdRealField> DivAssign<Rotation<T, 3>> for UnitQuaternion<T>

source§

fn div_assign(&mut self, rhs: Rotation<T, 3>)

Performs the /= operation. Read more
source§

impl<T: SimdRealField> DivAssign for UnitQuaternion<T>

source§

fn div_assign(&mut self, rhs: UnitQuaternion<T>)

Performs the /= operation. Read more
source§

impl<T> From<[Unit<Quaternion<<T as SimdValue>::Element>>; 16]> for UnitQuaternion<T>

source§

fn from(arr: [UnitQuaternion<T::Element>; 16]) -> Self

Converts to this type from the input type.
source§

impl<T> From<[Unit<Quaternion<<T as SimdValue>::Element>>; 2]> for UnitQuaternion<T>
where T: From<[<T as SimdValue>::Element; 2]> + Scalar + Copy + PrimitiveSimdValue, T::Element: Scalar + Copy,

source§

fn from(arr: [UnitQuaternion<T::Element>; 2]) -> Self

Converts to this type from the input type.
source§

impl<T> From<[Unit<Quaternion<<T as SimdValue>::Element>>; 4]> for UnitQuaternion<T>
where T: From<[<T as SimdValue>::Element; 4]> + Scalar + Copy + PrimitiveSimdValue, T::Element: Scalar + Copy,

source§

fn from(arr: [UnitQuaternion<T::Element>; 4]) -> Self

Converts to this type from the input type.
source§

impl<T> From<[Unit<Quaternion<<T as SimdValue>::Element>>; 8]> for UnitQuaternion<T>
where T: From<[<T as SimdValue>::Element; 8]> + Scalar + Copy + PrimitiveSimdValue, T::Element: Scalar + Copy,

source§

fn from(arr: [UnitQuaternion<T::Element>; 8]) -> Self

Converts to this type from the input type.
source§

impl<T: SimdRealField> From<Rotation<T, 3>> for UnitQuaternion<T>

source§

fn from(q: Rotation3<T>) -> Self

Converts to this type from the input type.
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Isometry<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Isometry<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<&'b Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>> for UnitQuaternion<T>

source§

type Output = Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Vector<T, U3, SB>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<&'b Matrix<T, Const<3>, Const<1>, SB>> for &'a UnitQuaternion<T>

source§

type Output = Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Vector<T, Const<3>, SB>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b OPoint<T, Const<3>>> for &'a UnitQuaternion<T>

source§

type Output = OPoint<T, Const<3>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Point3<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b OPoint<T, Const<3>>> for UnitQuaternion<T>

source§

type Output = OPoint<T, Const<3>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Point3<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Rotation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Rotation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Rotation<T, 3>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Rotation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Similarity<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Similarity<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T, C> Mul<&'b Transform<T, C, 3>> for &'a UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Transform<T, C, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T, C> Mul<&'b Transform<T, C, 3>> for UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Transform<T, C, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Translation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Translation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Translation<T, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: &'b Translation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Unit<DualQuaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b UnitDualQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Unit<DualQuaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b UnitDualQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<&'b Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Unit<Vector<T, U3, SB>>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<&'b Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>>> for UnitQuaternion<T>

source§

type Output = Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b Unit<Vector<T, U3, SB>>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, 'b, T: SimdRealField> Mul<&'b Unit<Quaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b UnitQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> Mul<&'b Unit<Quaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: &'b UnitQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Isometry<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<Isometry<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Isometry<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>> for &'a UnitQuaternion<T>

source§

type Output = Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Vector<T, U3, SB>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField, SB: Storage<T, Const<3>>> Mul<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>> for UnitQuaternion<T>

source§

type Output = Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Vector<T, U3, SB>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<OPoint<T, Const<3>>> for &'a UnitQuaternion<T>

source§

type Output = OPoint<T, Const<3>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Point3<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<OPoint<T, Const<3>>> for UnitQuaternion<T>

source§

type Output = OPoint<T, Const<3>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Point3<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Rotation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Rotation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<Rotation<T, 3>> for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Rotation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Similarity<T, Unit<Quaternion<T>>, 3>> for &'a UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<Similarity<T, Unit<Quaternion<T>>, 3>> for UnitQuaternion<T>

source§

type Output = Similarity<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Similarity<T, UnitQuaternion<T>, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T, C> Mul<Transform<T, C, 3>> for &'a UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Transform<T, C, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<T, C> Mul<Transform<T, C, 3>> for UnitQuaternion<T>

source§

type Output = Transform<T, <C as TCategoryMul<TAffine>>::Representative, 3>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Transform<T, C, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Translation<T, 3>> for &'a UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Translation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<Translation<T, 3>> for UnitQuaternion<T>

source§

type Output = Isometry<T, Unit<Quaternion<T>>, 3>

The resulting type after applying the * operator.
source§

fn mul(self, right: Translation<T, 3>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Unit<DualQuaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: UnitDualQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul<Unit<DualQuaternion<T>>> for UnitQuaternion<T>

source§

type Output = Unit<DualQuaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: UnitDualQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField, SB: Storage<T, Const<3>>> Mul<Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Unit<Vector<T, U3, SB>>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField, SB: Storage<T, Const<3>>> Mul<Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, SB>>> for UnitQuaternion<T>

source§

type Output = Unit<Matrix<T, Const<{ typenum::$D::USIZE }>, Const<1>, ArrayStorage<T, 3, 1>>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: Unit<Vector<T, U3, SB>>) -> Self::Output

Performs the * operation. Read more
source§

impl<'a, T: SimdRealField> Mul<Unit<Quaternion<T>>> for &'a UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: UnitQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<T: SimdRealField> Mul for UnitQuaternion<T>

source§

type Output = Unit<Quaternion<T>>

The resulting type after applying the * operator.
source§

fn mul(self, rhs: UnitQuaternion<T>) -> Self::Output

Performs the * operation. Read more
source§

impl<'b, T: SimdRealField> MulAssign<&'b Rotation<T, 3>> for UnitQuaternion<T>

source§

fn mul_assign(&mut self, rhs: &'b Rotation<T, 3>)

Performs the *= operation. Read more
source§

impl<'b, T: SimdRealField> MulAssign<&'b Unit<Quaternion<T>>> for UnitQuaternion<T>

source§

fn mul_assign(&mut self, rhs: &'b UnitQuaternion<T>)

Performs the *= operation. Read more
source§

impl<T: SimdRealField> MulAssign<Rotation<T, 3>> for UnitQuaternion<T>

source§

fn mul_assign(&mut self, rhs: Rotation<T, 3>)

Performs the *= operation. Read more
source§

impl<T: SimdRealField> MulAssign for UnitQuaternion<T>

source§

fn mul_assign(&mut self, rhs: UnitQuaternion<T>)

Performs the *= operation. Read more
source§

impl<T: SimdRealField> One for UnitQuaternion<T>

source§

fn one() -> Self

Returns the multiplicative identity element of Self, 1. Read more
source§

fn set_one(&mut self)

Sets self to the multiplicative identity element of Self, 1.
source§

fn is_one(&self) -> bool
where Self: PartialEq,

Returns true if self is equal to the multiplicative identity. Read more
source§

impl<T: Scalar + ClosedNeg + PartialEq> PartialEq for UnitQuaternion<T>

source§

fn eq(&self, rhs: &Self) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
source§

impl<T: RealField + RelativeEq<Epsilon = T>> RelativeEq for UnitQuaternion<T>

source§

fn default_max_relative() -> Self::Epsilon

The default relative tolerance for testing values that are far-apart. Read more
source§

fn relative_eq( &self, other: &Self, epsilon: Self::Epsilon, max_relative: Self::Epsilon, ) -> bool

A test for equality that uses a relative comparison if the values are far apart.
source§

fn relative_ne( &self, other: &Rhs, epsilon: Self::Epsilon, max_relative: Self::Epsilon, ) -> bool

The inverse of RelativeEq::relative_eq.
source§

impl<T: Scalar + SimdValue> SimdValue for UnitQuaternion<T>
where T::Element: Scalar,

source§

type Element = Unit<Quaternion<<T as SimdValue>::Element>>

The type of the elements of each lane of this SIMD value.
source§

type SimdBool = <T as SimdValue>::SimdBool

Type of the result of comparing two SIMD values like self.
source§

fn lanes() -> usize

The number of lanes of this SIMD value.
source§

fn splat(val: Self::Element) -> Self

Initializes an SIMD value with each lanes set to val.
source§

fn extract(&self, i: usize) -> Self::Element

Extracts the i-th lane of self. Read more
source§

unsafe fn extract_unchecked(&self, i: usize) -> Self::Element

Extracts the i-th lane of self without bound-checking.
source§

fn replace(&mut self, i: usize, val: Self::Element)

Replaces the i-th lane of self by val. Read more
source§

unsafe fn replace_unchecked(&mut self, i: usize, val: Self::Element)

Replaces the i-th lane of self by val without bound-checking.
source§

fn select(self, cond: Self::SimdBool, other: Self) -> Self

Merges self and other depending on the lanes of cond. Read more
source§

fn map_lanes(self, f: impl Fn(Self::Element) -> Self::Element) -> Self
where Self: Clone,

Applies a function to each lane of self. Read more
source§

fn zip_map_lanes( self, b: Self, f: impl Fn(Self::Element, Self::Element) -> Self::Element, ) -> Self
where Self: Clone,

Applies a function to each lane of self paired with the corresponding lane of b. Read more
source§

impl<T1, T2, R> SubsetOf<Isometry<T2, R, 3>> for UnitQuaternion<T1>
where T1: RealField, T2: RealField + SupersetOf<T1>, R: AbstractRotation<T2, 3> + SupersetOf<Self>,

source§

fn to_superset(&self) -> Isometry<T2, R, 3>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(iso: &Isometry<T2, R, 3>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(iso: &Isometry<T2, R, 3>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1: RealField, T2: RealField + SupersetOf<T1>> SubsetOf<Matrix<T2, Const<{ typenum::$D::USIZE }>, Const<{ typenum::$D::USIZE }>, ArrayStorage<T2, 4, 4>>> for UnitQuaternion<T1>

source§

fn to_superset(&self) -> Matrix4<T2>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(m: &Matrix4<T2>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(m: &Matrix4<T2>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1, T2> SubsetOf<Rotation<T2, 3>> for UnitQuaternion<T1>
where T1: RealField, T2: RealField + SupersetOf<T1>,

source§

fn to_superset(&self) -> Rotation3<T2>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(rot: &Rotation3<T2>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(rot: &Rotation3<T2>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1, T2, R> SubsetOf<Similarity<T2, R, 3>> for UnitQuaternion<T1>
where T1: RealField, T2: RealField + SupersetOf<T1>, R: AbstractRotation<T2, 3> + SupersetOf<Self>,

source§

fn to_superset(&self) -> Similarity<T2, R, 3>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(sim: &Similarity<T2, R, 3>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(sim: &Similarity<T2, R, 3>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1, T2, C> SubsetOf<Transform<T2, C, 3>> for UnitQuaternion<T1>

source§

fn to_superset(&self) -> Transform<T2, C, 3>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(t: &Transform<T2, C, 3>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(t: &Transform<T2, C, 3>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1, T2> SubsetOf<Unit<DualQuaternion<T2>>> for UnitQuaternion<T1>
where T1: RealField, T2: RealField + SupersetOf<T1>,

source§

fn to_superset(&self) -> UnitDualQuaternion<T2>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(dq: &UnitDualQuaternion<T2>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(dq: &UnitDualQuaternion<T2>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T1, T2> SubsetOf<Unit<Quaternion<T2>>> for UnitQuaternion<T1>
where T1: Scalar, T2: Scalar + SupersetOf<T1>,

source§

fn to_superset(&self) -> UnitQuaternion<T2>

The inclusion map: converts self to the equivalent element of its superset.
source§

fn is_in_subset(uq: &UnitQuaternion<T2>) -> bool

Checks if element is actually part of the subset Self (and can be converted to it).
source§

fn from_superset_unchecked(uq: &UnitQuaternion<T2>) -> Self

Use with care! Same as self.to_superset but without any property checks. Always succeeds.
source§

fn from_superset(element: &T) -> Option<Self>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
source§

impl<T: RealField + UlpsEq<Epsilon = T>> UlpsEq for UnitQuaternion<T>

source§

fn default_max_ulps() -> u32

The default ULPs to tolerate when testing values that are far-apart. Read more
source§

fn ulps_eq(&self, other: &Self, epsilon: Self::Epsilon, max_ulps: u32) -> bool

A test for equality that uses units in the last place (ULP) if the values are far apart.
source§

fn ulps_ne(&self, other: &Rhs, epsilon: Self::Epsilon, max_ulps: u32) -> bool

The inverse of UlpsEq::ulps_eq.
source§

impl<T: Scalar + ClosedNeg + Eq> Eq for UnitQuaternion<T>