procrustes.rotational
Rotational-Orthogonal Procrustes Module.
- procrustes.rotational.rotational(a: ndarray, b: ndarray, pad: bool = True, translate: bool = False, scale: bool = False, unpad_col: bool = False, unpad_row: bool = False, check_finite: bool = True, weight: ndarray | None = None, lapack_driver: str = 'gesvd') ProcrustesResult [source]
Perform rotational Procrustes.
Given a matrix \(\mathbf{A}_{m \times n}\) and a reference matrix \(\mathbf{B}_{m \times n}\), find the rotational transformation matrix \(\mathbf{R}_{n \times n}\) that makes \(\mathbf{A}\) as close as possible to \(\mathbf{B}\). In other words,
\[\underbrace{\min}_{\left\{\mathbf{R} \left| {\mathbf{R}^{-1} = {\mathbf{R}}^\dagger \atop \left| \mathbf{R} \right| = 1} \right. \right\}} \|\mathbf{A}\mathbf{R} - \mathbf{B}\|_{F}^2\]This Procrustes method requires the \(\mathbf{A}\) and \(\mathbf{B}\) matrices to have the same shape, which is gauranteed with the default
pad
argument for any given \(\mathbf{A}\) and \(\mathbf{B}\) matrices. In preparing the \(\mathbf{A}\) and \(\mathbf{B}\) matrices, the (optional) order of operations is: 1) unpad zero rows/columns, 2) translate the matrices to the origin, 3) weight entries of \(\mathbf{A}\), 4) scale the matrices to have unit norm, 5) pad matrices with zero rows/columns so they have the same shape.- Parameters:
a (ndarray) -- The 2D-array \(\mathbf{A}\) which is going to be transformed.
b (ndarray) -- The 2D-array \(\mathbf{B}\) representing the reference matrix.
pad (bool, optional) -- Add zero rows (at the bottom) and/or columns (to the right-hand side) of matrices \(\mathbf{A}\) and \(\mathbf{B}\) so that they have the same shape.
translate (bool, optional) -- If True, both arrays are centered at origin (columns of the arrays will have mean zero).
scale (bool, optional) -- If True, both arrays are normalized with respect to the Frobenius norm, i.e., \(\text{Tr}\left[\mathbf{A}^\dagger\mathbf{A}\right] = 1\) and \(\text{Tr}\left[\mathbf{B}^\dagger\mathbf{B}\right] = 1\).
unpad_col (bool, optional) -- If True, zero columns (with values less than 1.0e-8) on the right-hand side of the intial \(\mathbf{A}\) and \(\mathbf{B}\) matrices are removed.
unpad_row (bool, optional) -- If True, zero rows (with values less than 1.0e-8) at the bottom of the intial \(\mathbf{A}\) and \(\mathbf{B}\) matrices are removed.
check_finite (bool, optional) -- If True, convert the input to an array, checking for NaNs or Infs.
weight (ndarray, optional) -- The 1D-array representing the weights of each row of \(\mathbf{A}\). This defines the elements of the diagonal matrix \(\mathbf{W}\) that is multiplied by \(\mathbf{A}\) matrix, i.e., \(\mathbf{A} \rightarrow \mathbf{WA}\).
lapack_driver ({'gesvd', 'gesdd'}, optional) -- Whether to use the more efficient divide-and-conquer approach ('gesdd') or the more robust general rectangular approach ('gesvd') to compute the singular-value decomposition with scipy.linalg.svd.
- Returns:
res -- The Procrustes result represented as a class:utils.ProcrustesResult object.
- Return type:
Notes
The optimal rotational matrix is obtained by,
\[\mathbf{R}_{\text{opt}} = \arg \underbrace{\min}_{\left\{\mathbf{R} \left| {\mathbf{R}^{-1} = {\mathbf{R}}^\dagger \atop \left| \mathbf{R} \right| = 1} \right. \right\}} \|\mathbf{A}\mathbf{R} - \mathbf{B}\|_{F}^2 = \arg \underbrace{\max}_{\left\{\mathbf{R} \left| {\mathbf{R}^{-1} = {\mathbf{R}}^\dagger \atop \left| \mathbf{R} \right| = 1} \right. \right\}} \text{Tr}\left[\mathbf{R}^\dagger {\mathbf{A}}^\dagger \mathbf{B} \right]\]The solution is obtained by taking the singular value decomposition (SVD) of the \(\mathbf{A}^\dagger \mathbf{B}\) matrix,
\[\begin{split}\mathbf{A}^\dagger \mathbf{B} &= \tilde{\mathbf{U}} \tilde{\mathbf{\Sigma}} \tilde{\mathbf{V}}^{\dagger} \\ \mathbf{R}_{\text{opt}} &= \tilde{\mathbf{U}} \tilde{\mathbf{S}} \tilde{\mathbf{V}}^{\dagger}\end{split}\]where \(\tilde{\mathbf{S}}_{n \times m}\) is almost an identity matrix,
\[\begin{split}\tilde{\mathbf{S}}_{m \times n} \equiv \begin{bmatrix} 1 & 0 & \cdots & 0 & 0 \\ 0 & 1 & \ddots & \vdots &0 \\ 0 & \ddots &\ddots & 0 &\vdots \\ \vdots&0 & 0 & 1 &0 \\ 0 & 0 & 0 \cdots &0 &\operatorname{sgn} \left(\left|\mathbf{U}\mathbf{V}^\dagger\right|\right) \end{bmatrix}\end{split}\]in which the smallest singular value is replaced by
\[\begin{split}\operatorname{sgn} \left(\left|\tilde{\mathbf{U}} \tilde{\mathbf{V}}^\dagger\right|\right) = \begin{cases} +1 \qquad \left|\tilde{\mathbf{U}} \tilde{\mathbf{V}}^\dagger\right| \geq 0 \\ -1 \qquad \left|\tilde{\mathbf{U}} \tilde{\mathbf{V}}^\dagger\right| < 0 \end{cases}\end{split}\]Examples
>>> import numpy as np >>> array_a = np.array([[1.5, 7.4], [8.5, 4.5]]) >>> array_b = np.array([[6.29325035, 4.17193001, 0., 0,], ... [9.19238816, -2.82842712, 0., 0.], ... [0., 0., 0., 0.]]) >>> res = rotational(array_a,array_b,translate=False,scale=False) >>> res.t # rotational transformation array([[ 0.70710678, -0.70710678], [ 0.70710678, 0.70710678]]) >>> res.error # one-sided Procrustes error 1.483808210011695e-17