interpolative.py 30 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941
  1. # ******************************************************************************
  2. # Copyright (C) 2013 Kenneth L. Ho
  3. #
  4. # Redistribution and use in source and binary forms, with or without
  5. # modification, are permitted provided that the following conditions are met:
  6. #
  7. # Redistributions of source code must retain the above copyright notice, this
  8. # list of conditions and the following disclaimer. Redistributions in binary
  9. # form must reproduce the above copyright notice, this list of conditions and
  10. # the following disclaimer in the documentation and/or other materials
  11. # provided with the distribution.
  12. #
  13. # None of the names of the copyright holders may be used to endorse or
  14. # promote products derived from this software without specific prior written
  15. # permission.
  16. #
  17. # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  18. # AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  19. # IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  20. # ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
  21. # LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  22. # CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  23. # SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  24. # INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  25. # CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  26. # ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  27. # POSSIBILITY OF SUCH DAMAGE.
  28. # ******************************************************************************
  29. r"""
  30. ======================================================================
  31. Interpolative matrix decomposition (:mod:`scipy.linalg.interpolative`)
  32. ======================================================================
  33. .. versionadded:: 0.13
  34. .. versionchanged:: 1.15.0
  35. The underlying algorithms have been ported to Python from the original Fortran77
  36. code. See references below for more details.
  37. .. currentmodule:: scipy.linalg.interpolative
  38. An interpolative decomposition (ID) of a matrix :math:`A \in
  39. \mathbb{C}^{m \times n}` of rank :math:`k \leq \min \{ m, n \}` is a
  40. factorization
  41. .. math::
  42. A \Pi =
  43. \begin{bmatrix}
  44. A \Pi_{1} & A \Pi_{2}
  45. \end{bmatrix} =
  46. A \Pi_{1}
  47. \begin{bmatrix}
  48. I & T
  49. \end{bmatrix},
  50. where :math:`\Pi = [\Pi_{1}, \Pi_{2}]` is a permutation matrix with
  51. :math:`\Pi_{1} \in \{ 0, 1 \}^{n \times k}`, i.e., :math:`A \Pi_{2} =
  52. A \Pi_{1} T`. This can equivalently be written as :math:`A = BP`,
  53. where :math:`B = A \Pi_{1}` and :math:`P = [I, T] \Pi^{\mathsf{T}}`
  54. are the *skeleton* and *interpolation matrices*, respectively.
  55. If :math:`A` does not have exact rank :math:`k`, then there exists an
  56. approximation in the form of an ID such that :math:`A = BP + E`, where
  57. :math:`\| E \| \sim \sigma_{k + 1}` is on the order of the :math:`(k +
  58. 1)`-th largest singular value of :math:`A`. Note that :math:`\sigma_{k
  59. + 1}` is the best possible error for a rank-:math:`k` approximation
  60. and, in fact, is achieved by the singular value decomposition (SVD)
  61. :math:`A \approx U S V^{*}`, where :math:`U \in \mathbb{C}^{m \times
  62. k}` and :math:`V \in \mathbb{C}^{n \times k}` have orthonormal columns
  63. and :math:`S = \mathop{\mathrm{diag}} (\sigma_{i}) \in \mathbb{C}^{k
  64. \times k}` is diagonal with nonnegative entries. The principal
  65. advantages of using an ID over an SVD are that:
  66. - it is cheaper to construct;
  67. - it preserves the structure of :math:`A`; and
  68. - it is more efficient to compute with in light of the identity submatrix of :math:`P`.
  69. Routines
  70. ========
  71. Main functionality:
  72. .. autosummary::
  73. :toctree: generated/
  74. interp_decomp
  75. reconstruct_matrix_from_id
  76. reconstruct_interp_matrix
  77. reconstruct_skel_matrix
  78. id_to_svd
  79. svd
  80. estimate_spectral_norm
  81. estimate_spectral_norm_diff
  82. estimate_rank
  83. References
  84. ==========
  85. This module uses the algorithms found in ID software package [1]_ by Martinsson,
  86. Rokhlin, Shkolnisky, and Tygert, which is a Fortran library for computing IDs using
  87. various algorithms, including the rank-revealing QR approach of [2]_ and the more
  88. recent randomized methods described in [3]_, [4]_, and [5]_.
  89. We advise the user to consult also the documentation for the `ID package
  90. <http://tygert.com/software.html>`_.
  91. .. [1] P.G. Martinsson, V. Rokhlin, Y. Shkolnisky, M. Tygert. "ID: a
  92. software package for low-rank approximation of matrices via interpolative
  93. decompositions, version 0.2." http://tygert.com/id_doc.4.pdf.
  94. .. [2] H. Cheng, Z. Gimbutas, P.G. Martinsson, V. Rokhlin. "On the
  95. compression of low rank matrices." *SIAM J. Sci. Comput.* 26 (4): 1389--1404,
  96. 2005. :doi:`10.1137/030602678`.
  97. .. [3] E. Liberty, F. Woolfe, P.G. Martinsson, V. Rokhlin, M.
  98. Tygert. "Randomized algorithms for the low-rank approximation of matrices."
  99. *Proc. Natl. Acad. Sci. U.S.A.* 104 (51): 20167--20172, 2007.
  100. :doi:`10.1073/pnas.0709640104`.
  101. .. [4] P.G. Martinsson, V. Rokhlin, M. Tygert. "A randomized
  102. algorithm for the decomposition of matrices." *Appl. Comput. Harmon. Anal.* 30
  103. (1): 47--68, 2011. :doi:`10.1016/j.acha.2010.02.003`.
  104. .. [5] F. Woolfe, E. Liberty, V. Rokhlin, M. Tygert. "A fast
  105. randomized algorithm for the approximation of matrices." *Appl. Comput.
  106. Harmon. Anal.* 25 (3): 335--366, 2008. :doi:`10.1016/j.acha.2007.12.002`.
  107. Tutorial
  108. ========
  109. Initializing
  110. ------------
  111. The first step is to import :mod:`scipy.linalg.interpolative` by issuing the
  112. command:
  113. >>> import scipy.linalg.interpolative as sli
  114. Now let's build a matrix. For this, we consider a Hilbert matrix, which is well
  115. know to have low rank:
  116. >>> from scipy.linalg import hilbert
  117. >>> n = 1000
  118. >>> A = hilbert(n)
  119. We can also do this explicitly via:
  120. >>> import numpy as np
  121. >>> n = 1000
  122. >>> A = np.empty((n, n), order='F')
  123. >>> for j in range(n):
  124. ... for i in range(n):
  125. ... A[i,j] = 1. / (i + j + 1)
  126. Note the use of the flag ``order='F'`` in :func:`numpy.empty`. This
  127. instantiates the matrix in Fortran-contiguous order and is important for
  128. avoiding data copying when passing to the backend.
  129. We then define multiplication routines for the matrix by regarding it as a
  130. :class:`scipy.sparse.linalg.LinearOperator`:
  131. >>> from scipy.sparse.linalg import aslinearoperator
  132. >>> L = aslinearoperator(A)
  133. This automatically sets up methods describing the action of the matrix and its
  134. adjoint on a vector.
  135. Computing an ID
  136. ---------------
  137. We have several choices of algorithm to compute an ID. These fall largely
  138. according to two dichotomies:
  139. 1. how the matrix is represented, i.e., via its entries or via its action on a
  140. vector; and
  141. 2. whether to approximate it to a fixed relative precision or to a fixed rank.
  142. We step through each choice in turn below.
  143. In all cases, the ID is represented by three parameters:
  144. 1. a rank ``k``;
  145. 2. an index array ``idx``; and
  146. 3. interpolation coefficients ``proj``.
  147. The ID is specified by the relation
  148. ``np.dot(A[:,idx[:k]], proj) == A[:,idx[k:]]``.
  149. From matrix entries
  150. ...................
  151. We first consider a matrix given in terms of its entries.
  152. To compute an ID to a fixed precision, type:
  153. >>> eps = 1e-3
  154. >>> k, idx, proj = sli.interp_decomp(A, eps)
  155. where ``eps < 1`` is the desired precision.
  156. To compute an ID to a fixed rank, use:
  157. >>> idx, proj = sli.interp_decomp(A, k)
  158. where ``k >= 1`` is the desired rank.
  159. Both algorithms use random sampling and are usually faster than the
  160. corresponding older, deterministic algorithms, which can be accessed via the
  161. commands:
  162. >>> k, idx, proj = sli.interp_decomp(A, eps, rand=False)
  163. and:
  164. >>> idx, proj = sli.interp_decomp(A, k, rand=False)
  165. respectively.
  166. From matrix action
  167. ..................
  168. Now consider a matrix given in terms of its action on a vector as a
  169. :class:`scipy.sparse.linalg.LinearOperator`.
  170. To compute an ID to a fixed precision, type:
  171. >>> k, idx, proj = sli.interp_decomp(L, eps)
  172. To compute an ID to a fixed rank, use:
  173. >>> idx, proj = sli.interp_decomp(L, k)
  174. These algorithms are randomized.
  175. Reconstructing an ID
  176. --------------------
  177. The ID routines above do not output the skeleton and interpolation matrices
  178. explicitly but instead return the relevant information in a more compact (and
  179. sometimes more useful) form. To build these matrices, write:
  180. >>> B = sli.reconstruct_skel_matrix(A, k, idx)
  181. for the skeleton matrix and:
  182. >>> P = sli.reconstruct_interp_matrix(idx, proj)
  183. for the interpolation matrix. The ID approximation can then be computed as:
  184. >>> C = np.dot(B, P)
  185. This can also be constructed directly using:
  186. >>> C = sli.reconstruct_matrix_from_id(B, idx, proj)
  187. without having to first compute ``P``.
  188. Alternatively, this can be done explicitly as well using:
  189. >>> B = A[:,idx[:k]]
  190. >>> P = np.hstack([np.eye(k), proj])[:,np.argsort(idx)]
  191. >>> C = np.dot(B, P)
  192. Computing an SVD
  193. ----------------
  194. An ID can be converted to an SVD via the command:
  195. >>> U, S, V = sli.id_to_svd(B, idx, proj)
  196. The SVD approximation is then:
  197. >>> approx = U @ np.diag(S) @ V.conj().T
  198. The SVD can also be computed "fresh" by combining both the ID and conversion
  199. steps into one command. Following the various ID algorithms above, there are
  200. correspondingly various SVD algorithms that one can employ.
  201. From matrix entries
  202. ...................
  203. We consider first SVD algorithms for a matrix given in terms of its entries.
  204. To compute an SVD to a fixed precision, type:
  205. >>> U, S, V = sli.svd(A, eps)
  206. To compute an SVD to a fixed rank, use:
  207. >>> U, S, V = sli.svd(A, k)
  208. Both algorithms use random sampling; for the deterministic versions, issue the
  209. keyword ``rand=False`` as above.
  210. From matrix action
  211. ..................
  212. Now consider a matrix given in terms of its action on a vector.
  213. To compute an SVD to a fixed precision, type:
  214. >>> U, S, V = sli.svd(L, eps)
  215. To compute an SVD to a fixed rank, use:
  216. >>> U, S, V = sli.svd(L, k)
  217. Utility routines
  218. ----------------
  219. Several utility routines are also available.
  220. To estimate the spectral norm of a matrix, use:
  221. >>> snorm = sli.estimate_spectral_norm(A)
  222. This algorithm is based on the randomized power method and thus requires only
  223. matrix-vector products. The number of iterations to take can be set using the
  224. keyword ``its`` (default: ``its=20``). The matrix is interpreted as a
  225. :class:`scipy.sparse.linalg.LinearOperator`, but it is also valid to supply it
  226. as a :class:`numpy.ndarray`, in which case it is trivially converted using
  227. :func:`scipy.sparse.linalg.aslinearoperator`.
  228. The same algorithm can also estimate the spectral norm of the difference of two
  229. matrices ``A1`` and ``A2`` as follows:
  230. >>> A1, A2 = A**2, A
  231. >>> diff = sli.estimate_spectral_norm_diff(A1, A2)
  232. This is often useful for checking the accuracy of a matrix approximation.
  233. Some routines in :mod:`scipy.linalg.interpolative` require estimating the rank
  234. of a matrix as well. This can be done with either:
  235. >>> k = sli.estimate_rank(A, eps)
  236. or:
  237. >>> k = sli.estimate_rank(L, eps)
  238. depending on the representation. The parameter ``eps`` controls the definition
  239. of the numerical rank.
  240. Finally, the random number generation required for all randomized routines can
  241. be controlled via providing NumPy pseudo-random generators with a fixed seed. See
  242. :class:`numpy.random.Generator` and :func:`numpy.random.default_rng` for more details.
  243. Remarks
  244. -------
  245. The above functions all automatically detect the appropriate interface and work
  246. with both real and complex data types, passing input arguments to the proper
  247. backend routine.
  248. """
  249. import scipy.linalg._decomp_interpolative as _backend
  250. import numpy as np
  251. __all__ = [
  252. 'estimate_rank',
  253. 'estimate_spectral_norm',
  254. 'estimate_spectral_norm_diff',
  255. 'id_to_svd',
  256. 'interp_decomp',
  257. 'reconstruct_interp_matrix',
  258. 'reconstruct_matrix_from_id',
  259. 'reconstruct_skel_matrix',
  260. 'svd',
  261. ]
  262. _DTYPE_ERROR = ValueError("invalid input dtype (input must be float64 or complex128)")
  263. _TYPE_ERROR = TypeError("invalid input type (must be array or LinearOperator)")
  264. def _C_contiguous_copy(A):
  265. """
  266. Same as np.ascontiguousarray, but ensure a copy
  267. """
  268. A = np.asarray(A)
  269. if A.flags.c_contiguous:
  270. A = A.copy()
  271. else:
  272. A = np.ascontiguousarray(A)
  273. return A
  274. def _is_real(A):
  275. try:
  276. if A.dtype == np.complex128:
  277. return False
  278. elif A.dtype == np.float64:
  279. return True
  280. else:
  281. raise _DTYPE_ERROR
  282. except AttributeError as e:
  283. raise _TYPE_ERROR from e
  284. def interp_decomp(A, eps_or_k, rand=True, rng=None):
  285. """
  286. Compute ID of a matrix.
  287. An ID of a matrix `A` is a factorization defined by a rank `k`, a column
  288. index array `idx`, and interpolation coefficients `proj` such that::
  289. numpy.dot(A[:,idx[:k]], proj) = A[:,idx[k:]]
  290. The original matrix can then be reconstructed as::
  291. numpy.hstack([A[:,idx[:k]],
  292. numpy.dot(A[:,idx[:k]], proj)]
  293. )[:,numpy.argsort(idx)]
  294. or via the routine :func:`reconstruct_matrix_from_id`. This can
  295. equivalently be written as::
  296. numpy.dot(A[:,idx[:k]],
  297. numpy.hstack([numpy.eye(k), proj])
  298. )[:,np.argsort(idx)]
  299. in terms of the skeleton and interpolation matrices::
  300. B = A[:,idx[:k]]
  301. and::
  302. P = numpy.hstack([numpy.eye(k), proj])[:,np.argsort(idx)]
  303. respectively. See also :func:`reconstruct_interp_matrix` and
  304. :func:`reconstruct_skel_matrix`.
  305. The ID can be computed to any relative precision or rank (depending on the
  306. value of `eps_or_k`). If a precision is specified (`eps_or_k < 1`), then
  307. this function has the output signature::
  308. k, idx, proj = interp_decomp(A, eps_or_k)
  309. Otherwise, if a rank is specified (`eps_or_k >= 1`), then the output
  310. signature is::
  311. idx, proj = interp_decomp(A, eps_or_k)
  312. .. This function automatically detects the form of the input parameters
  313. and passes them to the appropriate backend. For details, see
  314. :func:`_backend.iddp_id`, :func:`_backend.iddp_aid`,
  315. :func:`_backend.iddp_rid`, :func:`_backend.iddr_id`,
  316. :func:`_backend.iddr_aid`, :func:`_backend.iddr_rid`,
  317. :func:`_backend.idzp_id`, :func:`_backend.idzp_aid`,
  318. :func:`_backend.idzp_rid`, :func:`_backend.idzr_id`,
  319. :func:`_backend.idzr_aid`, and :func:`_backend.idzr_rid`.
  320. Parameters
  321. ----------
  322. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator` with `rmatvec`
  323. Matrix to be factored
  324. eps_or_k : float or int
  325. Relative error (if ``eps_or_k < 1``) or rank (if ``eps_or_k >= 1``) of
  326. approximation.
  327. rand : bool, optional
  328. Whether to use random sampling if `A` is of type :class:`numpy.ndarray`
  329. (randomized algorithms are always used if `A` is of type
  330. :class:`scipy.sparse.linalg.LinearOperator`).
  331. rng : `numpy.random.Generator`, optional
  332. Pseudorandom number generator state. When `rng` is None, a new
  333. `numpy.random.Generator` is created using entropy from the
  334. operating system. Types other than `numpy.random.Generator` are
  335. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  336. If `rand` is ``False``, the argument is ignored.
  337. Returns
  338. -------
  339. k : int
  340. Rank required to achieve specified relative precision if
  341. ``eps_or_k < 1``.
  342. idx : :class:`numpy.ndarray`
  343. Column index array.
  344. proj : :class:`numpy.ndarray`
  345. Interpolation coefficients.
  346. """ # numpy/numpydoc#87 # noqa: E501
  347. from scipy.sparse.linalg import LinearOperator
  348. rng = np.random.default_rng(rng)
  349. real = _is_real(A)
  350. if isinstance(A, np.ndarray):
  351. A = _C_contiguous_copy(A)
  352. if eps_or_k < 1:
  353. eps = eps_or_k
  354. if rand:
  355. if real:
  356. k, idx, proj = _backend.iddp_aid(A, eps, rng=rng)
  357. else:
  358. k, idx, proj = _backend.idzp_aid(A, eps, rng=rng)
  359. else:
  360. if real:
  361. k, idx, proj = _backend.iddp_id(A, eps)
  362. else:
  363. k, idx, proj = _backend.idzp_id(A, eps)
  364. return k, idx, proj
  365. else:
  366. k = int(eps_or_k)
  367. if rand:
  368. if real:
  369. idx, proj = _backend.iddr_aid(A, k, rng=rng)
  370. else:
  371. idx, proj = _backend.idzr_aid(A, k, rng=rng)
  372. else:
  373. if real:
  374. idx, proj = _backend.iddr_id(A, k)
  375. else:
  376. idx, proj = _backend.idzr_id(A, k)
  377. return idx, proj
  378. elif isinstance(A, LinearOperator):
  379. if eps_or_k < 1:
  380. eps = eps_or_k
  381. if real:
  382. k, idx, proj = _backend.iddp_rid(A, eps, rng=rng)
  383. else:
  384. k, idx, proj = _backend.idzp_rid(A, eps, rng=rng)
  385. return k, idx, proj
  386. else:
  387. k = int(eps_or_k)
  388. if real:
  389. idx, proj = _backend.iddr_rid(A, k, rng=rng)
  390. else:
  391. idx, proj = _backend.idzr_rid(A, k, rng=rng)
  392. return idx, proj
  393. else:
  394. raise _TYPE_ERROR
  395. def reconstruct_matrix_from_id(B, idx, proj):
  396. """
  397. Reconstruct matrix from its ID.
  398. A matrix `A` with skeleton matrix `B` and ID indices and coefficients `idx`
  399. and `proj`, respectively, can be reconstructed as::
  400. numpy.hstack([B, numpy.dot(B, proj)])[:,numpy.argsort(idx)]
  401. See also :func:`reconstruct_interp_matrix` and
  402. :func:`reconstruct_skel_matrix`.
  403. .. This function automatically detects the matrix data type and calls the
  404. appropriate backend. For details, see :func:`_backend.idd_reconid` and
  405. :func:`_backend.idz_reconid`.
  406. Parameters
  407. ----------
  408. B : :class:`numpy.ndarray`
  409. Skeleton matrix.
  410. idx : :class:`numpy.ndarray`
  411. Column index array.
  412. proj : :class:`numpy.ndarray`
  413. Interpolation coefficients.
  414. Returns
  415. -------
  416. :class:`numpy.ndarray`
  417. Reconstructed matrix.
  418. """
  419. B = np.atleast_2d(_C_contiguous_copy(B))
  420. if _is_real(B):
  421. return _backend.idd_reconid(B, idx, proj)
  422. else:
  423. return _backend.idz_reconid(B, idx, proj)
  424. def reconstruct_interp_matrix(idx, proj):
  425. """
  426. Reconstruct interpolation matrix from ID.
  427. The interpolation matrix can be reconstructed from the ID indices and
  428. coefficients `idx` and `proj`, respectively, as::
  429. P = numpy.hstack([numpy.eye(proj.shape[0]), proj])[:,numpy.argsort(idx)]
  430. The original matrix can then be reconstructed from its skeleton matrix ``B``
  431. via ``A = B @ P``
  432. See also :func:`reconstruct_matrix_from_id` and
  433. :func:`reconstruct_skel_matrix`.
  434. .. This function automatically detects the matrix data type and calls the
  435. appropriate backend. For details, see :func:`_backend.idd_reconint` and
  436. :func:`_backend.idz_reconint`.
  437. Parameters
  438. ----------
  439. idx : :class:`numpy.ndarray`
  440. 1D column index array.
  441. proj : :class:`numpy.ndarray`
  442. Interpolation coefficients.
  443. Returns
  444. -------
  445. :class:`numpy.ndarray`
  446. Interpolation matrix.
  447. """
  448. n, krank = len(idx), proj.shape[0]
  449. if _is_real(proj):
  450. p = np.zeros([krank, n], dtype=np.float64)
  451. else:
  452. p = np.zeros([krank, n], dtype=np.complex128)
  453. for ci in range(krank):
  454. p[ci, idx[ci]] = 1.0
  455. p[:, idx[krank:]] = proj[:, :]
  456. return p
  457. def reconstruct_skel_matrix(A, k, idx):
  458. """
  459. Reconstruct skeleton matrix from ID.
  460. The skeleton matrix can be reconstructed from the original matrix `A` and its
  461. ID rank and indices `k` and `idx`, respectively, as::
  462. B = A[:,idx[:k]]
  463. The original matrix can then be reconstructed via::
  464. numpy.hstack([B, numpy.dot(B, proj)])[:,numpy.argsort(idx)]
  465. See also :func:`reconstruct_matrix_from_id` and
  466. :func:`reconstruct_interp_matrix`.
  467. .. This function automatically detects the matrix data type and calls the
  468. appropriate backend. For details, see :func:`_backend.idd_copycols` and
  469. :func:`_backend.idz_copycols`.
  470. Parameters
  471. ----------
  472. A : :class:`numpy.ndarray`
  473. Original matrix.
  474. k : int
  475. Rank of ID.
  476. idx : :class:`numpy.ndarray`
  477. Column index array.
  478. Returns
  479. -------
  480. :class:`numpy.ndarray`
  481. Skeleton matrix.
  482. """
  483. return A[:, idx[:k]]
  484. def id_to_svd(B, idx, proj):
  485. """
  486. Convert ID to SVD.
  487. The SVD reconstruction of a matrix with skeleton matrix `B` and ID indices and
  488. coefficients `idx` and `proj`, respectively, is::
  489. U, S, V = id_to_svd(B, idx, proj)
  490. A = numpy.dot(U, numpy.dot(numpy.diag(S), V.conj().T))
  491. See also :func:`svd`.
  492. .. This function automatically detects the matrix data type and calls the
  493. appropriate backend. For details, see :func:`_backend.idd_id2svd` and
  494. :func:`_backend.idz_id2svd`.
  495. Parameters
  496. ----------
  497. B : :class:`numpy.ndarray`
  498. Skeleton matrix.
  499. idx : :class:`numpy.ndarray`
  500. 1D column index array.
  501. proj : :class:`numpy.ndarray`
  502. Interpolation coefficients.
  503. Returns
  504. -------
  505. U : :class:`numpy.ndarray`
  506. Left singular vectors.
  507. S : :class:`numpy.ndarray`
  508. Singular values.
  509. V : :class:`numpy.ndarray`
  510. Right singular vectors.
  511. """
  512. B = _C_contiguous_copy(B)
  513. if _is_real(B):
  514. U, S, V = _backend.idd_id2svd(B, idx, proj)
  515. else:
  516. U, S, V = _backend.idz_id2svd(B, idx, proj)
  517. return U, S, V
  518. def estimate_spectral_norm(A, its=20, rng=None):
  519. """
  520. Estimate spectral norm of a matrix by the randomized power method.
  521. .. This function automatically detects the matrix data type and calls the
  522. appropriate backend. For details, see :func:`_backend.idd_snorm` and
  523. :func:`_backend.idz_snorm`.
  524. Parameters
  525. ----------
  526. A : :class:`scipy.sparse.linalg.LinearOperator`
  527. Matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with the
  528. `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  529. its : int, optional
  530. Number of power method iterations.
  531. rng : `numpy.random.Generator`, optional
  532. Pseudorandom number generator state. When `rng` is None, a new
  533. `numpy.random.Generator` is created using entropy from the
  534. operating system. Types other than `numpy.random.Generator` are
  535. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  536. If `rand` is ``False``, the argument is ignored.
  537. Returns
  538. -------
  539. float
  540. Spectral norm estimate.
  541. """
  542. from scipy.sparse.linalg import aslinearoperator
  543. rng = np.random.default_rng(rng)
  544. A = aslinearoperator(A)
  545. if _is_real(A):
  546. return _backend.idd_snorm(A, its=its, rng=rng)
  547. else:
  548. return _backend.idz_snorm(A, its=its, rng=rng)
  549. def estimate_spectral_norm_diff(A, B, its=20, rng=None):
  550. """
  551. Estimate spectral norm of the difference of two matrices by the randomized
  552. power method.
  553. .. This function automatically detects the matrix data type and calls the
  554. appropriate backend. For details, see :func:`_backend.idd_diffsnorm` and
  555. :func:`_backend.idz_diffsnorm`.
  556. Parameters
  557. ----------
  558. A : :class:`scipy.sparse.linalg.LinearOperator`
  559. First matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with the
  560. `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  561. B : :class:`scipy.sparse.linalg.LinearOperator`
  562. Second matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with
  563. the `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  564. its : int, optional
  565. Number of power method iterations.
  566. rng : `numpy.random.Generator`, optional
  567. Pseudorandom number generator state. When `rng` is None, a new
  568. `numpy.random.Generator` is created using entropy from the
  569. operating system. Types other than `numpy.random.Generator` are
  570. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  571. If `rand` is ``False``, the argument is ignored.
  572. Returns
  573. -------
  574. float
  575. Spectral norm estimate of matrix difference.
  576. """
  577. from scipy.sparse.linalg import aslinearoperator
  578. rng = np.random.default_rng(rng)
  579. A = aslinearoperator(A)
  580. B = aslinearoperator(B)
  581. if _is_real(A):
  582. return _backend.idd_diffsnorm(A, B, its=its, rng=rng)
  583. else:
  584. return _backend.idz_diffsnorm(A, B, its=its, rng=rng)
  585. def svd(A, eps_or_k, rand=True, rng=None):
  586. """
  587. Compute SVD of a matrix via an ID.
  588. An SVD of a matrix `A` is a factorization::
  589. A = U @ np.diag(S) @ V.conj().T
  590. where `U` and `V` have orthonormal columns and `S` is nonnegative.
  591. The SVD can be computed to any relative precision or rank (depending on the
  592. value of `eps_or_k`).
  593. See also :func:`interp_decomp` and :func:`id_to_svd`.
  594. .. This function automatically detects the form of the input parameters and
  595. passes them to the appropriate backend. For details, see
  596. :func:`_backend.iddp_svd`, :func:`_backend.iddp_asvd`,
  597. :func:`_backend.iddp_rsvd`, :func:`_backend.iddr_svd`,
  598. :func:`_backend.iddr_asvd`, :func:`_backend.iddr_rsvd`,
  599. :func:`_backend.idzp_svd`, :func:`_backend.idzp_asvd`,
  600. :func:`_backend.idzp_rsvd`, :func:`_backend.idzr_svd`,
  601. :func:`_backend.idzr_asvd`, and :func:`_backend.idzr_rsvd`.
  602. Parameters
  603. ----------
  604. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator`
  605. Matrix to be factored, given as either a :class:`numpy.ndarray` or a
  606. :class:`scipy.sparse.linalg.LinearOperator` with the `matvec` and
  607. `rmatvec` methods (to apply the matrix and its adjoint).
  608. eps_or_k : float or int
  609. Relative error (if ``eps_or_k < 1``) or rank (if ``eps_or_k >= 1``) of
  610. approximation.
  611. rand : bool, optional
  612. Whether to use random sampling if `A` is of type :class:`numpy.ndarray`
  613. (randomized algorithms are always used if `A` is of type
  614. :class:`scipy.sparse.linalg.LinearOperator`).
  615. rng : `numpy.random.Generator`, optional
  616. Pseudorandom number generator state. When `rng` is None, a new
  617. `numpy.random.Generator` is created using entropy from the
  618. operating system. Types other than `numpy.random.Generator` are
  619. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  620. If `rand` is ``False``, the argument is ignored.
  621. Returns
  622. -------
  623. U : :class:`numpy.ndarray`
  624. 2D array of left singular vectors.
  625. S : :class:`numpy.ndarray`
  626. 1D array of singular values.
  627. V : :class:`numpy.ndarray`
  628. 2D array right singular vectors.
  629. """
  630. from scipy.sparse.linalg import LinearOperator
  631. rng = np.random.default_rng(rng)
  632. real = _is_real(A)
  633. if isinstance(A, np.ndarray):
  634. A = _C_contiguous_copy(A)
  635. if eps_or_k < 1:
  636. eps = eps_or_k
  637. if rand:
  638. if real:
  639. U, S, V = _backend.iddp_asvd(A, eps, rng=rng)
  640. else:
  641. U, S, V = _backend.idzp_asvd(A, eps, rng=rng)
  642. else:
  643. if real:
  644. U, S, V = _backend.iddp_svd(A, eps)
  645. V = V.T.conj()
  646. else:
  647. U, S, V = _backend.idzp_svd(A, eps)
  648. V = V.T.conj()
  649. else:
  650. k = int(eps_or_k)
  651. if k > min(A.shape):
  652. raise ValueError(f"Approximation rank {k} exceeds min(A.shape) = "
  653. f" {min(A.shape)} ")
  654. if rand:
  655. if real:
  656. U, S, V = _backend.iddr_asvd(A, k, rng=rng)
  657. else:
  658. U, S, V = _backend.idzr_asvd(A, k, rng=rng)
  659. else:
  660. if real:
  661. U, S, V = _backend.iddr_svd(A, k)
  662. V = V.T.conj()
  663. else:
  664. U, S, V = _backend.idzr_svd(A, k)
  665. V = V.T.conj()
  666. elif isinstance(A, LinearOperator):
  667. if eps_or_k < 1:
  668. eps = eps_or_k
  669. if real:
  670. U, S, V = _backend.iddp_rsvd(A, eps, rng=rng)
  671. else:
  672. U, S, V = _backend.idzp_rsvd(A, eps, rng=rng)
  673. else:
  674. k = int(eps_or_k)
  675. if real:
  676. U, S, V = _backend.iddr_rsvd(A, k, rng=rng)
  677. else:
  678. U, S, V = _backend.idzr_rsvd(A, k, rng=rng)
  679. else:
  680. raise _TYPE_ERROR
  681. return U, S, V
  682. def estimate_rank(A, eps, rng=None):
  683. """
  684. Estimate matrix rank to a specified relative precision using randomized
  685. methods.
  686. The matrix `A` can be given as either a :class:`numpy.ndarray` or a
  687. :class:`scipy.sparse.linalg.LinearOperator`, with different algorithms used
  688. for each case. If `A` is of type :class:`numpy.ndarray`, then the output
  689. rank is typically about 8 higher than the actual numerical rank.
  690. .. This function automatically detects the form of the input parameters and
  691. passes them to the appropriate backend. For details,
  692. see :func:`_backend.idd_estrank`, :func:`_backend.idd_findrank`,
  693. :func:`_backend.idz_estrank`, and :func:`_backend.idz_findrank`.
  694. Parameters
  695. ----------
  696. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator`
  697. Matrix whose rank is to be estimated, given as either a
  698. :class:`numpy.ndarray` or a :class:`scipy.sparse.linalg.LinearOperator`
  699. with the `rmatvec` method (to apply the matrix adjoint).
  700. eps : float
  701. Relative error for numerical rank definition.
  702. rng : `numpy.random.Generator`, optional
  703. Pseudorandom number generator state. When `rng` is None, a new
  704. `numpy.random.Generator` is created using entropy from the
  705. operating system. Types other than `numpy.random.Generator` are
  706. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  707. If `rand` is ``False``, the argument is ignored.
  708. Returns
  709. -------
  710. int
  711. Estimated matrix rank.
  712. """
  713. from scipy.sparse.linalg import LinearOperator
  714. rng = np.random.default_rng(rng)
  715. real = _is_real(A)
  716. if isinstance(A, np.ndarray):
  717. A = _C_contiguous_copy(A)
  718. if real:
  719. rank, _ = _backend.idd_estrank(A, eps, rng=rng)
  720. else:
  721. rank, _ = _backend.idz_estrank(A, eps, rng=rng)
  722. if rank == 0:
  723. # special return value for nearly full rank
  724. rank = min(A.shape)
  725. return rank
  726. elif isinstance(A, LinearOperator):
  727. if real:
  728. return _backend.idd_findrank(A, eps, rng=rng)[0]
  729. else:
  730. return _backend.idz_findrank(A, eps, rng=rng)[0]
  731. else:
  732. raise _TYPE_ERROR