Matrix Functions

#include <Imath/ImathMatrixAlgo.h>

Functions that operate on matrices

template<class T>
bool Imath::extractScaling(const Matrix44<T> &mat, Vec3<T> &scl, bool exc = true)

Extract the scaling component of the given 4x4 matrix.

Parameters:
  • mat[in] The input matrix

  • scl[out] The extracted scale, i.e. the output value

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
Matrix44<T> Imath::sansScaling(const Matrix44<T> &mat, bool exc = true)

Return the given 4x4 matrix with scaling removed.

Parameters:
  • mat[in] The input matrix

  • exc[in] If true, throw an exception if the scaling in mat

template<class T>
bool Imath::removeScaling(Matrix44<T> &mat, bool exc = true)

Remove scaling from the given 4x4 matrix in place.

Return true if the scale could be successfully extracted, false if the matrix is degenerate.

Parameters:
  • mat[in] The matrix to operate on

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractScalingAndShear(const Matrix44<T> &mat, Vec3<T> &scl, Vec3<T> &shr, bool exc = true)

Extract the scaling and shear components of the given 4x4 matrix.

Return true if the scale could be successfully extracted, false if the matrix is degenerate.

Parameters:
  • mat[in] The input matrix

  • scl[out] The extracted scale

  • shr[out] The extracted shear

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
Matrix44<T> Imath::sansScalingAndShear(const Matrix44<T> &mat, bool exc = true)

Return the given 4x4 matrix with scaling and shear removed.

Parameters:
  • mat[in] The input matrix

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

template<class T>
void Imath::sansScalingAndShear(Matrix44<T> &result, const Matrix44<T> &mat, bool exc = true)

Extract scaling and shear from the given 4x4 matrix in-place.

Parameters:
  • result[inout] The output matrix

  • mat[in] The return value if result is degenerate

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

template<class T>
bool Imath::removeScalingAndShear(Matrix44<T> &mat, bool exc = true)

Remove scaling and shear from the given 4x4 matrix in place.

Parameters:
  • mat[inout] The matrix to operate on

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractAndRemoveScalingAndShear(Matrix44<T> &mat, Vec3<T> &scl, Vec3<T> &shr, bool exc = true)

Remove scaling and shear from the given 4x4 matrix in place, returning the extracted values.

Parameters:
  • mat[inout] The matrix to operate on

  • scl[out] The extracted scale

  • shr[out] The extracted shear

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
void Imath::extractEulerXYZ(const Matrix44<T> &mat, Vec3<T> &rot)

Extract the rotation from the given 4x4 matrix in the form of XYZ euler angles.

Parameters:
  • mat[in] The input matrix

  • rot[out] The extracted XYZ euler angle vector

template<class T>
void Imath::extractEulerZYX(const Matrix44<T> &mat, Vec3<T> &rot)

Extract the rotation from the given 4x4 matrix in the form of ZYX euler angles.

Parameters:
  • mat[in] The input matrix

  • rot[out] The extracted ZYX euler angle vector

template<class T>
Quat<T> Imath::extractQuat(const Matrix44<T> &mat)

Extract the rotation from the given 4x4 matrix in the form of a quaternion.

Parameters:

mat[in] The input matrix

Returns:

The extracted quaternion

template<class T>
bool Imath::extractSHRT(const Matrix44<T> &mat, Vec3<T> &s, Vec3<T> &h, Vec3<T> &r, Vec3<T> &t, bool exc, typename Euler<T>::Order rOrder)

Extract the scaling, shear, rotation, and translation components of the given 4x4 matrix.

The values are such that:

M = S * H * R * T

Parameters:
  • mat[in] The input matrix

  • s[out] The extracted scale

  • h[out] The extracted shear

  • r[out] The extracted rotation

  • t[out] The extracted translation

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

  • rOrder[in] The order with which to extract the rotation

Returns:

True if the values could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractSHRT(const Matrix44<T> &mat, Vec3<T> &s, Vec3<T> &h, Vec3<T> &r, Vec3<T> &t, bool exc = true)

Extract the scaling, shear, rotation, and translation components of the given 4x4 matrix.

Parameters:
  • mat[in] The input matrix

  • s[out] The extracted scale

  • h[out] The extracted shear

  • r[out] The extracted rotation, in XYZ euler angles

  • t[out] The extracted translation

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the values could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractSHRT(const Matrix44<T> &mat, Vec3<T> &s, Vec3<T> &h, Euler<T> &r, Vec3<T> &t, bool exc = true)

