quatica.decomp.qsvd
quatica.decomp.qsvd
Quaternion SVD Implementations for QuatIca
This module provides Q-SVD routines leveraging existing primitives in utils.py:
- Classical Q-SVD via real-block embedding and LAPACK
- Full Classical Q-SVD (no truncation)
- QR decomposition for quaternion matrices
- Randomized Q-SVD (rand_qsvd) using Gaussian sketching + power iterations
- Pass-efficient Q-SVD (pass_eff_qsvd) alternating a single QR per pass
All routines operate on quaternion arrays (numpy.quaternion) and reuse utilities: - real_expand(Q), real_contract(R, m, n) - quat_matmat(A, B), quat_hermitian(A) - quat_frobenius_norm(A)
References: - Ahmadi-Asl, S., Nobakht Kooshkghazi, M., & Leplat, V. (2025). Pass-efficient Randomized Algorithms for Low-rank Approximation of Quaternion Matrices. arXiv:2507.13731
Future Work: - Ma, R.-R., & Bai, Z.-J. (2018). A Structure-Preserving One-Sided Jacobi Method for Computing the SVD of a Quaternion Matrix. arXiv:1811.08671 This paper presents a more advanced, structure-preserving Q-SVD algorithm that will be implemented in future releases for improved accuracy and efficiency.
classical_qsvd(X_quat, R)
Classical Q-SVD via real-block embedding and LAPACK.
Parameters:
X_quat : numpy.ndarray with dtype=quaternion Input quaternion matrix of shape (m, n) R : int Target rank for truncation
Returns:
tuple : (U_quat, s, V_quat) U_quat : m×R quaternion matrix with orthonormal columns s : R-length array of singular values V_quat : n×R quaternion matrix with orthonormal columns
Notes:
Complexity: O((4m)(4n)min(4m,4n)), leverages optimized LAPACK. Handles degenerate singular values from real-block embedding by extracting every 4th singular value to obtain true quaternion singular values.
For reconstruction: X ≈ U @ diag(s) @ V^H (when R < min(m,n))
classical_qsvd_full(X_quat)
Full Classical Q-SVD via real-block embedding and LAPACK.
Parameters:
X_quat : numpy.ndarray with dtype=quaternion Input quaternion matrix of shape (m, n)
Returns:
tuple : (U_quat, s, V_quat) U_quat : m×m quaternion matrix with orthonormal columns s : min(m,n)-length array of singular values V_quat : n×n quaternion matrix with orthonormal columns
Notes:
Returns the FULL Q-SVD decomposition without truncation. For reconstruction: X = U @ Σ @ V^H where Σ is the diagonal matrix with s.
pass_eff_qsvd(X_quat, R, oversample=10, n_passes=2)
Pass-efficient Q-SVD alternating single multiplies + QR.
Parameters:
X_quat : numpy.ndarray with dtype=quaternion Input quaternion matrix of shape (m, n) R : int Target rank oversample : int, optional Oversampling parameter (default: 10) n_passes : int, optional Number of passes over the matrix (default: 2)
Returns:
tuple : (U_quat, s, V_quat) U_quat : m×R quaternion matrix with orthonormal columns s : R-length array of singular values V_quat : n×R quaternion matrix with orthonormal columns
Notes:
Advantage: fewer large passes over X, good cache behavior.
qr_qua(X_quat)
QR decomposition of quaternion matrix using real-block embedding.
Parameters:
X_quat : numpy.ndarray with dtype=quaternion Input quaternion matrix of shape (m, n)
Returns:
tuple : (Q_quat, R_quat) Q_quat : m×n quaternion matrix with orthonormal columns R_quat : n×n upper triangular quaternion matrix
Notes:
Uses real-block embedding + SciPy's QR + contraction back to quaternion. Produces perfect reconstruction and orthonormal Q matrix.
rand_qsvd(X_quat, R, oversample=10, n_iter=2)
Randomized Q-SVD using Gaussian sketching + power iterations.
Parameters:
X_quat : numpy.ndarray with dtype=quaternion Input quaternion matrix of shape (m, n) R : int Target rank oversample : int, optional Oversampling parameter (default: 10) n_iter : int, optional Number of power iterations (default: 2)
Returns:
tuple : (U_quat, s, V_quat) U_quat : m×R quaternion matrix with orthonormal columns s : R-length array of singular values V_quat : n×R quaternion matrix with orthonormal columns
Notes:
Cost: O(mn(R+P)) + O((R+P)²n) where P = oversample.