interpolative.py 30 KB

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