SciPy

scipy.linalg.svd

scipy.linalg.svd(a, full_matrices=True, compute_uv=True, overwrite_a=False, check_finite=True, lapack_driver='gesdd')[source]

Singular Value Decomposition.

Factorizes the matrix a into two unitary matrices U and Vh, and a 1-D array s of singular values (real, non-negative) such that a == U @ S @ Vh, where S is a suitably shaped matrix of zeros with main diagonal s.

Parameters:
a : (M, N) array_like

Matrix to decompose.

full_matrices : bool, optional

If True (default), U and Vh are of shape (M, M), (N, N). If False, the shapes are (M, K) and (K, N), where K = min(M, N).

compute_uv : bool, optional

Whether to compute also U and Vh in addition to s. Default is True.

overwrite_a : bool, optional

Whether to overwrite a; may improve performance. Default is False.

check_finite : bool, optional

Whether to check that the input matrix contains only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs.

lapack_driver : {‘gesdd’, ‘gesvd’}, optional

Whether to use the more efficient divide-and-conquer approach ('gesdd') or general rectangular approach ('gesvd') to compute the SVD. MATLAB and Octave use the 'gesvd' approach. Default is 'gesdd'.

New in version 0.18.

Returns:
U : ndarray

Unitary matrix having left singular vectors as columns. Of shape (M, M) or (M, K), depending on full_matrices.

s : ndarray

The singular values, sorted in non-increasing order. Of shape (K,), with K = min(M, N).

Vh : ndarray

Unitary matrix having right singular vectors as rows. Of shape (N, N) or (K, N) depending on full_matrices.

For ``compute_uv=False``, only ``s`` is returned.
Raises:
LinAlgError

If SVD computation does not converge.

See also

svdvals
Compute singular values of a matrix.
diagsvd
Construct the Sigma matrix, given the vector s.

Examples

>>> from scipy import linalg
>>> m, n = 9, 6
>>> a = np.random.randn(m, n) + 1.j*np.random.randn(m, n)
>>> U, s, Vh = linalg.svd(a)
>>> U.shape,  s.shape, Vh.shape
((9, 9), (6,), (6, 6))

Reconstruct the original matrix from the decomposition:

>>> sigma = np.zeros((m, n))
>>> for i in range(min(m, n)):
...     sigma[i, i] = s[i]
>>> a1 = np.dot(U, np.dot(sigma, Vh))
>>> np.allclose(a, a1)
True

Alternatively, use full_matrices=False (notice that the shape of U is then (m, n) instead of (m, m)):

>>> U, s, Vh = linalg.svd(a, full_matrices=False)
>>> U.shape, s.shape, Vh.shape
((9, 6), (6,), (6, 6))
>>> S = np.diag(s)
>>> np.allclose(a, np.dot(U, np.dot(S, Vh)))
True
>>> s2 = linalg.svd(a, compute_uv=False)
>>> np.allclose(s, s2)
True

Previous topic

scipy.linalg.lu_solve

Next topic

scipy.linalg.svdvals