test_interpolative.py 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232
  1. # ******************************************************************************
  2. # Copyright (C) 2013 Kenneth L. Ho
  3. # Redistribution and use in source and binary forms, with or without
  4. # modification, are permitted provided that the following conditions are met:
  5. #
  6. # Redistributions of source code must retain the above copyright notice, this
  7. # list of conditions and the following disclaimer. Redistributions in binary
  8. # form must reproduce the above copyright notice, this list of conditions and
  9. # the following disclaimer in the documentation and/or other materials
  10. # provided with the distribution.
  11. #
  12. # None of the names of the copyright holders may be used to endorse or
  13. # promote products derived from this software without specific prior written
  14. # permission.
  15. #
  16. # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  17. # AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  18. # IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  19. # ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
  20. # LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  21. # CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  22. # SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  23. # INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  24. # CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  25. # ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  26. # POSSIBILITY OF SUCH DAMAGE.
  27. # ******************************************************************************
  28. import scipy.linalg.interpolative as pymatrixid
  29. import numpy as np
  30. from scipy.linalg import hilbert, svdvals, norm
  31. from scipy.sparse.linalg import aslinearoperator
  32. from scipy.linalg.interpolative import interp_decomp
  33. from numpy.testing import (assert_, assert_allclose, assert_equal,
  34. assert_array_equal)
  35. import pytest
  36. from pytest import raises as assert_raises
  37. @pytest.fixture()
  38. def eps():
  39. yield 1e-12
  40. @pytest.fixture()
  41. def rng():
  42. rng = np.random.default_rng(1718313768084012)
  43. yield rng
  44. @pytest.fixture(params=[np.float64, np.complex128])
  45. def A(request):
  46. # construct Hilbert matrix
  47. # set parameters
  48. n = 300
  49. yield hilbert(n).astype(request.param)
  50. @pytest.fixture()
  51. def L(A):
  52. yield aslinearoperator(A)
  53. @pytest.fixture()
  54. def rank(A, eps):
  55. S = np.linalg.svd(A, compute_uv=False)
  56. try:
  57. rank = np.nonzero(S < eps)[0][0]
  58. except IndexError:
  59. rank = A.shape[0]
  60. return rank
  61. class TestInterpolativeDecomposition:
  62. @pytest.mark.parametrize(
  63. "rand,lin_op",
  64. [(False, False), (True, False), (True, True)])
  65. def test_real_id_fixed_precision(self, A, L, eps, rand, lin_op, rng):
  66. # Test ID routines on a Hilbert matrix.
  67. A_or_L = A if not lin_op else L
  68. k, idx, proj = pymatrixid.interp_decomp(A_or_L, eps, rand=rand, rng=rng)
  69. B = pymatrixid.reconstruct_matrix_from_id(A[:, idx[:k]], idx, proj)
  70. assert_allclose(A, B, rtol=eps, atol=1e-08)
  71. @pytest.mark.parametrize(
  72. "rand,lin_op",
  73. [(False, False), (True, False), (True, True)])
  74. def test_real_id_fixed_rank(self, A, L, eps, rank, rand, lin_op, rng):
  75. k = rank
  76. A_or_L = A if not lin_op else L
  77. idx, proj = pymatrixid.interp_decomp(A_or_L, k, rand=rand, rng=rng)
  78. B = pymatrixid.reconstruct_matrix_from_id(A[:, idx[:k]], idx, proj)
  79. assert_allclose(A, B, rtol=eps, atol=1e-08)
  80. @pytest.mark.parametrize("rand,lin_op", [(False, False)])
  81. def test_real_id_skel_and_interp_matrices(
  82. self, A, L, eps, rank, rand, lin_op, rng):
  83. k = rank
  84. A_or_L = A if not lin_op else L
  85. idx, proj = pymatrixid.interp_decomp(A_or_L, k, rand=rand, rng=rng)
  86. P = pymatrixid.reconstruct_interp_matrix(idx, proj)
  87. B = pymatrixid.reconstruct_skel_matrix(A, k, idx)
  88. assert_allclose(B, A[:, idx[:k]], rtol=eps, atol=1e-08)
  89. assert_allclose(B @ P, A, rtol=eps, atol=1e-08)
  90. @pytest.mark.parametrize(
  91. "rand,lin_op",
  92. [(False, False), (True, False), (True, True)])
  93. def test_svd_fixed_precision(self, A, L, eps, rand, lin_op, rng):
  94. A_or_L = A if not lin_op else L
  95. U, S, V = pymatrixid.svd(A_or_L, eps, rand=rand, rng=rng)
  96. B = U * S @ V.T.conj()
  97. assert_allclose(A, B, rtol=eps, atol=1e-08)
  98. @pytest.mark.parametrize(
  99. "rand,lin_op",
  100. [(False, False), (True, False), (True, True)])
  101. def test_svd_fixed_rank(self, A, L, eps, rank, rand, lin_op, rng):
  102. k = rank
  103. A_or_L = A if not lin_op else L
  104. U, S, V = pymatrixid.svd(A_or_L, k, rand=rand, rng=rng)
  105. B = U * S @ V.T.conj()
  106. assert_allclose(A, B, rtol=eps, atol=1e-08)
  107. def test_id_to_svd(self, A, eps, rank):
  108. k = rank
  109. idx, proj = pymatrixid.interp_decomp(A, k, rand=False)
  110. U, S, V = pymatrixid.id_to_svd(A[:, idx[:k]], idx, proj)
  111. B = U * S @ V.T.conj()
  112. assert_allclose(A, B, rtol=eps, atol=1e-08)
  113. def test_estimate_spectral_norm(self, A, rng):
  114. s = svdvals(A)
  115. norm_2_est = pymatrixid.estimate_spectral_norm(A, rng=rng)
  116. assert_allclose(norm_2_est, s[0], rtol=1e-6, atol=1e-8)
  117. def test_estimate_spectral_norm_diff(self, A, rng):
  118. B = A.copy()
  119. B[:, 0] *= 1.2
  120. s = svdvals(A - B)
  121. norm_2_est = pymatrixid.estimate_spectral_norm_diff(A, B, rng=rng)
  122. assert_allclose(norm_2_est, s[0], rtol=1e-6, atol=1e-8)
  123. def test_rank_estimates_array(self, A, rng):
  124. B = np.array([[1, 1, 0], [0, 0, 1], [0, 0, 1]], dtype=A.dtype)
  125. for M in [A, B]:
  126. rank_tol = 1e-9
  127. rank_np = np.linalg.matrix_rank(M, norm(M, 2) * rank_tol)
  128. rank_est = pymatrixid.estimate_rank(M, rank_tol, rng=rng)
  129. assert_(rank_est >= rank_np)
  130. assert_(rank_est <= rank_np + 10)
  131. def test_rank_estimates_lin_op(self, A, rng):
  132. B = np.array([[1, 1, 0], [0, 0, 1], [0, 0, 1]], dtype=A.dtype)
  133. for M in [A, B]:
  134. ML = aslinearoperator(M)
  135. rank_tol = 1e-9
  136. rank_np = np.linalg.matrix_rank(M, norm(M, 2) * rank_tol)
  137. rank_est = pymatrixid.estimate_rank(ML, rank_tol, rng=rng)
  138. assert_(rank_est >= rank_np - 4)
  139. assert_(rank_est <= rank_np + 4)
  140. def test_badcall(self):
  141. A = hilbert(5).astype(np.float32)
  142. with assert_raises(ValueError):
  143. pymatrixid.interp_decomp(A, 1e-6, rand=False)
  144. def test_rank_too_large(self):
  145. # svd(array, k) should not segfault
  146. a = np.ones((4, 3))
  147. with assert_raises(ValueError):
  148. pymatrixid.svd(a, 4)
  149. def test_full_rank(self):
  150. eps = 1.0e-12
  151. rng = np.random.default_rng(1234)
  152. # fixed precision
  153. A = rng.random((16, 8))
  154. k, idx, proj = pymatrixid.interp_decomp(A, eps)
  155. assert_equal(k, A.shape[1])
  156. P = pymatrixid.reconstruct_interp_matrix(idx, proj)
  157. B = pymatrixid.reconstruct_skel_matrix(A, k, idx)
  158. assert_allclose(A, B @ P)
  159. # fixed rank
  160. idx, proj = pymatrixid.interp_decomp(A, k)
  161. P = pymatrixid.reconstruct_interp_matrix(idx, proj)
  162. B = pymatrixid.reconstruct_skel_matrix(A, k, idx)
  163. assert_allclose(A, B @ P)
  164. @pytest.mark.parametrize("dtype", [np.float64, np.complex128])
  165. @pytest.mark.parametrize("rand", [True, False])
  166. @pytest.mark.parametrize("eps", [1, 0.1])
  167. def test_bug_9793(self, dtype, rand, eps):
  168. A = np.array([[-1, -1, -1, 0, 0, 0],
  169. [0, 0, 0, 1, 1, 1],
  170. [1, 0, 0, 1, 0, 0],
  171. [0, 1, 0, 0, 1, 0],
  172. [0, 0, 1, 0, 0, 1]],
  173. dtype=dtype, order="C")
  174. B = A.copy()
  175. interp_decomp(A.T, eps, rand=rand)
  176. assert_array_equal(A, B)
  177. def test_svd_aslinearoperator_shape_check(self):
  178. # See gh-issue #22451
  179. rng = np.random.default_rng(1744580941832515)
  180. x = rng.uniform(size=[7, 5])
  181. xl = aslinearoperator(x)
  182. u, s, v = pymatrixid.svd(xl, 3)
  183. assert_equal(u.shape, (7, 3))
  184. assert_equal(s.shape, (3,))
  185. assert_equal(v.shape, (5, 3))
  186. x = rng.uniform(size=[4, 9])
  187. xl = aslinearoperator(x)
  188. u, s, v = pymatrixid.svd(xl, 2)
  189. assert_equal(u.shape, (4, 2))
  190. assert_equal(s.shape, (2,))
  191. assert_equal(v.shape, (9, 2))