interpolative.py 31 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004
  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(n):
  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 ``numpy.random.rand`` with the given integer
  253. seed.
  254. To simply generate some random numbers, type:
  255. >>> sli.rand(n)
  256. where ``n`` is the number of random numbers to generate.
  257. Remarks
  258. -------
  259. The above functions all automatically detect the appropriate interface and work
  260. with both real and complex data types, passing input arguments to the proper
  261. backend routine.
  262. """
  263. import scipy.linalg._interpolative_backend as _backend
  264. import numpy as np
  265. import sys
  266. __all__ = [
  267. 'estimate_rank',
  268. 'estimate_spectral_norm',
  269. 'estimate_spectral_norm_diff',
  270. 'id_to_svd',
  271. 'interp_decomp',
  272. 'rand',
  273. 'reconstruct_interp_matrix',
  274. 'reconstruct_matrix_from_id',
  275. 'reconstruct_skel_matrix',
  276. 'seed',
  277. 'svd',
  278. ]
  279. _DTYPE_ERROR = ValueError("invalid input dtype (input must be float64 or complex128)")
  280. _TYPE_ERROR = TypeError("invalid input type (must be array or LinearOperator)")
  281. _32BIT_ERROR = ValueError("interpolative decomposition on 32-bit systems "
  282. "with complex128 is buggy")
  283. _IS_32BIT = (sys.maxsize < 2**32)
  284. def _is_real(A):
  285. try:
  286. if A.dtype == np.complex128:
  287. return False
  288. elif A.dtype == np.float64:
  289. return True
  290. else:
  291. raise _DTYPE_ERROR
  292. except AttributeError as e:
  293. raise _TYPE_ERROR from e
  294. def seed(seed=None):
  295. """
  296. Seed the internal random number generator used in this ID package.
  297. The generator is a lagged Fibonacci method with 55-element internal state.
  298. Parameters
  299. ----------
  300. seed : int, sequence, 'default', optional
  301. If 'default', the random seed is reset to a default value.
  302. If `seed` is a sequence containing 55 floating-point numbers
  303. in range [0,1], these are used to set the internal state of
  304. the generator.
  305. If the value is an integer, the internal state is obtained
  306. from `numpy.random.RandomState` (MT19937) with the integer
  307. used as the initial seed.
  308. If `seed` is omitted (None), ``numpy.random.rand`` is used to
  309. initialize the generator.
  310. """
  311. # For details, see :func:`_backend.id_srand`, :func:`_backend.id_srandi`,
  312. # and :func:`_backend.id_srando`.
  313. if isinstance(seed, str) and seed == 'default':
  314. _backend.id_srando()
  315. elif hasattr(seed, '__len__'):
  316. state = np.asfortranarray(seed, dtype=float)
  317. if state.shape != (55,):
  318. raise ValueError("invalid input size")
  319. elif state.min() < 0 or state.max() > 1:
  320. raise ValueError("values not in range [0,1]")
  321. _backend.id_srandi(state)
  322. elif seed is None:
  323. _backend.id_srandi(np.random.rand(55))
  324. else:
  325. rnd = np.random.RandomState(seed)
  326. _backend.id_srandi(rnd.rand(55))
  327. def rand(*shape):
  328. """
  329. Generate standard uniform pseudorandom numbers via a very efficient lagged
  330. Fibonacci method.
  331. This routine is used for all random number generation in this package and
  332. can affect ID and SVD results.
  333. Parameters
  334. ----------
  335. *shape
  336. Shape of output array
  337. """
  338. # For details, see :func:`_backend.id_srand`, and :func:`_backend.id_srando`.
  339. return _backend.id_srand(np.prod(shape)).reshape(shape)
  340. def interp_decomp(A, eps_or_k, rand=True):
  341. """
  342. Compute ID of a matrix.
  343. An ID of a matrix `A` is a factorization defined by a rank `k`, a column
  344. index array `idx`, and interpolation coefficients `proj` such that::
  345. numpy.dot(A[:,idx[:k]], proj) = A[:,idx[k:]]
  346. The original matrix can then be reconstructed as::
  347. numpy.hstack([A[:,idx[:k]],
  348. numpy.dot(A[:,idx[:k]], proj)]
  349. )[:,numpy.argsort(idx)]
  350. or via the routine :func:`reconstruct_matrix_from_id`. This can
  351. equivalently be written as::
  352. numpy.dot(A[:,idx[:k]],
  353. numpy.hstack([numpy.eye(k), proj])
  354. )[:,np.argsort(idx)]
  355. in terms of the skeleton and interpolation matrices::
  356. B = A[:,idx[:k]]
  357. and::
  358. P = numpy.hstack([numpy.eye(k), proj])[:,np.argsort(idx)]
  359. respectively. See also :func:`reconstruct_interp_matrix` and
  360. :func:`reconstruct_skel_matrix`.
  361. The ID can be computed to any relative precision or rank (depending on the
  362. value of `eps_or_k`). If a precision is specified (`eps_or_k < 1`), then
  363. this function has the output signature::
  364. k, idx, proj = interp_decomp(A, eps_or_k)
  365. Otherwise, if a rank is specified (`eps_or_k >= 1`), then the output
  366. signature is::
  367. idx, proj = interp_decomp(A, eps_or_k)
  368. .. This function automatically detects the form of the input parameters
  369. and passes them to the appropriate backend. For details, see
  370. :func:`_backend.iddp_id`, :func:`_backend.iddp_aid`,
  371. :func:`_backend.iddp_rid`, :func:`_backend.iddr_id`,
  372. :func:`_backend.iddr_aid`, :func:`_backend.iddr_rid`,
  373. :func:`_backend.idzp_id`, :func:`_backend.idzp_aid`,
  374. :func:`_backend.idzp_rid`, :func:`_backend.idzr_id`,
  375. :func:`_backend.idzr_aid`, and :func:`_backend.idzr_rid`.
  376. Parameters
  377. ----------
  378. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator` with `rmatvec`
  379. Matrix to be factored
  380. eps_or_k : float or int
  381. Relative error (if `eps_or_k < 1`) or rank (if `eps_or_k >= 1`) of
  382. approximation.
  383. rand : bool, optional
  384. Whether to use random sampling if `A` is of type :class:`numpy.ndarray`
  385. (randomized algorithms are always used if `A` is of type
  386. :class:`scipy.sparse.linalg.LinearOperator`).
  387. Returns
  388. -------
  389. k : int
  390. Rank required to achieve specified relative precision if
  391. `eps_or_k < 1`.
  392. idx : :class:`numpy.ndarray`
  393. Column index array.
  394. proj : :class:`numpy.ndarray`
  395. Interpolation coefficients.
  396. """
  397. from scipy.sparse.linalg import LinearOperator
  398. real = _is_real(A)
  399. if isinstance(A, np.ndarray):
  400. if eps_or_k < 1:
  401. eps = eps_or_k
  402. if rand:
  403. if real:
  404. k, idx, proj = _backend.iddp_aid(eps, A)
  405. else:
  406. if _IS_32BIT:
  407. raise _32BIT_ERROR
  408. k, idx, proj = _backend.idzp_aid(eps, A)
  409. else:
  410. if real:
  411. k, idx, proj = _backend.iddp_id(eps, A)
  412. else:
  413. k, idx, proj = _backend.idzp_id(eps, A)
  414. return k, idx - 1, proj
  415. else:
  416. k = int(eps_or_k)
  417. if rand:
  418. if real:
  419. idx, proj = _backend.iddr_aid(A, k)
  420. else:
  421. if _IS_32BIT:
  422. raise _32BIT_ERROR
  423. idx, proj = _backend.idzr_aid(A, k)
  424. else:
  425. if real:
  426. idx, proj = _backend.iddr_id(A, k)
  427. else:
  428. idx, proj = _backend.idzr_id(A, k)
  429. return idx - 1, proj
  430. elif isinstance(A, LinearOperator):
  431. m, n = A.shape
  432. matveca = A.rmatvec
  433. if eps_or_k < 1:
  434. eps = eps_or_k
  435. if real:
  436. k, idx, proj = _backend.iddp_rid(eps, m, n, matveca)
  437. else:
  438. if _IS_32BIT:
  439. raise _32BIT_ERROR
  440. k, idx, proj = _backend.idzp_rid(eps, m, n, matveca)
  441. return k, idx - 1, proj
  442. else:
  443. k = int(eps_or_k)
  444. if real:
  445. idx, proj = _backend.iddr_rid(m, n, matveca, k)
  446. else:
  447. if _IS_32BIT:
  448. raise _32BIT_ERROR
  449. idx, proj = _backend.idzr_rid(m, n, matveca, k)
  450. return idx - 1, proj
  451. else:
  452. raise _TYPE_ERROR
  453. def reconstruct_matrix_from_id(B, idx, proj):
  454. """
  455. Reconstruct matrix from its ID.
  456. A matrix `A` with skeleton matrix `B` and ID indices and coefficients `idx`
  457. and `proj`, respectively, can be reconstructed as::
  458. numpy.hstack([B, numpy.dot(B, proj)])[:,numpy.argsort(idx)]
  459. See also :func:`reconstruct_interp_matrix` and
  460. :func:`reconstruct_skel_matrix`.
  461. .. This function automatically detects the matrix data type and calls the
  462. appropriate backend. For details, see :func:`_backend.idd_reconid` and
  463. :func:`_backend.idz_reconid`.
  464. Parameters
  465. ----------
  466. B : :class:`numpy.ndarray`
  467. Skeleton matrix.
  468. idx : :class:`numpy.ndarray`
  469. Column index array.
  470. proj : :class:`numpy.ndarray`
  471. Interpolation coefficients.
  472. Returns
  473. -------
  474. :class:`numpy.ndarray`
  475. Reconstructed matrix.
  476. """
  477. if _is_real(B):
  478. return _backend.idd_reconid(B, idx + 1, proj)
  479. else:
  480. return _backend.idz_reconid(B, idx + 1, proj)
  481. def reconstruct_interp_matrix(idx, proj):
  482. """
  483. Reconstruct interpolation matrix from ID.
  484. The interpolation matrix can be reconstructed from the ID indices and
  485. coefficients `idx` and `proj`, respectively, as::
  486. P = numpy.hstack([numpy.eye(proj.shape[0]), proj])[:,numpy.argsort(idx)]
  487. The original matrix can then be reconstructed from its skeleton matrix `B`
  488. via::
  489. numpy.dot(B, P)
  490. See also :func:`reconstruct_matrix_from_id` and
  491. :func:`reconstruct_skel_matrix`.
  492. .. This function automatically detects the matrix data type and calls the
  493. appropriate backend. For details, see :func:`_backend.idd_reconint` and
  494. :func:`_backend.idz_reconint`.
  495. Parameters
  496. ----------
  497. idx : :class:`numpy.ndarray`
  498. Column index array.
  499. proj : :class:`numpy.ndarray`
  500. Interpolation coefficients.
  501. Returns
  502. -------
  503. :class:`numpy.ndarray`
  504. Interpolation matrix.
  505. """
  506. if _is_real(proj):
  507. return _backend.idd_reconint(idx + 1, proj)
  508. else:
  509. return _backend.idz_reconint(idx + 1, proj)
  510. def reconstruct_skel_matrix(A, k, idx):
  511. """
  512. Reconstruct skeleton matrix from ID.
  513. The skeleton matrix can be reconstructed from the original matrix `A` and its
  514. ID rank and indices `k` and `idx`, respectively, as::
  515. B = A[:,idx[:k]]
  516. The original matrix can then be reconstructed via::
  517. numpy.hstack([B, numpy.dot(B, proj)])[:,numpy.argsort(idx)]
  518. See also :func:`reconstruct_matrix_from_id` and
  519. :func:`reconstruct_interp_matrix`.
  520. .. This function automatically detects the matrix data type and calls the
  521. appropriate backend. For details, see :func:`_backend.idd_copycols` and
  522. :func:`_backend.idz_copycols`.
  523. Parameters
  524. ----------
  525. A : :class:`numpy.ndarray`
  526. Original matrix.
  527. k : int
  528. Rank of ID.
  529. idx : :class:`numpy.ndarray`
  530. Column index array.
  531. Returns
  532. -------
  533. :class:`numpy.ndarray`
  534. Skeleton matrix.
  535. """
  536. if _is_real(A):
  537. return _backend.idd_copycols(A, k, idx + 1)
  538. else:
  539. return _backend.idz_copycols(A, k, idx + 1)
  540. def id_to_svd(B, idx, proj):
  541. """
  542. Convert ID to SVD.
  543. The SVD reconstruction of a matrix with skeleton matrix `B` and ID indices and
  544. coefficients `idx` and `proj`, respectively, is::
  545. U, S, V = id_to_svd(B, idx, proj)
  546. A = numpy.dot(U, numpy.dot(numpy.diag(S), V.conj().T))
  547. See also :func:`svd`.
  548. .. This function automatically detects the matrix data type and calls the
  549. appropriate backend. For details, see :func:`_backend.idd_id2svd` and
  550. :func:`_backend.idz_id2svd`.
  551. Parameters
  552. ----------
  553. B : :class:`numpy.ndarray`
  554. Skeleton matrix.
  555. idx : :class:`numpy.ndarray`
  556. Column index array.
  557. proj : :class:`numpy.ndarray`
  558. Interpolation coefficients.
  559. Returns
  560. -------
  561. U : :class:`numpy.ndarray`
  562. Left singular vectors.
  563. S : :class:`numpy.ndarray`
  564. Singular values.
  565. V : :class:`numpy.ndarray`
  566. Right singular vectors.
  567. """
  568. if _is_real(B):
  569. U, V, S = _backend.idd_id2svd(B, idx + 1, proj)
  570. else:
  571. U, V, S = _backend.idz_id2svd(B, idx + 1, proj)
  572. return U, S, V
  573. def estimate_spectral_norm(A, its=20):
  574. """
  575. Estimate spectral norm of a matrix by the randomized power method.
  576. .. This function automatically detects the matrix data type and calls the
  577. appropriate backend. For details, see :func:`_backend.idd_snorm` and
  578. :func:`_backend.idz_snorm`.
  579. Parameters
  580. ----------
  581. A : :class:`scipy.sparse.linalg.LinearOperator`
  582. Matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with the
  583. `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  584. its : int, optional
  585. Number of power method iterations.
  586. Returns
  587. -------
  588. float
  589. Spectral norm estimate.
  590. """
  591. from scipy.sparse.linalg import aslinearoperator
  592. A = aslinearoperator(A)
  593. m, n = A.shape
  594. matvec = lambda x: A. matvec(x)
  595. matveca = lambda x: A.rmatvec(x)
  596. if _is_real(A):
  597. return _backend.idd_snorm(m, n, matveca, matvec, its=its)
  598. else:
  599. return _backend.idz_snorm(m, n, matveca, matvec, its=its)
  600. def estimate_spectral_norm_diff(A, B, its=20):
  601. """
  602. Estimate spectral norm of the difference of two matrices by the randomized
  603. power method.
  604. .. This function automatically detects the matrix data type and calls the
  605. appropriate backend. For details, see :func:`_backend.idd_diffsnorm` and
  606. :func:`_backend.idz_diffsnorm`.
  607. Parameters
  608. ----------
  609. A : :class:`scipy.sparse.linalg.LinearOperator`
  610. First matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with the
  611. `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  612. B : :class:`scipy.sparse.linalg.LinearOperator`
  613. Second matrix given as a :class:`scipy.sparse.linalg.LinearOperator` with
  614. the `matvec` and `rmatvec` methods (to apply the matrix and its adjoint).
  615. its : int, optional
  616. Number of power method iterations.
  617. Returns
  618. -------
  619. float
  620. Spectral norm estimate of matrix difference.
  621. """
  622. from scipy.sparse.linalg import aslinearoperator
  623. A = aslinearoperator(A)
  624. B = aslinearoperator(B)
  625. m, n = A.shape
  626. matvec1 = lambda x: A. matvec(x)
  627. matveca1 = lambda x: A.rmatvec(x)
  628. matvec2 = lambda x: B. matvec(x)
  629. matveca2 = lambda x: B.rmatvec(x)
  630. if _is_real(A):
  631. return _backend.idd_diffsnorm(
  632. m, n, matveca1, matveca2, matvec1, matvec2, its=its)
  633. else:
  634. return _backend.idz_diffsnorm(
  635. m, n, matveca1, matveca2, matvec1, matvec2, its=its)
  636. def svd(A, eps_or_k, rand=True):
  637. """
  638. Compute SVD of a matrix via an ID.
  639. An SVD of a matrix `A` is a factorization::
  640. A = numpy.dot(U, numpy.dot(numpy.diag(S), V.conj().T))
  641. where `U` and `V` have orthonormal columns and `S` is nonnegative.
  642. The SVD can be computed to any relative precision or rank (depending on the
  643. value of `eps_or_k`).
  644. See also :func:`interp_decomp` and :func:`id_to_svd`.
  645. .. This function automatically detects the form of the input parameters and
  646. passes them to the appropriate backend. For details, see
  647. :func:`_backend.iddp_svd`, :func:`_backend.iddp_asvd`,
  648. :func:`_backend.iddp_rsvd`, :func:`_backend.iddr_svd`,
  649. :func:`_backend.iddr_asvd`, :func:`_backend.iddr_rsvd`,
  650. :func:`_backend.idzp_svd`, :func:`_backend.idzp_asvd`,
  651. :func:`_backend.idzp_rsvd`, :func:`_backend.idzr_svd`,
  652. :func:`_backend.idzr_asvd`, and :func:`_backend.idzr_rsvd`.
  653. Parameters
  654. ----------
  655. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator`
  656. Matrix to be factored, given as either a :class:`numpy.ndarray` or a
  657. :class:`scipy.sparse.linalg.LinearOperator` with the `matvec` and
  658. `rmatvec` methods (to apply the matrix and its adjoint).
  659. eps_or_k : float or int
  660. Relative error (if `eps_or_k < 1`) or rank (if `eps_or_k >= 1`) of
  661. approximation.
  662. rand : bool, optional
  663. Whether to use random sampling if `A` is of type :class:`numpy.ndarray`
  664. (randomized algorithms are always used if `A` is of type
  665. :class:`scipy.sparse.linalg.LinearOperator`).
  666. Returns
  667. -------
  668. U : :class:`numpy.ndarray`
  669. Left singular vectors.
  670. S : :class:`numpy.ndarray`
  671. Singular values.
  672. V : :class:`numpy.ndarray`
  673. Right singular vectors.
  674. """
  675. from scipy.sparse.linalg import LinearOperator
  676. real = _is_real(A)
  677. if isinstance(A, np.ndarray):
  678. if eps_or_k < 1:
  679. eps = eps_or_k
  680. if rand:
  681. if real:
  682. U, V, S = _backend.iddp_asvd(eps, A)
  683. else:
  684. if _IS_32BIT:
  685. raise _32BIT_ERROR
  686. U, V, S = _backend.idzp_asvd(eps, A)
  687. else:
  688. if real:
  689. U, V, S = _backend.iddp_svd(eps, A)
  690. else:
  691. U, V, S = _backend.idzp_svd(eps, A)
  692. else:
  693. k = int(eps_or_k)
  694. if k > min(A.shape):
  695. raise ValueError("Approximation rank %s exceeds min(A.shape) = "
  696. " %s " % (k, min(A.shape)))
  697. if rand:
  698. if real:
  699. U, V, S = _backend.iddr_asvd(A, k)
  700. else:
  701. if _IS_32BIT:
  702. raise _32BIT_ERROR
  703. U, V, S = _backend.idzr_asvd(A, k)
  704. else:
  705. if real:
  706. U, V, S = _backend.iddr_svd(A, k)
  707. else:
  708. U, V, S = _backend.idzr_svd(A, k)
  709. elif isinstance(A, LinearOperator):
  710. m, n = A.shape
  711. matvec = lambda x: A.matvec(x)
  712. matveca = lambda x: A.rmatvec(x)
  713. if eps_or_k < 1:
  714. eps = eps_or_k
  715. if real:
  716. U, V, S = _backend.iddp_rsvd(eps, m, n, matveca, matvec)
  717. else:
  718. if _IS_32BIT:
  719. raise _32BIT_ERROR
  720. U, V, S = _backend.idzp_rsvd(eps, m, n, matveca, matvec)
  721. else:
  722. k = int(eps_or_k)
  723. if real:
  724. U, V, S = _backend.iddr_rsvd(m, n, matveca, matvec, k)
  725. else:
  726. if _IS_32BIT:
  727. raise _32BIT_ERROR
  728. U, V, S = _backend.idzr_rsvd(m, n, matveca, matvec, k)
  729. else:
  730. raise _TYPE_ERROR
  731. return U, S, V
  732. def estimate_rank(A, eps):
  733. """
  734. Estimate matrix rank to a specified relative precision using randomized
  735. methods.
  736. The matrix `A` can be given as either a :class:`numpy.ndarray` or a
  737. :class:`scipy.sparse.linalg.LinearOperator`, with different algorithms used
  738. for each case. If `A` is of type :class:`numpy.ndarray`, then the output
  739. rank is typically about 8 higher than the actual numerical rank.
  740. .. This function automatically detects the form of the input parameters and
  741. passes them to the appropriate backend. For details,
  742. see :func:`_backend.idd_estrank`, :func:`_backend.idd_findrank`,
  743. :func:`_backend.idz_estrank`, and :func:`_backend.idz_findrank`.
  744. Parameters
  745. ----------
  746. A : :class:`numpy.ndarray` or :class:`scipy.sparse.linalg.LinearOperator`
  747. Matrix whose rank is to be estimated, given as either a
  748. :class:`numpy.ndarray` or a :class:`scipy.sparse.linalg.LinearOperator`
  749. with the `rmatvec` method (to apply the matrix adjoint).
  750. eps : float
  751. Relative error for numerical rank definition.
  752. Returns
  753. -------
  754. int
  755. Estimated matrix rank.
  756. """
  757. from scipy.sparse.linalg import LinearOperator
  758. real = _is_real(A)
  759. if isinstance(A, np.ndarray):
  760. if real:
  761. rank = _backend.idd_estrank(eps, A)
  762. else:
  763. rank = _backend.idz_estrank(eps, A)
  764. if rank == 0:
  765. # special return value for nearly full rank
  766. rank = min(A.shape)
  767. return rank
  768. elif isinstance(A, LinearOperator):
  769. m, n = A.shape
  770. matveca = A.rmatvec
  771. if real:
  772. return _backend.idd_findrank(eps, m, n, matveca)
  773. else:
  774. return _backend.idz_findrank(eps, m, n, matveca)
  775. else:
  776. raise _TYPE_ERROR