test_linsolve.py 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799
  1. import sys
  2. import threading
  3. import numpy as np
  4. from numpy import array, finfo, arange, eye, all, unique, ones, dot
  5. import numpy.random as random
  6. from numpy.testing import (
  7. assert_array_almost_equal, assert_almost_equal,
  8. assert_equal, assert_array_equal, assert_, assert_allclose,
  9. assert_warns, suppress_warnings)
  10. import pytest
  11. from pytest import raises as assert_raises
  12. import scipy.linalg
  13. from scipy.linalg import norm, inv
  14. from scipy.sparse import (spdiags, SparseEfficiencyWarning, csc_matrix,
  15. csr_matrix, identity, isspmatrix, dok_matrix, lil_matrix, bsr_matrix)
  16. from scipy.sparse.linalg import SuperLU
  17. from scipy.sparse.linalg._dsolve import (spsolve, use_solver, splu, spilu,
  18. MatrixRankWarning, _superlu, spsolve_triangular, factorized)
  19. import scipy.sparse
  20. from scipy._lib._testutils import check_free_memory
  21. sup_sparse_efficiency = suppress_warnings()
  22. sup_sparse_efficiency.filter(SparseEfficiencyWarning)
  23. # scikits.umfpack is not a SciPy dependency but it is optionally used in
  24. # dsolve, so check whether it's available
  25. try:
  26. import scikits.umfpack as umfpack
  27. has_umfpack = True
  28. except ImportError:
  29. has_umfpack = False
  30. def toarray(a):
  31. if isspmatrix(a):
  32. return a.toarray()
  33. else:
  34. return a
  35. def setup_bug_8278():
  36. N = 2 ** 6
  37. h = 1/N
  38. Ah1D = scipy.sparse.diags([-1, 2, -1], [-1, 0, 1],
  39. shape=(N-1, N-1))/(h**2)
  40. eyeN = scipy.sparse.eye(N - 1)
  41. A = (scipy.sparse.kron(eyeN, scipy.sparse.kron(eyeN, Ah1D))
  42. + scipy.sparse.kron(eyeN, scipy.sparse.kron(Ah1D, eyeN))
  43. + scipy.sparse.kron(Ah1D, scipy.sparse.kron(eyeN, eyeN)))
  44. b = np.random.rand((N-1)**3)
  45. return A, b
  46. class TestFactorized:
  47. def setup_method(self):
  48. n = 5
  49. d = arange(n) + 1
  50. self.n = n
  51. self.A = spdiags((d, 2*d, d[::-1]), (-3, 0, 5), n, n).tocsc()
  52. random.seed(1234)
  53. def _check_singular(self):
  54. A = csc_matrix((5,5), dtype='d')
  55. b = ones(5)
  56. assert_array_almost_equal(0. * b, factorized(A)(b))
  57. def _check_non_singular(self):
  58. # Make a diagonal dominant, to make sure it is not singular
  59. n = 5
  60. a = csc_matrix(random.rand(n, n))
  61. b = ones(n)
  62. expected = splu(a).solve(b)
  63. assert_array_almost_equal(factorized(a)(b), expected)
  64. def test_singular_without_umfpack(self):
  65. use_solver(useUmfpack=False)
  66. with assert_raises(RuntimeError, match="Factor is exactly singular"):
  67. self._check_singular()
  68. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  69. def test_singular_with_umfpack(self):
  70. use_solver(useUmfpack=True)
  71. with suppress_warnings() as sup:
  72. sup.filter(RuntimeWarning, "divide by zero encountered in double_scalars")
  73. assert_warns(umfpack.UmfpackWarning, self._check_singular)
  74. def test_non_singular_without_umfpack(self):
  75. use_solver(useUmfpack=False)
  76. self._check_non_singular()
  77. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  78. def test_non_singular_with_umfpack(self):
  79. use_solver(useUmfpack=True)
  80. self._check_non_singular()
  81. def test_cannot_factorize_nonsquare_matrix_without_umfpack(self):
  82. use_solver(useUmfpack=False)
  83. msg = "can only factor square matrices"
  84. with assert_raises(ValueError, match=msg):
  85. factorized(self.A[:, :4])
  86. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  87. def test_factorizes_nonsquare_matrix_with_umfpack(self):
  88. use_solver(useUmfpack=True)
  89. # does not raise
  90. factorized(self.A[:,:4])
  91. def test_call_with_incorrectly_sized_matrix_without_umfpack(self):
  92. use_solver(useUmfpack=False)
  93. solve = factorized(self.A)
  94. b = random.rand(4)
  95. B = random.rand(4, 3)
  96. BB = random.rand(self.n, 3, 9)
  97. with assert_raises(ValueError, match="is of incompatible size"):
  98. solve(b)
  99. with assert_raises(ValueError, match="is of incompatible size"):
  100. solve(B)
  101. with assert_raises(ValueError,
  102. match="object too deep for desired array"):
  103. solve(BB)
  104. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  105. def test_call_with_incorrectly_sized_matrix_with_umfpack(self):
  106. use_solver(useUmfpack=True)
  107. solve = factorized(self.A)
  108. b = random.rand(4)
  109. B = random.rand(4, 3)
  110. BB = random.rand(self.n, 3, 9)
  111. # does not raise
  112. solve(b)
  113. msg = "object too deep for desired array"
  114. with assert_raises(ValueError, match=msg):
  115. solve(B)
  116. with assert_raises(ValueError, match=msg):
  117. solve(BB)
  118. def test_call_with_cast_to_complex_without_umfpack(self):
  119. use_solver(useUmfpack=False)
  120. solve = factorized(self.A)
  121. b = random.rand(4)
  122. for t in [np.complex64, np.complex128]:
  123. with assert_raises(TypeError, match="Cannot cast array data"):
  124. solve(b.astype(t))
  125. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  126. def test_call_with_cast_to_complex_with_umfpack(self):
  127. use_solver(useUmfpack=True)
  128. solve = factorized(self.A)
  129. b = random.rand(4)
  130. for t in [np.complex64, np.complex128]:
  131. assert_warns(np.ComplexWarning, solve, b.astype(t))
  132. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  133. def test_assume_sorted_indices_flag(self):
  134. # a sparse matrix with unsorted indices
  135. unsorted_inds = np.array([2, 0, 1, 0])
  136. data = np.array([10, 16, 5, 0.4])
  137. indptr = np.array([0, 1, 2, 4])
  138. A = csc_matrix((data, unsorted_inds, indptr), (3, 3))
  139. b = ones(3)
  140. # should raise when incorrectly assuming indices are sorted
  141. use_solver(useUmfpack=True, assumeSortedIndices=True)
  142. with assert_raises(RuntimeError,
  143. match="UMFPACK_ERROR_invalid_matrix"):
  144. factorized(A)
  145. # should sort indices and succeed when not assuming indices are sorted
  146. use_solver(useUmfpack=True, assumeSortedIndices=False)
  147. expected = splu(A.copy()).solve(b)
  148. assert_equal(A.has_sorted_indices, 0)
  149. assert_array_almost_equal(factorized(A)(b), expected)
  150. @pytest.mark.slow
  151. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  152. def test_bug_8278(self):
  153. check_free_memory(8000)
  154. use_solver(useUmfpack=True)
  155. A, b = setup_bug_8278()
  156. A = A.tocsc()
  157. f = factorized(A)
  158. x = f(b)
  159. assert_array_almost_equal(A @ x, b)
  160. class TestLinsolve:
  161. def setup_method(self):
  162. use_solver(useUmfpack=False)
  163. def test_singular(self):
  164. A = csc_matrix((5,5), dtype='d')
  165. b = array([1, 2, 3, 4, 5],dtype='d')
  166. with suppress_warnings() as sup:
  167. sup.filter(MatrixRankWarning, "Matrix is exactly singular")
  168. x = spsolve(A, b)
  169. assert_(not np.isfinite(x).any())
  170. def test_singular_gh_3312(self):
  171. # "Bad" test case that leads SuperLU to call LAPACK with invalid
  172. # arguments. Check that it fails moderately gracefully.
  173. ij = np.array([(17, 0), (17, 6), (17, 12), (10, 13)], dtype=np.int32)
  174. v = np.array([0.284213, 0.94933781, 0.15767017, 0.38797296])
  175. A = csc_matrix((v, ij.T), shape=(20, 20))
  176. b = np.arange(20)
  177. try:
  178. # should either raise a runtime error or return value
  179. # appropriate for singular input (which yields the warning)
  180. with suppress_warnings() as sup:
  181. sup.filter(MatrixRankWarning, "Matrix is exactly singular")
  182. x = spsolve(A, b)
  183. assert not np.isfinite(x).any()
  184. except RuntimeError:
  185. pass
  186. def test_twodiags(self):
  187. A = spdiags([[1, 2, 3, 4, 5], [6, 5, 8, 9, 10]], [0, 1], 5, 5)
  188. b = array([1, 2, 3, 4, 5])
  189. # condition number of A
  190. cond_A = norm(A.toarray(), 2) * norm(inv(A.toarray()), 2)
  191. for t in ['f','d','F','D']:
  192. eps = finfo(t).eps # floating point epsilon
  193. b = b.astype(t)
  194. for format in ['csc','csr']:
  195. Asp = A.astype(t).asformat(format)
  196. x = spsolve(Asp,b)
  197. assert_(norm(b - Asp@x) < 10 * cond_A * eps)
  198. def test_bvector_smoketest(self):
  199. Adense = array([[0., 1., 1.],
  200. [1., 0., 1.],
  201. [0., 0., 1.]])
  202. As = csc_matrix(Adense)
  203. random.seed(1234)
  204. x = random.randn(3)
  205. b = As@x
  206. x2 = spsolve(As, b)
  207. assert_array_almost_equal(x, x2)
  208. def test_bmatrix_smoketest(self):
  209. Adense = array([[0., 1., 1.],
  210. [1., 0., 1.],
  211. [0., 0., 1.]])
  212. As = csc_matrix(Adense)
  213. random.seed(1234)
  214. x = random.randn(3, 4)
  215. Bdense = As.dot(x)
  216. Bs = csc_matrix(Bdense)
  217. x2 = spsolve(As, Bs)
  218. assert_array_almost_equal(x, x2.toarray())
  219. @sup_sparse_efficiency
  220. def test_non_square(self):
  221. # A is not square.
  222. A = ones((3, 4))
  223. b = ones((4, 1))
  224. assert_raises(ValueError, spsolve, A, b)
  225. # A2 and b2 have incompatible shapes.
  226. A2 = csc_matrix(eye(3))
  227. b2 = array([1.0, 2.0])
  228. assert_raises(ValueError, spsolve, A2, b2)
  229. @sup_sparse_efficiency
  230. def test_example_comparison(self):
  231. row = array([0,0,1,2,2,2])
  232. col = array([0,2,2,0,1,2])
  233. data = array([1,2,3,-4,5,6])
  234. sM = csr_matrix((data,(row,col)), shape=(3,3), dtype=float)
  235. M = sM.toarray()
  236. row = array([0,0,1,1,0,0])
  237. col = array([0,2,1,1,0,0])
  238. data = array([1,1,1,1,1,1])
  239. sN = csr_matrix((data, (row,col)), shape=(3,3), dtype=float)
  240. N = sN.toarray()
  241. sX = spsolve(sM, sN)
  242. X = scipy.linalg.solve(M, N)
  243. assert_array_almost_equal(X, sX.toarray())
  244. @sup_sparse_efficiency
  245. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  246. def test_shape_compatibility(self):
  247. use_solver(useUmfpack=True)
  248. A = csc_matrix([[1., 0], [0, 2]])
  249. bs = [
  250. [1, 6],
  251. array([1, 6]),
  252. [[1], [6]],
  253. array([[1], [6]]),
  254. csc_matrix([[1], [6]]),
  255. csr_matrix([[1], [6]]),
  256. dok_matrix([[1], [6]]),
  257. bsr_matrix([[1], [6]]),
  258. array([[1., 2., 3.], [6., 8., 10.]]),
  259. csc_matrix([[1., 2., 3.], [6., 8., 10.]]),
  260. csr_matrix([[1., 2., 3.], [6., 8., 10.]]),
  261. dok_matrix([[1., 2., 3.], [6., 8., 10.]]),
  262. bsr_matrix([[1., 2., 3.], [6., 8., 10.]]),
  263. ]
  264. for b in bs:
  265. x = np.linalg.solve(A.toarray(), toarray(b))
  266. for spmattype in [csc_matrix, csr_matrix, dok_matrix, lil_matrix]:
  267. x1 = spsolve(spmattype(A), b, use_umfpack=True)
  268. x2 = spsolve(spmattype(A), b, use_umfpack=False)
  269. # check solution
  270. if x.ndim == 2 and x.shape[1] == 1:
  271. # interprets also these as "vectors"
  272. x = x.ravel()
  273. assert_array_almost_equal(toarray(x1), x, err_msg=repr((b, spmattype, 1)))
  274. assert_array_almost_equal(toarray(x2), x, err_msg=repr((b, spmattype, 2)))
  275. # dense vs. sparse output ("vectors" are always dense)
  276. if isspmatrix(b) and x.ndim > 1:
  277. assert_(isspmatrix(x1), repr((b, spmattype, 1)))
  278. assert_(isspmatrix(x2), repr((b, spmattype, 2)))
  279. else:
  280. assert_(isinstance(x1, np.ndarray), repr((b, spmattype, 1)))
  281. assert_(isinstance(x2, np.ndarray), repr((b, spmattype, 2)))
  282. # check output shape
  283. if x.ndim == 1:
  284. # "vector"
  285. assert_equal(x1.shape, (A.shape[1],))
  286. assert_equal(x2.shape, (A.shape[1],))
  287. else:
  288. # "matrix"
  289. assert_equal(x1.shape, x.shape)
  290. assert_equal(x2.shape, x.shape)
  291. A = csc_matrix((3, 3))
  292. b = csc_matrix((1, 3))
  293. assert_raises(ValueError, spsolve, A, b)
  294. @sup_sparse_efficiency
  295. def test_ndarray_support(self):
  296. A = array([[1., 2.], [2., 0.]])
  297. x = array([[1., 1.], [0.5, -0.5]])
  298. b = array([[2., 0.], [2., 2.]])
  299. assert_array_almost_equal(x, spsolve(A, b))
  300. def test_gssv_badinput(self):
  301. N = 10
  302. d = arange(N) + 1.0
  303. A = spdiags((d, 2*d, d[::-1]), (-3, 0, 5), N, N)
  304. for spmatrix in (csc_matrix, csr_matrix):
  305. A = spmatrix(A)
  306. b = np.arange(N)
  307. def not_c_contig(x):
  308. return x.repeat(2)[::2]
  309. def not_1dim(x):
  310. return x[:,None]
  311. def bad_type(x):
  312. return x.astype(bool)
  313. def too_short(x):
  314. return x[:-1]
  315. badops = [not_c_contig, not_1dim, bad_type, too_short]
  316. for badop in badops:
  317. msg = "%r %r" % (spmatrix, badop)
  318. # Not C-contiguous
  319. assert_raises((ValueError, TypeError), _superlu.gssv,
  320. N, A.nnz, badop(A.data), A.indices, A.indptr,
  321. b, int(spmatrix == csc_matrix), err_msg=msg)
  322. assert_raises((ValueError, TypeError), _superlu.gssv,
  323. N, A.nnz, A.data, badop(A.indices), A.indptr,
  324. b, int(spmatrix == csc_matrix), err_msg=msg)
  325. assert_raises((ValueError, TypeError), _superlu.gssv,
  326. N, A.nnz, A.data, A.indices, badop(A.indptr),
  327. b, int(spmatrix == csc_matrix), err_msg=msg)
  328. def test_sparsity_preservation(self):
  329. ident = csc_matrix([
  330. [1, 0, 0],
  331. [0, 1, 0],
  332. [0, 0, 1]])
  333. b = csc_matrix([
  334. [0, 1],
  335. [1, 0],
  336. [0, 0]])
  337. x = spsolve(ident, b)
  338. assert_equal(ident.nnz, 3)
  339. assert_equal(b.nnz, 2)
  340. assert_equal(x.nnz, 2)
  341. assert_allclose(x.A, b.A, atol=1e-12, rtol=1e-12)
  342. def test_dtype_cast(self):
  343. A_real = scipy.sparse.csr_matrix([[1, 2, 0],
  344. [0, 0, 3],
  345. [4, 0, 5]])
  346. A_complex = scipy.sparse.csr_matrix([[1, 2, 0],
  347. [0, 0, 3],
  348. [4, 0, 5 + 1j]])
  349. b_real = np.array([1,1,1])
  350. b_complex = np.array([1,1,1]) + 1j*np.array([1,1,1])
  351. x = spsolve(A_real, b_real)
  352. assert_(np.issubdtype(x.dtype, np.floating))
  353. x = spsolve(A_real, b_complex)
  354. assert_(np.issubdtype(x.dtype, np.complexfloating))
  355. x = spsolve(A_complex, b_real)
  356. assert_(np.issubdtype(x.dtype, np.complexfloating))
  357. x = spsolve(A_complex, b_complex)
  358. assert_(np.issubdtype(x.dtype, np.complexfloating))
  359. @pytest.mark.slow
  360. @pytest.mark.skipif(not has_umfpack, reason="umfpack not available")
  361. def test_bug_8278(self):
  362. check_free_memory(8000)
  363. use_solver(useUmfpack=True)
  364. A, b = setup_bug_8278()
  365. x = spsolve(A, b)
  366. assert_array_almost_equal(A @ x, b)
  367. class TestSplu:
  368. def setup_method(self):
  369. use_solver(useUmfpack=False)
  370. n = 40
  371. d = arange(n) + 1
  372. self.n = n
  373. self.A = spdiags((d, 2*d, d[::-1]), (-3, 0, 5), n, n)
  374. random.seed(1234)
  375. def _smoketest(self, spxlu, check, dtype):
  376. if np.issubdtype(dtype, np.complexfloating):
  377. A = self.A + 1j*self.A.T
  378. else:
  379. A = self.A
  380. A = A.astype(dtype)
  381. lu = spxlu(A)
  382. rng = random.RandomState(1234)
  383. # Input shapes
  384. for k in [None, 1, 2, self.n, self.n+2]:
  385. msg = "k=%r" % (k,)
  386. if k is None:
  387. b = rng.rand(self.n)
  388. else:
  389. b = rng.rand(self.n, k)
  390. if np.issubdtype(dtype, np.complexfloating):
  391. b = b + 1j*rng.rand(*b.shape)
  392. b = b.astype(dtype)
  393. x = lu.solve(b)
  394. check(A, b, x, msg)
  395. x = lu.solve(b, 'T')
  396. check(A.T, b, x, msg)
  397. x = lu.solve(b, 'H')
  398. check(A.T.conj(), b, x, msg)
  399. @sup_sparse_efficiency
  400. def test_splu_smoketest(self):
  401. self._internal_test_splu_smoketest()
  402. def _internal_test_splu_smoketest(self):
  403. # Check that splu works at all
  404. def check(A, b, x, msg=""):
  405. eps = np.finfo(A.dtype).eps
  406. r = A @ x
  407. assert_(abs(r - b).max() < 1e3*eps, msg)
  408. self._smoketest(splu, check, np.float32)
  409. self._smoketest(splu, check, np.float64)
  410. self._smoketest(splu, check, np.complex64)
  411. self._smoketest(splu, check, np.complex128)
  412. @sup_sparse_efficiency
  413. def test_spilu_smoketest(self):
  414. self._internal_test_spilu_smoketest()
  415. def _internal_test_spilu_smoketest(self):
  416. errors = []
  417. def check(A, b, x, msg=""):
  418. r = A @ x
  419. err = abs(r - b).max()
  420. assert_(err < 1e-2, msg)
  421. if b.dtype in (np.float64, np.complex128):
  422. errors.append(err)
  423. self._smoketest(spilu, check, np.float32)
  424. self._smoketest(spilu, check, np.float64)
  425. self._smoketest(spilu, check, np.complex64)
  426. self._smoketest(spilu, check, np.complex128)
  427. assert_(max(errors) > 1e-5)
  428. @sup_sparse_efficiency
  429. def test_spilu_drop_rule(self):
  430. # Test passing in the drop_rule argument to spilu.
  431. A = identity(2)
  432. rules = [
  433. b'basic,area'.decode('ascii'), # unicode
  434. b'basic,area', # ascii
  435. [b'basic', b'area'.decode('ascii')]
  436. ]
  437. for rule in rules:
  438. # Argument should be accepted
  439. assert_(isinstance(spilu(A, drop_rule=rule), SuperLU))
  440. def test_splu_nnz0(self):
  441. A = csc_matrix((5,5), dtype='d')
  442. assert_raises(RuntimeError, splu, A)
  443. def test_spilu_nnz0(self):
  444. A = csc_matrix((5,5), dtype='d')
  445. assert_raises(RuntimeError, spilu, A)
  446. def test_splu_basic(self):
  447. # Test basic splu functionality.
  448. n = 30
  449. rng = random.RandomState(12)
  450. a = rng.rand(n, n)
  451. a[a < 0.95] = 0
  452. # First test with a singular matrix
  453. a[:, 0] = 0
  454. a_ = csc_matrix(a)
  455. # Matrix is exactly singular
  456. assert_raises(RuntimeError, splu, a_)
  457. # Make a diagonal dominant, to make sure it is not singular
  458. a += 4*eye(n)
  459. a_ = csc_matrix(a)
  460. lu = splu(a_)
  461. b = ones(n)
  462. x = lu.solve(b)
  463. assert_almost_equal(dot(a, x), b)
  464. def test_splu_perm(self):
  465. # Test the permutation vectors exposed by splu.
  466. n = 30
  467. a = random.random((n, n))
  468. a[a < 0.95] = 0
  469. # Make a diagonal dominant, to make sure it is not singular
  470. a += 4*eye(n)
  471. a_ = csc_matrix(a)
  472. lu = splu(a_)
  473. # Check that the permutation indices do belong to [0, n-1].
  474. for perm in (lu.perm_r, lu.perm_c):
  475. assert_(all(perm > -1))
  476. assert_(all(perm < n))
  477. assert_equal(len(unique(perm)), len(perm))
  478. # Now make a symmetric, and test that the two permutation vectors are
  479. # the same
  480. # Note: a += a.T relies on undefined behavior.
  481. a = a + a.T
  482. a_ = csc_matrix(a)
  483. lu = splu(a_)
  484. assert_array_equal(lu.perm_r, lu.perm_c)
  485. @pytest.mark.parametrize("splu_fun, rtol", [(splu, 1e-7), (spilu, 1e-1)])
  486. def test_natural_permc(self, splu_fun, rtol):
  487. # Test that the "NATURAL" permc_spec does not permute the matrix
  488. np.random.seed(42)
  489. n = 500
  490. p = 0.01
  491. A = scipy.sparse.random(n, n, p)
  492. x = np.random.rand(n)
  493. # Make A diagonal dominant to make sure it is not singular
  494. A += (n+1)*scipy.sparse.identity(n)
  495. A_ = csc_matrix(A)
  496. b = A_ @ x
  497. # without permc_spec, permutation is not identity
  498. lu = splu_fun(A_)
  499. assert_(np.any(lu.perm_c != np.arange(n)))
  500. # with permc_spec="NATURAL", permutation is identity
  501. lu = splu_fun(A_, permc_spec="NATURAL")
  502. assert_array_equal(lu.perm_c, np.arange(n))
  503. # Also, lu decomposition is valid
  504. x2 = lu.solve(b)
  505. assert_allclose(x, x2, rtol=rtol)
  506. @pytest.mark.skipif(not hasattr(sys, 'getrefcount'), reason="no sys.getrefcount")
  507. def test_lu_refcount(self):
  508. # Test that we are keeping track of the reference count with splu.
  509. n = 30
  510. a = random.random((n, n))
  511. a[a < 0.95] = 0
  512. # Make a diagonal dominant, to make sure it is not singular
  513. a += 4*eye(n)
  514. a_ = csc_matrix(a)
  515. lu = splu(a_)
  516. # And now test that we don't have a refcount bug
  517. rc = sys.getrefcount(lu)
  518. for attr in ('perm_r', 'perm_c'):
  519. perm = getattr(lu, attr)
  520. assert_equal(sys.getrefcount(lu), rc + 1)
  521. del perm
  522. assert_equal(sys.getrefcount(lu), rc)
  523. def test_bad_inputs(self):
  524. A = self.A.tocsc()
  525. assert_raises(ValueError, splu, A[:,:4])
  526. assert_raises(ValueError, spilu, A[:,:4])
  527. for lu in [splu(A), spilu(A)]:
  528. b = random.rand(42)
  529. B = random.rand(42, 3)
  530. BB = random.rand(self.n, 3, 9)
  531. assert_raises(ValueError, lu.solve, b)
  532. assert_raises(ValueError, lu.solve, B)
  533. assert_raises(ValueError, lu.solve, BB)
  534. assert_raises(TypeError, lu.solve,
  535. b.astype(np.complex64))
  536. assert_raises(TypeError, lu.solve,
  537. b.astype(np.complex128))
  538. @sup_sparse_efficiency
  539. def test_superlu_dlamch_i386_nan(self):
  540. # SuperLU 4.3 calls some functions returning floats without
  541. # declaring them. On i386@linux call convention, this fails to
  542. # clear floating point registers after call. As a result, NaN
  543. # can appear in the next floating point operation made.
  544. #
  545. # Here's a test case that triggered the issue.
  546. n = 8
  547. d = np.arange(n) + 1
  548. A = spdiags((d, 2*d, d[::-1]), (-3, 0, 5), n, n)
  549. A = A.astype(np.float32)
  550. spilu(A)
  551. A = A + 1j*A
  552. B = A.A
  553. assert_(not np.isnan(B).any())
  554. @sup_sparse_efficiency
  555. def test_lu_attr(self):
  556. def check(dtype, complex_2=False):
  557. A = self.A.astype(dtype)
  558. if complex_2:
  559. A = A + 1j*A.T
  560. n = A.shape[0]
  561. lu = splu(A)
  562. # Check that the decomposition is as advertized
  563. Pc = np.zeros((n, n))
  564. Pc[np.arange(n), lu.perm_c] = 1
  565. Pr = np.zeros((n, n))
  566. Pr[lu.perm_r, np.arange(n)] = 1
  567. Ad = A.toarray()
  568. lhs = Pr.dot(Ad).dot(Pc)
  569. rhs = (lu.L @ lu.U).toarray()
  570. eps = np.finfo(dtype).eps
  571. assert_allclose(lhs, rhs, atol=100*eps)
  572. check(np.float32)
  573. check(np.float64)
  574. check(np.complex64)
  575. check(np.complex128)
  576. check(np.complex64, True)
  577. check(np.complex128, True)
  578. @pytest.mark.slow
  579. @sup_sparse_efficiency
  580. def test_threads_parallel(self):
  581. oks = []
  582. def worker():
  583. try:
  584. self.test_splu_basic()
  585. self._internal_test_splu_smoketest()
  586. self._internal_test_spilu_smoketest()
  587. oks.append(True)
  588. except Exception:
  589. pass
  590. threads = [threading.Thread(target=worker)
  591. for k in range(20)]
  592. for t in threads:
  593. t.start()
  594. for t in threads:
  595. t.join()
  596. assert_equal(len(oks), 20)
  597. class TestSpsolveTriangular:
  598. def setup_method(self):
  599. use_solver(useUmfpack=False)
  600. def test_zero_diagonal(self):
  601. n = 5
  602. rng = np.random.default_rng(43876432987)
  603. A = rng.standard_normal((n, n))
  604. b = np.arange(n)
  605. A = scipy.sparse.tril(A, k=0, format='csr')
  606. x = spsolve_triangular(A, b, unit_diagonal=True, lower=True)
  607. A.setdiag(1)
  608. assert_allclose(A.dot(x), b)
  609. # Regression test from gh-15199
  610. A = np.array([[0, 0, 0], [1, 0, 0], [1, 1, 0]], dtype=np.float64)
  611. b = np.array([1., 2., 3.])
  612. with suppress_warnings() as sup:
  613. sup.filter(SparseEfficiencyWarning, "CSR matrix format is")
  614. spsolve_triangular(A, b, unit_diagonal=True)
  615. def test_singular(self):
  616. n = 5
  617. A = csr_matrix((n, n))
  618. b = np.arange(n)
  619. for lower in (True, False):
  620. assert_raises(scipy.linalg.LinAlgError, spsolve_triangular, A, b, lower=lower)
  621. @sup_sparse_efficiency
  622. def test_bad_shape(self):
  623. # A is not square.
  624. A = np.zeros((3, 4))
  625. b = ones((4, 1))
  626. assert_raises(ValueError, spsolve_triangular, A, b)
  627. # A2 and b2 have incompatible shapes.
  628. A2 = csr_matrix(eye(3))
  629. b2 = array([1.0, 2.0])
  630. assert_raises(ValueError, spsolve_triangular, A2, b2)
  631. @sup_sparse_efficiency
  632. def test_input_types(self):
  633. A = array([[1., 0.], [1., 2.]])
  634. b = array([[2., 0.], [2., 2.]])
  635. for matrix_type in (array, csc_matrix, csr_matrix):
  636. x = spsolve_triangular(matrix_type(A), b, lower=True)
  637. assert_array_almost_equal(A.dot(x), b)
  638. @pytest.mark.slow
  639. @pytest.mark.timeout(120) # prerelease_deps_coverage_64bit_blas job
  640. @sup_sparse_efficiency
  641. def test_random(self):
  642. def random_triangle_matrix(n, lower=True):
  643. A = scipy.sparse.random(n, n, density=0.1, format='coo')
  644. if lower:
  645. A = scipy.sparse.tril(A)
  646. else:
  647. A = scipy.sparse.triu(A)
  648. A = A.tocsr(copy=False)
  649. for i in range(n):
  650. A[i, i] = np.random.rand() + 1
  651. return A
  652. np.random.seed(1234)
  653. for lower in (True, False):
  654. for n in (10, 10**2, 10**3):
  655. A = random_triangle_matrix(n, lower=lower)
  656. for m in (1, 10):
  657. for b in (np.random.rand(n, m),
  658. np.random.randint(-9, 9, (n, m)),
  659. np.random.randint(-9, 9, (n, m)) +
  660. np.random.randint(-9, 9, (n, m)) * 1j):
  661. x = spsolve_triangular(A, b, lower=lower)
  662. assert_array_almost_equal(A.dot(x), b)
  663. x = spsolve_triangular(A, b, lower=lower,
  664. unit_diagonal=True)
  665. A.setdiag(1)
  666. assert_array_almost_equal(A.dot(x), b)