_decomp_svd.py 17 KB

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