t3toolbox.backend.linalg.truncated_svd#
- t3toolbox.backend.linalg.truncated_svd(A: NDArray, min_rank: int = None, max_rank: int = None, rtol: float = None, atol: float = None, use_jax: bool = False) t3toolbox.backend.common.typ.Tuple[NDArray, NDArray, NDArray]#
Compute (truncated) singular value decomposition of matrix.
A = U @ diag(ss) @ Vt Equality may be approximate if truncation is used.
- Parameters:
A (NDArray) – Matrix. shape=(N, M)
min_rank (int) – Minimum rank for truncation. Should have 1 <= min_rank <= max_rank <= minimum(N, M).
min_rank – Maximum rank for truncation. Should have 1 <= min_rank <= max_rank <= minimum(N, M).
rtol (float) – Relative tolerance for truncation. Remove singular values satisfying sigma < maximum(atol, rtol*sigma1).
atol (float) – Absolute tolerance for truncation. Remove singular values satisfying sigma < maximum(atol, rtol*sigma1).
xnp – Linear algebra backend. Default: np (numpy)
- Returns:
U (NDArray) – Left singular vectors. shape=(N, k). U.T @ U = identity matrix
ss (NDArray) – Singular values. Non-negative. shape=(k,).
Vt (NDArray) – Right singular vectors. shape=(k, M) Vt @ Vt.T = identity matrix
Examples
>>> import numpy as np >>> import t3toolbox.dense as dense >>> A = np.random.randn(55,70) >>> U, ss, Vt = dense.truncated_svd(A) >>> A2 = np.einsum('ix,x,xj->ij', U, ss, Vt) >>> print(np.linalg.norm(A - A2)) 1.0428742517412705e-13 >>> rank = len(ss) >>> print(np.linalg.norm(U.T @ U - np.eye(rank))) 1.1907994177245428e-14 >>> print(np.linalg.norm(Vt @ Vt.T - np.eye(rank))) 1.1027751835566194e-14
>>> import numpy as np >>> import t3toolbox.dense as dense >>> A = np.random.randn(55, 70) @ np.diag(1.0 / np.arange(1,71)**2) # Create matrix with spectral decay >>> U, ss, Vt = dense.truncated_svd(A, rtol=1e-2) # Truncated SVD with relative tolerance 1e-2 >>> A2 = np.einsum('ix,x,xj->ij', U, ss, Vt) >>> truncated_rank = len(ss) >>> print(truncated_rank) 10 >>> relerr_num = np.linalg.norm(A - A2, 2) # Check error in induced 2-norm >>> relerr_den = np.linalg.norm(A, 2) >>> print(relerr_num / relerr_den) # should be just less than rtol=1e-2 0.008530627920514714 >>> U, ss, Vt = dense.truncated_svd(A, atol=1e-2) # Truncated SVD with absolute tolerance 1e-2 >>> A2 = np.einsum('ix,x,xj->ij', U, ss, Vt) >>> truncated_rank = len(ss) >>> print(truncated_rank) 24 >>> err = np.linalg.norm(A - A2, 2) # Check error in induced 2-norm >>> print(err) # should be just less than atol=1e-2 0.00882416786402483