test_array_api.py 7.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339
  1. import pytest
  2. import numpy as np
  3. import numpy.testing as npt
  4. import scipy.sparse
  5. import scipy.sparse.linalg as spla
  6. sparray_types = ('bsr', 'coo', 'csc', 'csr', 'dia', 'dok', 'lil')
  7. sparray_classes = [
  8. getattr(scipy.sparse, f'{T}_array') for T in sparray_types
  9. ]
  10. A = np.array([
  11. [0, 1, 2, 0],
  12. [2, 0, 0, 3],
  13. [1, 4, 0, 0]
  14. ])
  15. B = np.array([
  16. [0, 1],
  17. [2, 0]
  18. ])
  19. X = np.array([
  20. [1, 0, 0, 1],
  21. [2, 1, 2, 0],
  22. [0, 2, 1, 0],
  23. [0, 0, 1, 2]
  24. ], dtype=float)
  25. sparrays = [sparray(A) for sparray in sparray_classes]
  26. square_sparrays = [sparray(B) for sparray in sparray_classes]
  27. eig_sparrays = [sparray(X) for sparray in sparray_classes]
  28. parametrize_sparrays = pytest.mark.parametrize(
  29. "A", sparrays, ids=sparray_types
  30. )
  31. parametrize_square_sparrays = pytest.mark.parametrize(
  32. "B", square_sparrays, ids=sparray_types
  33. )
  34. parametrize_eig_sparrays = pytest.mark.parametrize(
  35. "X", eig_sparrays, ids=sparray_types
  36. )
  37. @parametrize_sparrays
  38. def test_sum(A):
  39. assert not isinstance(A.sum(axis=0), np.matrix), \
  40. "Expected array, got matrix"
  41. assert A.sum(axis=0).shape == (4,)
  42. assert A.sum(axis=1).shape == (3,)
  43. @parametrize_sparrays
  44. def test_mean(A):
  45. assert not isinstance(A.mean(axis=1), np.matrix), \
  46. "Expected array, got matrix"
  47. @parametrize_sparrays
  48. def test_todense(A):
  49. assert not isinstance(A.todense(), np.matrix), \
  50. "Expected array, got matrix"
  51. @parametrize_sparrays
  52. def test_indexing(A):
  53. if A.__class__.__name__[:3] in ('dia', 'coo', 'bsr'):
  54. return
  55. with pytest.raises(NotImplementedError):
  56. A[1, :]
  57. with pytest.raises(NotImplementedError):
  58. A[:, 1]
  59. with pytest.raises(NotImplementedError):
  60. A[1, [1, 2]]
  61. with pytest.raises(NotImplementedError):
  62. A[[1, 2], 1]
  63. assert A[[0]]._is_array, "Expected sparse array, got sparse matrix"
  64. assert A[1, [[1, 2]]]._is_array, "Expected ndarray, got sparse array"
  65. assert A[[[1, 2]], 1]._is_array, "Expected ndarray, got sparse array"
  66. assert A[:, [1, 2]]._is_array, "Expected sparse array, got something else"
  67. @parametrize_sparrays
  68. def test_dense_addition(A):
  69. X = np.random.random(A.shape)
  70. assert not isinstance(A + X, np.matrix), "Expected array, got matrix"
  71. @parametrize_sparrays
  72. def test_sparse_addition(A):
  73. assert (A + A)._is_array, "Expected array, got matrix"
  74. @parametrize_sparrays
  75. def test_elementwise_mul(A):
  76. assert np.all((A * A).todense() == A.power(2).todense())
  77. @parametrize_sparrays
  78. def test_elementwise_rmul(A):
  79. with pytest.raises(TypeError):
  80. None * A
  81. with pytest.raises(ValueError):
  82. np.eye(3) * scipy.sparse.csr_array(np.arange(6).reshape(2, 3))
  83. assert np.all((2 * A) == (A.todense() * 2))
  84. assert np.all((A.todense() * A) == (A.todense() ** 2))
  85. @parametrize_sparrays
  86. def test_matmul(A):
  87. assert np.all((A @ A.T).todense() == A.dot(A.T).todense())
  88. @parametrize_square_sparrays
  89. def test_pow(B):
  90. assert (B**0)._is_array, "Expected array, got matrix"
  91. assert (B**2)._is_array, "Expected array, got matrix"
  92. @parametrize_sparrays
  93. def test_sparse_divide(A):
  94. assert isinstance(A / A, np.ndarray)
  95. @parametrize_sparrays
  96. def test_dense_divide(A):
  97. assert (A / 2)._is_array, "Expected array, got matrix"
  98. @parametrize_sparrays
  99. def test_no_A_attr(A):
  100. with pytest.warns(np.VisibleDeprecationWarning):
  101. A.A
  102. @parametrize_sparrays
  103. def test_no_H_attr(A):
  104. with pytest.warns(np.VisibleDeprecationWarning):
  105. A.H
  106. @parametrize_sparrays
  107. def test_getrow_getcol(A):
  108. assert A.getcol(0)._is_array
  109. assert A.getrow(0)._is_array
  110. @parametrize_sparrays
  111. def test_docstr(A):
  112. if A.__doc__ is None:
  113. return
  114. docstr = A.__doc__.lower()
  115. for phrase in ('matrix', 'matrices'):
  116. assert phrase not in docstr
  117. # -- linalg --
  118. @parametrize_sparrays
  119. def test_as_linearoperator(A):
  120. L = spla.aslinearoperator(A)
  121. npt.assert_allclose(L * [1, 2, 3, 4], A @ [1, 2, 3, 4])
  122. @parametrize_square_sparrays
  123. def test_inv(B):
  124. if B.__class__.__name__[:3] != 'csc':
  125. return
  126. C = spla.inv(B)
  127. assert C._is_array
  128. npt.assert_allclose(C.todense(), np.linalg.inv(B.todense()))
  129. @parametrize_square_sparrays
  130. def test_expm(B):
  131. if B.__class__.__name__[:3] != 'csc':
  132. return
  133. Bmat = scipy.sparse.csc_matrix(B)
  134. C = spla.expm(B)
  135. assert C._is_array
  136. npt.assert_allclose(
  137. C.todense(),
  138. spla.expm(Bmat).todense()
  139. )
  140. @parametrize_square_sparrays
  141. def test_expm_multiply(B):
  142. if B.__class__.__name__[:3] != 'csc':
  143. return
  144. npt.assert_allclose(
  145. spla.expm_multiply(B, np.array([1, 2])),
  146. spla.expm(B) @ [1, 2]
  147. )
  148. @parametrize_sparrays
  149. def test_norm(A):
  150. C = spla.norm(A)
  151. npt.assert_allclose(C, np.linalg.norm(A.todense()))
  152. @parametrize_square_sparrays
  153. def test_onenormest(B):
  154. C = spla.onenormest(B)
  155. npt.assert_allclose(C, np.linalg.norm(B.todense(), 1))
  156. @parametrize_square_sparrays
  157. def test_spsolve(B):
  158. if B.__class__.__name__[:3] not in ('csc', 'csr'):
  159. return
  160. npt.assert_allclose(
  161. spla.spsolve(B, [1, 2]),
  162. np.linalg.solve(B.todense(), [1, 2])
  163. )
  164. def test_spsolve_triangular():
  165. X = scipy.sparse.csr_array([
  166. [1, 0, 0, 0],
  167. [2, 1, 0, 0],
  168. [3, 2, 1, 0],
  169. [4, 3, 2, 1],
  170. ])
  171. spla.spsolve_triangular(X, [1, 2, 3, 4])
  172. @parametrize_square_sparrays
  173. def test_factorized(B):
  174. if B.__class__.__name__[:3] != 'csc':
  175. return
  176. LU = spla.factorized(B)
  177. npt.assert_allclose(
  178. LU(np.array([1, 2])),
  179. np.linalg.solve(B.todense(), [1, 2])
  180. )
  181. @parametrize_square_sparrays
  182. @pytest.mark.parametrize(
  183. "solver",
  184. ["bicg", "bicgstab", "cg", "cgs", "gmres", "lgmres", "minres", "qmr",
  185. "gcrotmk", "tfqmr"]
  186. )
  187. def test_solvers(B, solver):
  188. if solver == "minres":
  189. kwargs = {}
  190. else:
  191. kwargs = {'atol': 1e-5}
  192. x, info = getattr(spla, solver)(B, np.array([1, 2]), **kwargs)
  193. assert info >= 0 # no errors, even if perhaps did not converge fully
  194. npt.assert_allclose(x, [1, 1], atol=1e-1)
  195. @parametrize_sparrays
  196. @pytest.mark.parametrize(
  197. "solver",
  198. ["lsqr", "lsmr"]
  199. )
  200. def test_lstsqr(A, solver):
  201. x, *_ = getattr(spla, solver)(A, [1, 2, 3])
  202. npt.assert_allclose(A @ x, [1, 2, 3])
  203. @parametrize_eig_sparrays
  204. def test_eigs(X):
  205. e, v = spla.eigs(X, k=1)
  206. npt.assert_allclose(
  207. X @ v,
  208. e[0] * v
  209. )
  210. @parametrize_eig_sparrays
  211. def test_eigsh(X):
  212. X = X + X.T
  213. e, v = spla.eigsh(X, k=1)
  214. npt.assert_allclose(
  215. X @ v,
  216. e[0] * v
  217. )
  218. @parametrize_eig_sparrays
  219. def test_svds(X):
  220. u, s, vh = spla.svds(X, k=3)
  221. u2, s2, vh2 = np.linalg.svd(X.todense())
  222. s = np.sort(s)
  223. s2 = np.sort(s2[:3])
  224. npt.assert_allclose(s, s2, atol=1e-3)
  225. def test_splu():
  226. X = scipy.sparse.csc_array([
  227. [1, 0, 0, 0],
  228. [2, 1, 0, 0],
  229. [3, 2, 1, 0],
  230. [4, 3, 2, 1],
  231. ])
  232. LU = spla.splu(X)
  233. npt.assert_allclose(LU.solve(np.array([1, 2, 3, 4])), [1, 0, 0, 0])
  234. def test_spilu():
  235. X = scipy.sparse.csc_array([
  236. [1, 0, 0, 0],
  237. [2, 1, 0, 0],
  238. [3, 2, 1, 0],
  239. [4, 3, 2, 1],
  240. ])
  241. LU = spla.spilu(X)
  242. npt.assert_allclose(LU.solve(np.array([1, 2, 3, 4])), [1, 0, 0, 0])
  243. @parametrize_sparrays
  244. def test_power_operator(A):
  245. # https://github.com/scipy/scipy/issues/15948
  246. npt.assert_equal((A**2).todense(), (A.todense())**2)