bundles / scipy latest / scipy / sparse / linalg / _eigen / _svds / svds
function
scipy.sparse.linalg._eigen._svds:svds
Signature
def svds ( A , k = 6 , ncv = None , tol = 0 , which = LM , v0 = None , maxiter = None , return_singular_vectors = True , solver = arpack , rng = None , options = None , * , random_state = None ) Summary
Partial singular value decomposition of a sparse matrix.
Extended Summary
Compute the largest or smallest k singular values and corresponding singular vectors of a sparse matrix A. The order in which the singular values are returned is not guaranteed.
In the descriptions below, let M, N = A.shape.
Parameters
A: ndarray, sparse matrix, or LinearOperatorMatrix to decompose of a floating point numeric dtype.
k: int, default: 6Number of singular values and singular vectors to compute. Must satisfy
1 <= k <= kmax, wherekmax=min(M, N)forsolver='propack'andkmax=min(M, N) - 1otherwise.ncv: int, optionalWhen
solver='arpack', this is the number of Lanczos vectors generated. See'arpack' <sparse.linalg.svds-arpack>for details. Whensolver='lobpcg'orsolver='propack', this parameter is ignored.tol: float, optionalTolerance for singular values. Zero (default) means machine precision.
which: {'LM', 'SM'}Which
ksingular values to find: either the largest magnitude ('LM') or smallest magnitude ('SM') singular values.v0: ndarray, optionalThe starting vector for iteration; see method-specific documentation (
'arpack' <sparse.linalg.svds-arpack>,'lobpcg' <sparse.linalg.svds-lobpcg>), or'propack' <sparse.linalg.svds-propack>for details.maxiter: int, optionalMaximum number of iterations; see method-specific documentation (
'arpack' <sparse.linalg.svds-arpack>,'lobpcg' <sparse.linalg.svds-lobpcg>), or'propack' <sparse.linalg.svds-propack>for details.return_singular_vectors: {True, False, "u", "vh"}Singular values are always computed and returned; this parameter controls the computation and return of singular vectors.
True: return singular vectors.False: do not return singular vectors."u": ifM <= N, compute only the left singular vectors and returnNonefor the right singular vectors. Otherwise, compute all singular vectors."vh": ifM > N, compute only the right singular vectors and returnNonefor the left singular vectors. Otherwise, compute all singular vectors.
If
solver='propack', the option is respected regardless of the matrix shape.solver: {'arpack', 'propack', 'lobpcg'}, optionalThe solver used.
'arpack' <sparse.linalg.svds-arpack>,'lobpcg' <sparse.linalg.svds-lobpcg>, and'propack' <sparse.linalg.svds-propack>are supported. Default:'arpack'.rng: {None, int, `numpy.random.Generator`}, optionalIf
rngis passed by keyword, types other than numpy.random.Generator are passed to numpy.random.default_rng to instantiate aGenerator. Ifrngis already aGeneratorinstance, then the provided instance is used. Specifyrngfor repeatable function behavior.If this argument is passed by position or
random_stateis passed by keyword, legacy behavior for the argumentrandom_stateapplies:If
random_stateis None (or numpy.random), the numpy.random.RandomState singleton is used.If
random_stateis an int, a newRandomStateinstance is used, seeded withrandom_state.If
random_stateis already aGeneratororRandomStateinstance then that instance is used.
options: dict, optionalA dictionary of solver-specific options. No solver-specific options are currently supported; this parameter is reserved for future use.
Returns
u: ndarray, shape=(M, k)Unitary matrix having left singular vectors as columns.
s: ndarray, shape=(k,)The singular values.
vh: ndarray, shape=(k, N)Unitary matrix having right singular vectors as rows.
Notes
When ARPACK or LOBPCG is selected as the method, the singular values and singular vectors are computed as the eigenvalues and eigenvectors of the corresponding Gram matrix, either A.conj().T @ A or A @ A.conj().T, depending on which one is computationally cheaper. It is then followed by the Rayleigh-Ritz method as postprocessing; see Using the normal matrix, in Rayleigh-Ritz method, (2022, Nov. 19), Wikipedia, https://w.wiki/4zms.
Alternatively, the PROPACK solver can be called.
Dtype of A is mapped, if possible, to float(32,64) and complex(64,128) for internal computations.
Examples
Construct a matrix `A` from singular values and vectors.import numpy as np from scipy import sparse, linalg, stats from scipy.sparse.linalg import svds, aslinearoperator, LinearOperator✓
rng = np.random.default_rng(258265244568965474821194062361901728911) orthogonal = stats.ortho_group.rvs(10, random_state=rng) s = [1e-3, 1, 2, 3, 4] # non-zero singular values u = orthogonal[:, :5] # left singular vectors vT = orthogonal[:, 5:].T # right singular vectors A = u @ np.diag(s) @ vT✓
u4, s4, vT4 = svds(A, k=4) A4 = u4 @ np.diag(s4) @ vT4 np.allclose(A4, A, atol=1e-3)✓
u5, s5, vT5 = svds(A, k=5) A5 = u5 @ np.diag(s5) @ vT5 np.allclose(A5, A)✓
np.allclose(s5, s)
✓(np.allclose(np.abs(u5), np.abs(u)) and np.allclose(np.abs(vT5), np.abs(vT)))✓
(np.allclose(u5.T @ u5, np.eye(5)) and np.allclose(vT5 @ vT5.T, np.eye(5)))✓
rng = np.random.default_rng(178686584221410808734965903901790843963) s = [1, 1 + 1e-6] # non-zero singular values u, _ = np.linalg.qr(rng.standard_normal((99, 2))) v, _ = np.linalg.qr(rng.standard_normal((99, 2))) vT = v.T A = u @ np.diag(s) @ vT A = A.astype(np.float32) u2, s2, vT2 = svds(A, k=2, rng=rng) np.allclose(s2, s)✓
(linalg.subspace_angles(u2[:, :1], u[:, :1]) + linalg.subspace_angles(u2[:, 1:], u[:, 1:])) (linalg.subspace_angles(vT2[:1, :].T, vT[:1, :].T) + linalg.subspace_angles(vT2[1:, :].T, vT[1:, :].T))✗
linalg.subspace_angles(u2, u).sum() < 1e-6
✗linalg.subspace_angles(vT2.T, vT.T).sum() < 1e-6
✗rng = np.random.default_rng(0) X_dense = rng.random(size=(100, 100)) X_dense[:, 2 * np.arange(50)] = 0 X = sparse.csr_array(X_dense) _, singular_values, _ = svds(X, k=5, rng=rng) print(singular_values)✓
rng = np.random.default_rng(102524723947864966825913730119128190974) G = sparse.random_array((8, 9), density=0.5, rng=rng) Glo = aslinearoperator(G) _, singular_values_svds, _ = svds(Glo, k=5, rng=rng) _, singular_values_svd, _ = linalg.svd(G.toarray()) np.allclose(singular_values_svds, singular_values_svd[-4::-1])✓
diff0 = lambda a: np.diff(a, axis=0)
✓n = 5 # The dimension of the space. M_from_diff0 = diff0(np.eye(n)) print(M_from_diff0.astype(int))✓
M = - np.eye(n - 1, n, dtype=int) np.fill_diagonal(M[:,1:], 1) np.allclose(M, M_from_diff0)✓
print(M.T)
✓print(M.T @ M)
✓print(M @ M.T)
✓def diff0t(a): if a.ndim == 1: a = a[:,np.newaxis] # Turn 1D into 2D array d = np.zeros((a.shape[0] + 1, a.shape[1]), dtype=a.dtype) d[0, :] = - a[0, :] d[1:-1, :] = a[0:-1, :] - a[1:, :] d[-1, :] = a[-1, :] return d✓
np.allclose(M.T, diff0t(np.eye(n-1)))
✓def diff0_func_aslo_def(n): return LinearOperator(matvec=diff0, matmat=diff0, rmatvec=diff0t, rmatmat=diff0t, shape=(n - 1, n)) diff0_func_aslo = diff0_func_aslo_def(n) diff0_matrix_aslo = aslinearoperator(M_from_diff0)✓
np.allclose(diff0_func_aslo(np.eye(n)), diff0_matrix_aslo(np.eye(n))) np.allclose(diff0_func_aslo.T(np.eye(n-1)), diff0_matrix_aslo.T(np.eye(n-1)))✓
n = 100 diff0_func_aslo = diff0_func_aslo_def(n) u, s, vT = svds(diff0_func_aslo, k=3, which='SM')✓
se = 2. * np.sin(np.pi * np.arange(1, 4) / (2. * n)) ue = np.sqrt(2 / n) * np.sin(np.pi * np.outer(np.arange(1, n), np.arange(1, 4)) / n) np.allclose(s, se, atol=1e-3) np.allclose(np.abs(u), np.abs(ue), atol=1e-6)✓
Aliases
-
scipy.sparse.linalg.svds