Extract the scaling, shear, rotation, and translation components of the given 4x4 matrix.

Parameters:
  • mat[in] The input matrix

  • s[out] The extracted scale

  • h[out] The extracted shear

  • r[out] The extracted rotation, in Euler angles

  • t[out] The extracted translation

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the values could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::checkForZeroScaleInRow(const T &scl, const Vec3<T> &row, bool exc = true)

Return true if the given scale can be removed from the given row matrix, false if scl is small enough that the operation would overflow.

If exc is true, throw an exception on overflow.

template<class T>
Matrix44<T> Imath::outerProduct(const Vec4<T> &a, const Vec4<T> &b)

Return the 4x4 outer product two 4-vectors.

template<class T>
Matrix44<T> Imath::rotationMatrix(const Vec3<T> &fromDirection, const Vec3<T> &toDirection)

Return a 4x4 matrix that rotates the vector fromDirection to toDirection

template<class T>
Matrix44<T> Imath::rotationMatrixWithUpDir(const Vec3<T> &fromDir, const Vec3<T> &toDir, const Vec3<T> &upDir)

Return a 4x4 matrix that rotates the fromDir vector so that it points towards `toDir1.

You may also specify that you want the up vector to be pointing in a certain direction 1upDir`.

template<class T>
void Imath::alignZAxisWithTargetDir(Matrix44<T> &result, Vec3<T> targetDir, Vec3<T> upDir)

Construct a 4x4 matrix that rotates the z-axis so that it points towards targetDir.

You must also specify that you want the up vector to be pointing in a certain direction upDir.

Notes: The following degenerate cases are handled: (a) when the directions given by toDir and upDir are parallel or opposite (the direction vectors must have a non-zero cross product); (b) when any of the given direction vectors have zero length

Parameters:
  • result[out] The output matrix

  • targetDir[in] The target direction vector

  • upDir[in] The up direction vector

template<class T>
Matrix44<T> Imath::computeLocalFrame(const Vec3<T> &p, const Vec3<T> &xDir, const Vec3<T> &normal)

Compute an orthonormal direct 4x4 frame from a position, an x axis direction and a normal to the y axis.

If the x axis and normal are perpendicular, then the normal will have the same direction as the z axis.

Parameters:
  • p[in] The position of the frame

  • xDir[in] The x axis direction of the frame

  • normal[in] A normal to the y axis of the frame

Returns:

The orthonormal frame

template<class T>
Matrix44<T> Imath::addOffset(const Matrix44<T> &inMat, const Vec3<T> &tOffset, const Vec3<T> &rOffset, const Vec3<T> &sOffset, const Vec3<T> &ref)

Add a translate/rotate/scale offset to a 4x4 input frame and put it in another frame of reference.

Parameters:
  • inMat[in] Input frame

  • tOffset[in] Translation offset

  • rOffset[in] Rotation offset in degrees

  • sOffset[in] Scale offset

  • ref[in] Frame of reference

Returns:

The offsetted frame

template<class T>
Matrix44<T> Imath::computeRSMatrix(bool keepRotateA, bool keepScaleA, const Matrix44<T> &A, const Matrix44<T> &B)

Compute 4x4 translate/rotate/scale matrix from A with the rotate/scale of B.

Parameters:
  • keepRotateA[in] If true, keep rotate from matrix A, use B otherwise

  • keepScaleA[in] If true, keep scale from matrix A, use B otherwise

  • A[in] Matrix A

  • B[in] Matrix B

Returns:

Matrix A with tweaked rotation/scale

template<class T>
bool Imath::extractScaling(const Matrix33<T> &mat, Vec2<T> &scl, bool exc = true)

Extract the scaling component of the given 3x3 matrix.

Parameters:
  • mat[in] The input matrix

  • scl[out] The extracted scale, i.e. the output value

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
Matrix33<T> Imath::sansScaling(const Matrix33<T> &mat, bool exc = true)

Return the given 3x3 matrix with scaling removed.

Parameters:
  • mat[in] The input matrix

  • exc[in] If true, throw an exception if the scaling in mat

template<class T>
bool Imath::removeScaling(Matrix33<T> &mat, bool exc = true)

Remove scaling from the given 3x3 matrix in place.

Return true if the scale could be successfully extracted, false if the matrix is degenerate.

