_decomp_svd.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503
  1. """SVD decomposition functions."""
  2. import numpy
  3. from numpy import zeros, r_, diag, dot, arccos, arcsin, where, clip
  4. # Local imports.
  5. from ._misc import LinAlgError, _datacopied
  6. from .lapack import get_lapack_funcs, _compute_lwork
  7. from ._decomp import _asarray_validated
  8. __all__ = ['svd', 'svdvals', 'diagsvd', 'orth', 'subspace_angles', 'null_space']
  9. def svd(a, full_matrices=True, compute_uv=True, overwrite_a=False,
  10. check_finite=True, lapack_driver='gesdd'):
  11. """
  12. Singular Value Decomposition.
  13. Factorizes the matrix `a` into two unitary matrices ``U`` and ``Vh``, and
  14. a 1-D array ``s`` of singular values (real, non-negative) such that
  15. ``a == U @ S @ Vh``, where ``S`` is a suitably shaped matrix of zeros with
  16. main diagonal ``s``.
  17. Parameters
  18. ----------
  19. a : (M, N) array_like
  20. Matrix to decompose.
  21. full_matrices : bool, optional
  22. If True (default), `U` and `Vh` are of shape ``(M, M)``, ``(N, N)``.
  23. If False, the shapes are ``(M, K)`` and ``(K, N)``, where
  24. ``K = min(M, N)``.
  25. compute_uv : bool, optional
  26. Whether to compute also ``U`` and ``Vh`` in addition to ``s``.
  27. Default is True.
  28. overwrite_a : bool, optional
  29. Whether to overwrite `a`; may improve performance.
  30. Default is False.
  31. check_finite : bool, optional
  32. Whether to check that the input matrix contains only finite numbers.
  33. Disabling may give a performance gain, but may result in problems
  34. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  35. lapack_driver : {'gesdd', 'gesvd'}, optional
  36. Whether to use the more efficient divide-and-conquer approach
  37. (``'gesdd'``) or general rectangular approach (``'gesvd'``)
  38. to compute the SVD. MATLAB and Octave use the ``'gesvd'`` approach.
  39. Default is ``'gesdd'``.
  40. .. versionadded:: 0.18
  41. Returns
  42. -------
  43. U : ndarray
  44. Unitary matrix having left singular vectors as columns.
  45. Of shape ``(M, M)`` or ``(M, K)``, depending on `full_matrices`.
  46. s : ndarray
  47. The singular values, sorted in non-increasing order.
  48. Of shape (K,), with ``K = min(M, N)``.
  49. Vh : ndarray
  50. Unitary matrix having right singular vectors as rows.
  51. Of shape ``(N, N)`` or ``(K, N)`` depending on `full_matrices`.
  52. For ``compute_uv=False``, only ``s`` is returned.
  53. Raises
  54. ------
  55. LinAlgError
  56. If SVD computation does not converge.
  57. See Also
  58. --------
  59. svdvals : Compute singular values of a matrix.
  60. diagsvd : Construct the Sigma matrix, given the vector s.
  61. Examples
  62. --------
  63. >>> import numpy as np
  64. >>> from scipy import linalg
  65. >>> rng = np.random.default_rng()
  66. >>> m, n = 9, 6
  67. >>> a = rng.standard_normal((m, n)) + 1.j*rng.standard_normal((m, n))
  68. >>> U, s, Vh = linalg.svd(a)
  69. >>> U.shape, s.shape, Vh.shape
  70. ((9, 9), (6,), (6, 6))
  71. Reconstruct the original matrix from the decomposition:
  72. >>> sigma = np.zeros((m, n))
  73. >>> for i in range(min(m, n)):
  74. ... sigma[i, i] = s[i]
  75. >>> a1 = np.dot(U, np.dot(sigma, Vh))
  76. >>> np.allclose(a, a1)
  77. True
  78. Alternatively, use ``full_matrices=False`` (notice that the shape of
  79. ``U`` is then ``(m, n)`` instead of ``(m, m)``):
  80. >>> U, s, Vh = linalg.svd(a, full_matrices=False)
  81. >>> U.shape, s.shape, Vh.shape
  82. ((9, 6), (6,), (6, 6))
  83. >>> S = np.diag(s)
  84. >>> np.allclose(a, np.dot(U, np.dot(S, Vh)))
  85. True
  86. >>> s2 = linalg.svd(a, compute_uv=False)
  87. >>> np.allclose(s, s2)
  88. True
  89. """
  90. a1 = _asarray_validated(a, check_finite=check_finite)
  91. if len(a1.shape) != 2:
  92. raise ValueError('expected matrix')
  93. m, n = a1.shape
  94. overwrite_a = overwrite_a or (_datacopied(a1, a))
  95. if not isinstance(lapack_driver, str):
  96. raise TypeError('lapack_driver must be a string')
  97. if lapack_driver not in ('gesdd', 'gesvd'):
  98. raise ValueError('lapack_driver must be "gesdd" or "gesvd", not "%s"'
  99. % (lapack_driver,))
  100. funcs = (lapack_driver, lapack_driver + '_lwork')
  101. gesXd, gesXd_lwork = get_lapack_funcs(funcs, (a1,), ilp64='preferred')
  102. # compute optimal lwork
  103. lwork = _compute_lwork(gesXd_lwork, a1.shape[0], a1.shape[1],
  104. compute_uv=compute_uv, full_matrices=full_matrices)
  105. # perform decomposition
  106. u, s, v, info = gesXd(a1, compute_uv=compute_uv, lwork=lwork,
  107. full_matrices=full_matrices, overwrite_a=overwrite_a)
  108. if info > 0:
  109. raise LinAlgError("SVD did not converge")
  110. if info < 0:
  111. raise ValueError('illegal value in %dth argument of internal gesdd'
  112. % -info)
  113. if compute_uv:
  114. return u, s, v
  115. else:
  116. return s
  117. def svdvals(a, overwrite_a=False, check_finite=True):
  118. """
  119. Compute singular values of a matrix.
  120. Parameters
  121. ----------
  122. a : (M, N) array_like
  123. Matrix to decompose.
  124. overwrite_a : bool, optional
  125. Whether to overwrite `a`; may improve performance.
  126. Default is False.
  127. check_finite : bool, optional
  128. Whether to check that the input matrix contains only finite numbers.
  129. Disabling may give a performance gain, but may result in problems
  130. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  131. Returns
  132. -------
  133. s : (min(M, N),) ndarray
  134. The singular values, sorted in decreasing order.
  135. Raises
  136. ------
  137. LinAlgError
  138. If SVD computation does not converge.
  139. See Also
  140. --------
  141. svd : Compute the full singular value decomposition of a matrix.
  142. diagsvd : Construct the Sigma matrix, given the vector s.
  143. Notes
  144. -----
  145. ``svdvals(a)`` only differs from ``svd(a, compute_uv=False)`` by its
  146. handling of the edge case of empty ``a``, where it returns an
  147. empty sequence:
  148. >>> import numpy as np
  149. >>> a = np.empty((0, 2))
  150. >>> from scipy.linalg import svdvals
  151. >>> svdvals(a)
  152. array([], dtype=float64)
  153. Examples
  154. --------
  155. >>> import numpy as np
  156. >>> from scipy.linalg import svdvals
  157. >>> m = np.array([[1.0, 0.0],
  158. ... [2.0, 3.0],
  159. ... [1.0, 1.0],
  160. ... [0.0, 2.0],
  161. ... [1.0, 0.0]])
  162. >>> svdvals(m)
  163. array([ 4.28091555, 1.63516424])
  164. We can verify the maximum singular value of `m` by computing the maximum
  165. length of `m.dot(u)` over all the unit vectors `u` in the (x,y) plane.
  166. We approximate "all" the unit vectors with a large sample. Because
  167. of linearity, we only need the unit vectors with angles in [0, pi].
  168. >>> t = np.linspace(0, np.pi, 2000)
  169. >>> u = np.array([np.cos(t), np.sin(t)])
  170. >>> np.linalg.norm(m.dot(u), axis=0).max()
  171. 4.2809152422538475
  172. `p` is a projection matrix with rank 1. With exact arithmetic,
  173. its singular values would be [1, 0, 0, 0].
  174. >>> v = np.array([0.1, 0.3, 0.9, 0.3])
  175. >>> p = np.outer(v, v)
  176. >>> svdvals(p)
  177. array([ 1.00000000e+00, 2.02021698e-17, 1.56692500e-17,
  178. 8.15115104e-34])
  179. The singular values of an orthogonal matrix are all 1. Here, we
  180. create a random orthogonal matrix by using the `rvs()` method of
  181. `scipy.stats.ortho_group`.
  182. >>> from scipy.stats import ortho_group
  183. >>> orth = ortho_group.rvs(4)
  184. >>> svdvals(orth)
  185. array([ 1., 1., 1., 1.])
  186. """
  187. a = _asarray_validated(a, check_finite=check_finite)
  188. if a.size:
  189. return svd(a, compute_uv=0, overwrite_a=overwrite_a,
  190. check_finite=False)
  191. elif len(a.shape) != 2:
  192. raise ValueError('expected matrix')
  193. else:
  194. return numpy.empty(0)
  195. def diagsvd(s, M, N):
  196. """
  197. Construct the sigma matrix in SVD from singular values and size M, N.
  198. Parameters
  199. ----------
  200. s : (M,) or (N,) array_like
  201. Singular values
  202. M : int
  203. Size of the matrix whose singular values are `s`.
  204. N : int
  205. Size of the matrix whose singular values are `s`.
  206. Returns
  207. -------
  208. S : (M, N) ndarray
  209. The S-matrix in the singular value decomposition
  210. See Also
  211. --------
  212. svd : Singular value decomposition of a matrix
  213. svdvals : Compute singular values of a matrix.
  214. Examples
  215. --------
  216. >>> import numpy as np
  217. >>> from scipy.linalg import diagsvd
  218. >>> vals = np.array([1, 2, 3]) # The array representing the computed svd
  219. >>> diagsvd(vals, 3, 4)
  220. array([[1, 0, 0, 0],
  221. [0, 2, 0, 0],
  222. [0, 0, 3, 0]])
  223. >>> diagsvd(vals, 4, 3)
  224. array([[1, 0, 0],
  225. [0, 2, 0],
  226. [0, 0, 3],
  227. [0, 0, 0]])
  228. """
  229. part = diag(s)
  230. typ = part.dtype.char
  231. MorN = len(s)
  232. if MorN == M:
  233. return r_['-1', part, zeros((M, N-M), typ)]
  234. elif MorN == N:
  235. return r_[part, zeros((M-N, N), typ)]
  236. else:
  237. raise ValueError("Length of s must be M or N.")
  238. # Orthonormal decomposition
  239. def orth(A, rcond=None):
  240. """
  241. Construct an orthonormal basis for the range of A using SVD
  242. Parameters
  243. ----------
  244. A : (M, N) array_like
  245. Input array
  246. rcond : float, optional
  247. Relative condition number. Singular values ``s`` smaller than
  248. ``rcond * max(s)`` are considered zero.
  249. Default: floating point eps * max(M,N).
  250. Returns
  251. -------
  252. Q : (M, K) ndarray
  253. Orthonormal basis for the range of A.
  254. K = effective rank of A, as determined by rcond
  255. See Also
  256. --------
  257. svd : Singular value decomposition of a matrix
  258. null_space : Matrix null space
  259. Examples
  260. --------
  261. >>> import numpy as np
  262. >>> from scipy.linalg import orth
  263. >>> A = np.array([[2, 0, 0], [0, 5, 0]]) # rank 2 array
  264. >>> orth(A)
  265. array([[0., 1.],
  266. [1., 0.]])
  267. >>> orth(A.T)
  268. array([[0., 1.],
  269. [1., 0.],
  270. [0., 0.]])
  271. """
  272. u, s, vh = svd(A, full_matrices=False)
  273. M, N = u.shape[0], vh.shape[1]
  274. if rcond is None:
  275. rcond = numpy.finfo(s.dtype).eps * max(M, N)
  276. tol = numpy.amax(s) * rcond
  277. num = numpy.sum(s > tol, dtype=int)
  278. Q = u[:, :num]
  279. return Q
  280. def null_space(A, rcond=None):
  281. """
  282. Construct an orthonormal basis for the null space of A using SVD
  283. Parameters
  284. ----------
  285. A : (M, N) array_like
  286. Input array
  287. rcond : float, optional
  288. Relative condition number. Singular values ``s`` smaller than
  289. ``rcond * max(s)`` are considered zero.
  290. Default: floating point eps * max(M,N).
  291. Returns
  292. -------
  293. Z : (N, K) ndarray
  294. Orthonormal basis for the null space of A.
  295. K = dimension of effective null space, as determined by rcond
  296. See Also
  297. --------
  298. svd : Singular value decomposition of a matrix
  299. orth : Matrix range
  300. Examples
  301. --------
  302. 1-D null space:
  303. >>> import numpy as np
  304. >>> from scipy.linalg import null_space
  305. >>> A = np.array([[1, 1], [1, 1]])
  306. >>> ns = null_space(A)
  307. >>> ns * np.sign(ns[0,0]) # Remove the sign ambiguity of the vector
  308. array([[ 0.70710678],
  309. [-0.70710678]])
  310. 2-D null space:
  311. >>> from numpy.random import default_rng
  312. >>> rng = default_rng()
  313. >>> B = rng.random((3, 5))
  314. >>> Z = null_space(B)
  315. >>> Z.shape
  316. (5, 2)
  317. >>> np.allclose(B.dot(Z), 0)
  318. True
  319. The basis vectors are orthonormal (up to rounding error):
  320. >>> Z.T.dot(Z)
  321. array([[ 1.00000000e+00, 6.92087741e-17],
  322. [ 6.92087741e-17, 1.00000000e+00]])
  323. """
  324. u, s, vh = svd(A, full_matrices=True)
  325. M, N = u.shape[0], vh.shape[1]
  326. if rcond is None:
  327. rcond = numpy.finfo(s.dtype).eps * max(M, N)
  328. tol = numpy.amax(s) * rcond
  329. num = numpy.sum(s > tol, dtype=int)
  330. Q = vh[num:,:].T.conj()
  331. return Q
  332. def subspace_angles(A, B):
  333. r"""
  334. Compute the subspace angles between two matrices.
  335. Parameters
  336. ----------
  337. A : (M, N) array_like
  338. The first input array.
  339. B : (M, K) array_like
  340. The second input array.
  341. Returns
  342. -------
  343. angles : ndarray, shape (min(N, K),)
  344. The subspace angles between the column spaces of `A` and `B` in
  345. descending order.
  346. See Also
  347. --------
  348. orth
  349. svd
  350. Notes
  351. -----
  352. This computes the subspace angles according to the formula
  353. provided in [1]_. For equivalence with MATLAB and Octave behavior,
  354. use ``angles[0]``.
  355. .. versionadded:: 1.0
  356. References
  357. ----------
  358. .. [1] Knyazev A, Argentati M (2002) Principal Angles between Subspaces
  359. in an A-Based Scalar Product: Algorithms and Perturbation
  360. Estimates. SIAM J. Sci. Comput. 23:2008-2040.
  361. Examples
  362. --------
  363. An Hadamard matrix, which has orthogonal columns, so we expect that
  364. the suspace angle to be :math:`\frac{\pi}{2}`:
  365. >>> import numpy as np
  366. >>> from scipy.linalg import hadamard, subspace_angles
  367. >>> rng = np.random.default_rng()
  368. >>> H = hadamard(4)
  369. >>> print(H)
  370. [[ 1 1 1 1]
  371. [ 1 -1 1 -1]
  372. [ 1 1 -1 -1]
  373. [ 1 -1 -1 1]]
  374. >>> np.rad2deg(subspace_angles(H[:, :2], H[:, 2:]))
  375. array([ 90., 90.])
  376. And the subspace angle of a matrix to itself should be zero:
  377. >>> subspace_angles(H[:, :2], H[:, :2]) <= 2 * np.finfo(float).eps
  378. array([ True, True], dtype=bool)
  379. The angles between non-orthogonal subspaces are in between these extremes:
  380. >>> x = rng.standard_normal((4, 3))
  381. >>> np.rad2deg(subspace_angles(x[:, :2], x[:, [2]]))
  382. array([ 55.832]) # random
  383. """
  384. # Steps here omit the U and V calculation steps from the paper
  385. # 1. Compute orthonormal bases of column-spaces
  386. A = _asarray_validated(A, check_finite=True)
  387. if len(A.shape) != 2:
  388. raise ValueError('expected 2D array, got shape %s' % (A.shape,))
  389. QA = orth(A)
  390. del A
  391. B = _asarray_validated(B, check_finite=True)
  392. if len(B.shape) != 2:
  393. raise ValueError('expected 2D array, got shape %s' % (B.shape,))
  394. if len(B) != len(QA):
  395. raise ValueError('A and B must have the same number of rows, got '
  396. '%s and %s' % (QA.shape[0], B.shape[0]))
  397. QB = orth(B)
  398. del B
  399. # 2. Compute SVD for cosine
  400. QA_H_QB = dot(QA.T.conj(), QB)
  401. sigma = svdvals(QA_H_QB)
  402. # 3. Compute matrix B
  403. if QA.shape[1] >= QB.shape[1]:
  404. B = QB - dot(QA, QA_H_QB)
  405. else:
  406. B = QA - dot(QB, QA_H_QB.T.conj())
  407. del QA, QB, QA_H_QB
  408. # 4. Compute SVD for sine
  409. mask = sigma ** 2 >= 0.5
  410. if mask.any():
  411. mu_arcsin = arcsin(clip(svdvals(B, overwrite_a=True), -1., 1.))
  412. else:
  413. mu_arcsin = 0.
  414. # 5. Compute the principal angles
  415. # with reverse ordering of sigma because smallest sigma belongs to largest
  416. # angle theta
  417. theta = where(mask, mu_arcsin, arccos(clip(sigma[::-1], -1., 1.)))
  418. return theta