SciPy

scipy.sparse.linalg.svds

scipy.sparse.linalg.svds(A, k=6, ncv=None, tol=0, which='LM', v0=None, maxiter=None, return_singular_vectors=True)[source]

Compute the largest k singular values/vectors for a sparse matrix.

Parameters
A{sparse matrix, LinearOperator}

Array to compute the SVD on, of shape (M, N)

kint, optional

Number of singular values and vectors to compute. Must be 1 <= k < min(A.shape).

ncvint, optional

The number of Lanczos vectors generated ncv must be greater than k+1 and smaller than n; it is recommended that ncv > 2*k Default: min(n, max(2*k + 1, 20))

tolfloat, optional

Tolerance for singular values. Zero (default) means machine precision.

whichstr, [‘LM’ | ‘SM’], optional

Which k singular values to find:

  • ‘LM’ : largest singular values

  • ‘SM’ : smallest singular values

New in version 0.12.0.

v0ndarray, optional

Starting vector for iteration, of length min(A.shape). Should be an (approximate) left singular vector if N > M and a right singular vector otherwise. Default: random

New in version 0.12.0.

maxiterint, optional

Maximum number of iterations.

New in version 0.12.0.

return_singular_vectorsbool or str, optional
  • True: return singular vectors (True) in addition to singular values.

New in version 0.12.0.

  • “u”: only return the u matrix, without computing vh (if N > M).

  • “vh”: only return the vh matrix, without computing u (if N <= M).

New in version 0.16.0.

Returns
undarray, shape=(M, k)

Unitary matrix having left singular vectors as columns. If return_singular_vectors is “vh”, this variable is not computed, and None is returned instead.

sndarray, shape=(k,)

The singular values.

vtndarray, shape=(k, N)

Unitary matrix having right singular vectors as rows. If return_singular_vectors is “u”, this variable is not computed, and None is returned instead.

Notes

This is a naive implementation using ARPACK as an eigensolver on A.H * A or A * A.H, depending on which one is more efficient.

Examples

>>> from scipy.sparse import csc_matrix
>>> from scipy.sparse.linalg import svds, eigs
>>> A = csc_matrix([[1, 0, 0], [5, 0, 2], [0, -1, 0], [0, 0, 3]], dtype=float)
>>> u, s, vt = svds(A, k=2)
>>> s
array([ 2.75193379,  5.6059665 ])
>>> np.sqrt(eigs(A.dot(A.T), k=2)[0]).real
array([ 5.6059665 ,  2.75193379])

Previous topic

scipy.sparse.linalg.lobpcg

Next topic

scipy.sparse.linalg.splu