Parameters:
  • mat[in] The matrix to operate on

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractScalingAndShear(const Matrix33<T> &mat, Vec2<T> &scl, T &shr, bool exc = true)

Extract the scaling and shear components of the given 3x3 matrix.

Return true if the scale could be successfully extracted, false if the matrix is degenerate.

Parameters:
  • mat[in] The input matrix

  • scl[out] The extracted scale

  • shr[out] The extracted shear

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
Matrix33<T> Imath::sansScalingAndShear(const Matrix33<T> &mat, bool exc = true)

Return the given 3x3 matrix with scaling and shear removed.

Parameters:
  • mat[in] The input matrix

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

template<class T>
bool Imath::removeScalingAndShear(Matrix33<T> &mat, bool exc = true)

Remove scaling and shear from the given 3x3e matrix in place.

Parameters:
  • mat[inout] The matrix to operate on

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::extractAndRemoveScalingAndShear(Matrix33<T> &mat, Vec2<T> &scl, T &shr, bool exc = true)

Remove scaling and shear from the given 3x3 matrix in place, returning the extracted values.

Parameters:
  • mat[inout] The matrix to operate on

  • scl[out] The extracted scale

  • shr[out] The extracted shear

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the scale could be extracted, false if the matrix is degenerate.

template<class T>
void Imath::extractEuler(const Matrix22<T> &mat, T &rot)

Extract the rotation from the given 2x2 matrix.

Parameters:
  • mat[in] The input matrix

  • rot[out] The extracted rotation value

template<class T>
void Imath::extractEuler(const Matrix33<T> &mat, T &rot)

Extract the rotation from the given 3x3 matrix.

Parameters:
  • mat[in] The input matrix

  • rot[out] The extracted rotation value

template<class T>
bool Imath::extractSHRT(const Matrix33<T> &mat, Vec2<T> &s, T &h, T &r, Vec2<T> &t, bool exc = true)

Extract the scaling, shear, rotation, and translation components of the given 3x3 matrix.

The values are such that:

M = S * H * R * T

Parameters:
  • mat[in] The input matrix

  • s[out] The extracted scale

  • h[out] The extracted shear

  • r[out] The extracted rotation

  • t[out] The extracted translation

  • exc[in] If true, throw an exception if the scaling in mat is very close to zero.

Returns:

True if the values could be extracted, false if the matrix is degenerate.

template<class T>
bool Imath::checkForZeroScaleInRow(const T &scl, const Vec2<T> &row, bool exc = true)

Return true if the given scale can be removed from the given row matrix, false if scl is small enough that the operation would overflow.

If exc is true, throw an exception on overflow.

template<class T>
Matrix33<T> Imath::outerProduct(const Vec3<T> &a, const Vec3<T> &b)

Return the 3xe outer product two 3-vectors.

template<typename T>
M44d Imath::procrustesRotationAndTranslation(const Vec3<T> *A, const Vec3<T> *B, const T *weights, const size_t numPoints, const bool doScaling = false)

Computes the translation and rotation that brings the ‘from’ points as close as possible to the ‘to’ points under the Frobenius norm.

To be more specific, let x be the matrix of ‘from’ points and y be the matrix of ‘to’ points, we want to find the matrix A of the form [ R t ] [ 0 1 ] that minimizes || (A*x - y)^T * W * (A*x - y) ||_F If doScaling is true, then a uniform scale is allowed also.

Parameters:
  • A – From points

  • B – To points

  • weights – Per-point weights

  • numPoints – The number of points in A, B, and weights (must be equal)

  • doScaling – If true, include a scaling transformation

Returns:

The procrustes transformation

template<typename T>
M44d Imath::procrustesRotationAndTranslation(const Vec3<T> *A, const Vec3<T> *B, const size_t numPoints, const bool doScaling = false)

Computes the translation and rotation that brings the ‘from’ points as close as possible to the ‘to’ points under the Frobenius norm.

To be more specific, let x be the matrix of ‘from’ points and y be the matrix of ‘to’ points, we want to find the matrix A of the form [ R t ] [ 0 1 ] that minimizes || (A*x - y)^T * W * (A*x - y) ||_F If doScaling is true, then a uniform scale is allowed also.

Parameters:
  • A – From points

  • B – To points

  • numPoints – The number of points in A and B (must be equal)

  • doScaling – If true, include a scaling transformation

Returns:

The procrustes transformation

