1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171
|
"""SVD decomposition functions."""
import numpy
from numpy import asarray_chkfinite, zeros, r_, diag
from scipy.linalg import calc_lwork
# Local imports.
from misc import LinAlgError, _datacopied
from lapack import get_lapack_funcs
__all__ = ['svd', 'svdvals', 'diagsvd', 'orth']
def svd(a, full_matrices=True, compute_uv=True, overwrite_a=False):
"""Singular Value Decomposition.
Factorizes the matrix a into two unitary matrices U and Vh and
an 1d-array s of singular values (real, non-negative) such that
a == U S Vh if S is an suitably shaped matrix of zeros whose
main diagonal is s.
Parameters
----------
a : array, shape (M, N)
Matrix to decompose
full_matrices : boolean
If true, U, Vh are shaped (M,M), (N,N)
If false, the shapes are (M,K), (K,N) where K = min(M,N)
compute_uv : boolean
Whether to compute also U, Vh in addition to s (Default: true)
overwrite_a : boolean
Whether data in a is overwritten (may improve performance)
Returns
-------
U: array, shape (M,M) or (M,K) depending on full_matrices
s: array, shape (K,)
The singular values, sorted so that s[i] >= s[i+1]. K = min(M, N)
Vh: array, 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
Examples
--------
>>> from scipy import random, linalg, allclose, dot
>>> a = random.randn(9, 6) + 1j*random.randn(9, 6)
>>> U, s, Vh = linalg.svd(a)
>>> U.shape, Vh.shape, s.shape
((9, 9), (6, 6), (6,))
>>> U, s, Vh = linalg.svd(a, full_matrices=False)
>>> U.shape, Vh.shape, s.shape
((9, 6), (6, 6), (6,))
>>> S = linalg.diagsvd(s, 6, 6)
>>> allclose(a, dot(U, dot(S, Vh)))
True
>>> s2 = linalg.svd(a, compute_uv=False)
>>> allclose(s, s2)
True
See also
--------
svdvals : return singular values of a matrix
diagsvd : return the Sigma matrix, given the vector s
"""
a1 = asarray_chkfinite(a)
if len(a1.shape) != 2:
raise ValueError('expected matrix')
m,n = a1.shape
overwrite_a = overwrite_a or (_datacopied(a1, a))
gesdd, = get_lapack_funcs(('gesdd',), (a1,))
if gesdd.module_name[:7] == 'flapack':
lwork = calc_lwork.gesdd(gesdd.prefix, m, n, compute_uv)[1]
u,s,v,info = gesdd(a1,compute_uv = compute_uv, lwork = lwork,
full_matrices=full_matrices, overwrite_a = overwrite_a)
else: # 'clapack'
raise NotImplementedError('calling gesdd from %s' % gesdd.module_name)
if info > 0:
raise LinAlgError("SVD did not converge")
if info < 0:
raise ValueError('illegal value in %d-th argument of internal gesdd'
% -info)
if compute_uv:
return u, s, v
else:
return s
def svdvals(a, overwrite_a=False):
"""Compute singular values of a matrix.
Parameters
----------
a : array, shape (M, N)
Matrix to decompose
overwrite_a : boolean
Whether data in a is overwritten (may improve performance)
Returns
-------
s: array, shape (K,)
The singular values, sorted so that s[i] >= s[i+1]. K = min(M, N)
Raises LinAlgError if SVD computation does not converge
See also
--------
svd : return the full singular value decomposition of a matrix
diagsvd : return the Sigma matrix, given the vector s
"""
return svd(a, compute_uv=0, overwrite_a=overwrite_a)
def diagsvd(s, M, N):
"""Construct the sigma matrix in SVD from singular values and size M,N.
Parameters
----------
s : array, shape (M,) or (N,)
Singular values
M : integer
N : integer
Size of the matrix whose singular values are s
Returns
-------
S : array, shape (M, N)
The S-matrix in the singular value decomposition
"""
part = diag(s)
typ = part.dtype.char
MorN = len(s)
if MorN == M:
return r_['-1', part, zeros((M, N-M), typ)]
elif MorN == N:
return r_[part, zeros((M-N,N), typ)]
else:
raise ValueError("Length of s must be M or N.")
# Orthonormal decomposition
def orth(A):
"""Construct an orthonormal basis for the range of A using SVD
Parameters
----------
A : array, shape (M, N)
Returns
-------
Q : array, shape (M, K)
Orthonormal basis for the range of A.
K = effective rank of A, as determined by automatic cutoff
See also
--------
svd : Singular value decomposition of a matrix
"""
u, s, vh = svd(A)
M, N = A.shape
eps = numpy.finfo(float).eps
tol = max(M,N) * numpy.amax(s) * eps
num = numpy.sum(s > tol, dtype=int)
Q = u[:,:num]
return Q
|