procrustes.generic

Generic Procrustes Module.

procrustes.generic.generic(a: numpy.ndarray, b: numpy.ndarray, pad: bool = True, translate: bool = False, scale: bool = False, unpad_col: bool = False, unpad_row: bool = False, check_finite: bool = True, weight: Optional[numpy.ndarray] = None, use_svd: bool = False) procrustes.utils.ProcrustesResult[source]

Perform generic one-sided Procrustes.

Given matrix \(\mathbf{A}_{m \times n}\) and a reference matrix \(\mathbf{B}_{m \times n}\), find the transformation matrix \(\mathbf{T}_{n \times n}\) that makes \(\mathbf{AT}\) as close as possible to \(\mathbf{B}\). In other words,

\[\underbrace{\text{min}}_{\mathbf{T}} \quad \|\mathbf{A} \mathbf{T} - \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}\).

  • use_svd (bool, optional) -- If True, the (Moore-Penrose) pseudo-inverse is computed by singular-value decomposition (SVD) including all 'large' singular values (using scipy.linalg.pinv2). If False, the the (Moore-Penrose) pseudo-inverse is computed by least-squares solver (using scipy.linalg.pinv). The least-squares implementation is less efficient, but more robust, than the SVD implementation.

Returns

res -- The Procrustes result represented as a class:utils.ProcrustesResult object.

Return type

ProcrustesResult

Notes

The optimal transformation matrix is obtained by solving the least-squares equations,

\[\mathbf{X}_\text{opt} = {(\mathbf{A}^{\top}\mathbf{A})}^{-1} \mathbf{A}^{\top} \mathbf{B}\]

If \(m < n\), the transformation matrix \(\mathbf{T}_\text{opt}\) is not unique, because the system of equations is underdetermined (i.e., there are fewer equations than unknowns).