template<typename T>
void Imath::jacobiSVD(const Matrix33<T> &A, Matrix33<T> &U, Vec3<T> &S, Matrix33<T> &V, const T tol = std::numeric_limits<T>::epsilon(), const bool forcePositiveDeterminant = false)

Compute the SVD of a 3x3 matrix using Jacobi transformations.

This method should be quite accurate (competitive with LAPACK) even for poorly conditioned matrices, and because it has been written specifically for the 3x3/4x4 case it is much faster than calling out to LAPACK.

The SVD of a 3x3/4x4 matrix A is defined as follows: A = U * S * V^T where S is the diagonal matrix of singular values and both U and V are orthonormal. By convention, the entries S are all positive and sorted from the largest to the smallest. However, some uses of this function may require that the matrix U*V^T have positive determinant; in this case, we may make the smallest singular value negative to ensure that this is satisfied.

Currently only available for single- and double-precision matrices.

template<typename T>
void Imath::jacobiSVD(const Matrix44<T> &A, Matrix44<T> &U, Vec4<T> &S, Matrix44<T> &V, const T tol = std::numeric_limits<T>::epsilon(), const bool forcePositiveDeterminant = false)

Compute the SVD of a 3x3 matrix using Jacobi transformations.

This method should be quite accurate (competitive with LAPACK) even for poorly conditioned matrices, and because it has been written specifically for the 3x3/4x4 case it is much faster than calling out to LAPACK.

The SVD of a 3x3/4x4 matrix A is defined as follows: A = U * S * V^T where S is the diagonal matrix of singular values and both U and V are orthonormal. By convention, the entries S are all positive and sorted from the largest to the smallest. However, some uses of this function may require that the matrix U*V^T have positive determinant; in this case, we may make the smallest singular value negative to ensure that this is satisfied.

Currently only available for single- and double-precision matrices.

template<typename T>
void Imath::jacobiEigenSolver(Matrix33<T> &A, Vec3<T> &S, Matrix33<T> &V, const T tol)

Compute the eigenvalues (S) and the eigenvectors (V) of a real symmetric matrix using Jacobi transformation, using a given tolerance tol.

Jacobi transformation of a 3x3/4x4 matrix A outputs S and V: A = V * S * V^T where V is orthonormal and S is the diagonal matrix of eigenvalues. Input matrix A must be symmetric. A is also modified during the computation so that upper diagonal entries of A become zero.

template<typename T>
inline void Imath::jacobiEigenSolver(Matrix33<T> &A, Vec3<T> &S, Matrix33<T> &V)

Compute the eigenvalues (S) and the eigenvectors (V) of a real symmetric matrix using Jacobi transformation.

Jacobi transformation of a 3x3/4x4 matrix A outputs S and V: A = V * S * V^T where V is orthonormal and S is the diagonal matrix of eigenvalues. Input matrix A must be symmetric. A is also modified during the computation so that upper diagonal entries of A become zero.

template<typename T>
void Imath::jacobiEigenSolver(Matrix44<T> &A, Vec4<T> &S, Matrix44<T> &V, const T tol)

Compute the eigenvalues (S) and the eigenvectors (V) of a real symmetric matrix using Jacobi transformation, using a given tolerance tol.

Jacobi transformation of a 3x3/4x4 matrix A outputs S and V: A = V * S * V^T where V is orthonormal and S is the diagonal matrix of eigenvalues. Input matrix A must be symmetric. A is also modified during the computation so that upper diagonal entries of A become zero.

template<typename T>
inline void Imath::jacobiEigenSolver(Matrix44<T> &A, Vec4<T> &S, Matrix44<T> &V)

Compute the eigenvalues (S) and the eigenvectors (V) of a real symmetric matrix using Jacobi transformation.

Jacobi transformation of a 3x3/4x4 matrix A outputs S and V: A = V * S * V^T where V is orthonormal and S is the diagonal matrix of eigenvalues. Input matrix A must be symmetric. A is also modified during the computation so that upper diagonal entries of A become zero.

template<typename TM, typename TV>
void Imath::maxEigenVector(TM &A, TV &S)

Compute a eigenvector corresponding to the abs max eigenvalue of a real symmetric matrix using Jacobi transformation.

template<typename TM, typename TV>
void Imath::minEigenVector(TM &A, TV &S)

Compute a eigenvector corresponding to the abs min eigenvalue of a real symmetric matrix using Jacobi transformation.