test_lgmres.py 6.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211
  1. """Tests for the linalg._isolve.lgmres module
  2. """
  3. from numpy.testing import (assert_, assert_allclose, assert_equal,
  4. suppress_warnings)
  5. import pytest
  6. from platform import python_implementation
  7. import numpy as np
  8. from numpy import zeros, array, allclose
  9. from scipy.linalg import norm
  10. from scipy.sparse import csr_matrix, eye, rand
  11. from scipy.sparse.linalg._interface import LinearOperator
  12. from scipy.sparse.linalg import splu
  13. from scipy.sparse.linalg._isolve import lgmres, gmres
  14. Am = csr_matrix(array([[-2, 1, 0, 0, 0, 9],
  15. [1, -2, 1, 0, 5, 0],
  16. [0, 1, -2, 1, 0, 0],
  17. [0, 0, 1, -2, 1, 0],
  18. [0, 3, 0, 1, -2, 1],
  19. [1, 0, 0, 0, 1, -2]]))
  20. b = array([1, 2, 3, 4, 5, 6])
  21. count = [0]
  22. def matvec(v):
  23. count[0] += 1
  24. return Am@v
  25. A = LinearOperator(matvec=matvec, shape=Am.shape, dtype=Am.dtype)
  26. def do_solve(**kw):
  27. count[0] = 0
  28. with suppress_warnings() as sup:
  29. sup.filter(DeprecationWarning, ".*called without specifying.*")
  30. x0, flag = lgmres(A, b, x0=zeros(A.shape[0]),
  31. inner_m=6, tol=1e-14, **kw)
  32. count_0 = count[0]
  33. assert_(allclose(A@x0, b, rtol=1e-12, atol=1e-12), norm(A@x0-b))
  34. return x0, count_0
  35. class TestLGMRES:
  36. def test_preconditioner(self):
  37. # Check that preconditioning works
  38. pc = splu(Am.tocsc())
  39. M = LinearOperator(matvec=pc.solve, shape=A.shape, dtype=A.dtype)
  40. x0, count_0 = do_solve()
  41. x1, count_1 = do_solve(M=M)
  42. assert_(count_1 == 3)
  43. assert_(count_1 < count_0/2)
  44. assert_(allclose(x1, x0, rtol=1e-14))
  45. def test_outer_v(self):
  46. # Check that the augmentation vectors behave as expected
  47. outer_v = []
  48. x0, count_0 = do_solve(outer_k=6, outer_v=outer_v)
  49. assert_(len(outer_v) > 0)
  50. assert_(len(outer_v) <= 6)
  51. x1, count_1 = do_solve(outer_k=6, outer_v=outer_v,
  52. prepend_outer_v=True)
  53. assert_(count_1 == 2, count_1)
  54. assert_(count_1 < count_0/2)
  55. assert_(allclose(x1, x0, rtol=1e-14))
  56. # ---
  57. outer_v = []
  58. x0, count_0 = do_solve(outer_k=6, outer_v=outer_v,
  59. store_outer_Av=False)
  60. assert_(array([v[1] is None for v in outer_v]).all())
  61. assert_(len(outer_v) > 0)
  62. assert_(len(outer_v) <= 6)
  63. x1, count_1 = do_solve(outer_k=6, outer_v=outer_v,
  64. prepend_outer_v=True)
  65. assert_(count_1 == 3, count_1)
  66. assert_(count_1 < count_0/2)
  67. assert_(allclose(x1, x0, rtol=1e-14))
  68. @pytest.mark.skipif(python_implementation() == 'PyPy',
  69. reason="Fails on PyPy CI runs. See #9507")
  70. def test_arnoldi(self):
  71. np.random.seed(1234)
  72. A = eye(2000) + rand(2000, 2000, density=5e-4)
  73. b = np.random.rand(2000)
  74. # The inner arnoldi should be equivalent to gmres
  75. with suppress_warnings() as sup:
  76. sup.filter(DeprecationWarning, ".*called without specifying.*")
  77. x0, flag0 = lgmres(A, b, x0=zeros(A.shape[0]),
  78. inner_m=15, maxiter=1)
  79. x1, flag1 = gmres(A, b, x0=zeros(A.shape[0]),
  80. restart=15, maxiter=1)
  81. assert_equal(flag0, 1)
  82. assert_equal(flag1, 1)
  83. norm = np.linalg.norm(A.dot(x0) - b)
  84. assert_(norm > 1e-4)
  85. assert_allclose(x0, x1)
  86. def test_cornercase(self):
  87. np.random.seed(1234)
  88. # Rounding error may prevent convergence with tol=0 --- ensure
  89. # that the return values in this case are correct, and no
  90. # exceptions are raised
  91. for n in [3, 5, 10, 100]:
  92. A = 2*eye(n)
  93. with suppress_warnings() as sup:
  94. sup.filter(DeprecationWarning, ".*called without specifying.*")
  95. b = np.ones(n)
  96. x, info = lgmres(A, b, maxiter=10)
  97. assert_equal(info, 0)
  98. assert_allclose(A.dot(x) - b, 0, atol=1e-14)
  99. x, info = lgmres(A, b, tol=0, maxiter=10)
  100. if info == 0:
  101. assert_allclose(A.dot(x) - b, 0, atol=1e-14)
  102. b = np.random.rand(n)
  103. x, info = lgmres(A, b, maxiter=10)
  104. assert_equal(info, 0)
  105. assert_allclose(A.dot(x) - b, 0, atol=1e-14)
  106. x, info = lgmres(A, b, tol=0, maxiter=10)
  107. if info == 0:
  108. assert_allclose(A.dot(x) - b, 0, atol=1e-14)
  109. def test_nans(self):
  110. A = eye(3, format='lil')
  111. A[1, 1] = np.nan
  112. b = np.ones(3)
  113. with suppress_warnings() as sup:
  114. sup.filter(DeprecationWarning, ".*called without specifying.*")
  115. x, info = lgmres(A, b, tol=0, maxiter=10)
  116. assert_equal(info, 1)
  117. def test_breakdown_with_outer_v(self):
  118. A = np.array([[1, 2], [3, 4]], dtype=float)
  119. b = np.array([1, 2])
  120. x = np.linalg.solve(A, b)
  121. v0 = np.array([1, 0])
  122. # The inner iteration should converge to the correct solution,
  123. # since it's in the outer vector list
  124. with suppress_warnings() as sup:
  125. sup.filter(DeprecationWarning, ".*called without specifying.*")
  126. xp, info = lgmres(A, b, outer_v=[(v0, None), (x, None)], maxiter=1)
  127. assert_allclose(xp, x, atol=1e-12)
  128. def test_breakdown_underdetermined(self):
  129. # Should find LSQ solution in the Krylov span in one inner
  130. # iteration, despite solver breakdown from nilpotent A.
  131. A = np.array([[0, 1, 1, 1],
  132. [0, 0, 1, 1],
  133. [0, 0, 0, 1],
  134. [0, 0, 0, 0]], dtype=float)
  135. bs = [
  136. np.array([1, 1, 1, 1]),
  137. np.array([1, 1, 1, 0]),
  138. np.array([1, 1, 0, 0]),
  139. np.array([1, 0, 0, 0]),
  140. ]
  141. for b in bs:
  142. with suppress_warnings() as sup:
  143. sup.filter(DeprecationWarning, ".*called without specifying.*")
  144. xp, info = lgmres(A, b, maxiter=1)
  145. resp = np.linalg.norm(A.dot(xp) - b)
  146. K = np.c_[b, A.dot(b), A.dot(A.dot(b)), A.dot(A.dot(A.dot(b)))]
  147. y, _, _, _ = np.linalg.lstsq(A.dot(K), b, rcond=-1)
  148. x = K.dot(y)
  149. res = np.linalg.norm(A.dot(x) - b)
  150. assert_allclose(resp, res, err_msg=repr(b))
  151. def test_denormals(self):
  152. # Check that no warnings are emitted if the matrix contains
  153. # numbers for which 1/x has no float representation, and that
  154. # the solver behaves properly.
  155. A = np.array([[1, 2], [3, 4]], dtype=float)
  156. A *= 100 * np.nextafter(0, 1)
  157. b = np.array([1, 1])
  158. with suppress_warnings() as sup:
  159. sup.filter(DeprecationWarning, ".*called without specifying.*")
  160. xp, info = lgmres(A, b)
  161. if info == 0:
  162. assert_allclose(A.dot(xp), b)