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
andVh
, and a 1-D arrays
of singular values (real, non-negative) such thata == U @ S @ Vh
, whereS
is a suitably shaped matrix of zeros with main diagonals
.- Parameters
- a(M, N) array_like
Matrix to decompose.
- full_matricesbool, optional
If True (default), U and Vh are of shape
(M, M)
,(N, N)
. If False, the shapes are(M, K)
and(K, N)
, whereK = min(M, N)
.- compute_uvbool, optional
Whether to compute also
U
andVh
in addition tos
. Default is True.- overwrite_abool, optional
Whether to overwrite a; may improve performance. Default is False.
- check_finitebool, 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
- Undarray
Unitary matrix having left singular vectors as columns. Of shape
(M, M)
or(M, K)
, depending on full_matrices.- sndarray
The singular values, sorted in non-increasing order. Of shape (K,), with
K = min(M, N)
.- Vhndarray
Unitary matrix having right singular vectors as rows. Of shape
(N, N)
or(K, N)
depending on full_matrices.- For
compute_uv=False
, onlys
is returned.
- Raises
- LinAlgError
If SVD computation does not converge.
See also
Examples
>>> from scipy import linalg >>> from numpy.random import default_rng >>> rng = default_rng() >>> m, n = 9, 6 >>> a = rng.standard_normal((m, n)) + 1.j*rng.standard_normal((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 ofU
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