test_linsolve.py 24 KB

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