test_lsq_linear.py 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290
  1. import warnings
  2. import pytest
  3. import numpy as np
  4. from numpy.linalg import lstsq
  5. from numpy.testing import assert_allclose, assert_equal, assert_
  6. from scipy.sparse import random_array, coo_array
  7. from scipy.sparse.linalg import aslinearoperator
  8. from scipy.optimize import lsq_linear
  9. from scipy.optimize._minimize import Bounds
  10. A = np.array([
  11. [0.171, -0.057],
  12. [-0.049, -0.248],
  13. [-0.166, 0.054],
  14. ])
  15. b = np.array([0.074, 1.014, -0.383])
  16. class BaseMixin:
  17. def setup_method(self):
  18. self.rnd = np.random.default_rng(1092892)
  19. def test_dense_no_bounds(self):
  20. for lsq_solver in self.lsq_solvers:
  21. res = lsq_linear(A, b, method=self.method, lsq_solver=lsq_solver)
  22. assert_allclose(res.x, lstsq(A, b, rcond=-1)[0])
  23. assert_allclose(res.x, res.unbounded_sol[0])
  24. def test_dense_bounds(self):
  25. # Solutions for comparison are taken from MATLAB.
  26. lb = np.array([-1, -10])
  27. ub = np.array([1, 0])
  28. unbounded_sol = lstsq(A, b, rcond=-1)[0]
  29. for lsq_solver in self.lsq_solvers:
  30. res = lsq_linear(A, b, (lb, ub), method=self.method,
  31. lsq_solver=lsq_solver)
  32. assert_allclose(res.x, lstsq(A, b, rcond=-1)[0])
  33. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  34. lb = np.array([0.0, -np.inf])
  35. for lsq_solver in self.lsq_solvers:
  36. res = lsq_linear(A, b, (lb, np.inf), method=self.method,
  37. lsq_solver=lsq_solver)
  38. assert_allclose(res.x, np.array([0.0, -4.084174437334673]),
  39. atol=1e-6)
  40. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  41. lb = np.array([-1, 0])
  42. for lsq_solver in self.lsq_solvers:
  43. res = lsq_linear(A, b, (lb, np.inf), method=self.method,
  44. lsq_solver=lsq_solver)
  45. assert_allclose(res.x, np.array([0.448427311733504, 0]),
  46. atol=1e-15)
  47. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  48. ub = np.array([np.inf, -5])
  49. for lsq_solver in self.lsq_solvers:
  50. res = lsq_linear(A, b, (-np.inf, ub), method=self.method,
  51. lsq_solver=lsq_solver)
  52. assert_allclose(res.x, np.array([-0.105560998682388, -5]))
  53. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  54. ub = np.array([-1, np.inf])
  55. for lsq_solver in self.lsq_solvers:
  56. res = lsq_linear(A, b, (-np.inf, ub), method=self.method,
  57. lsq_solver=lsq_solver)
  58. assert_allclose(res.x, np.array([-1, -4.181102129483254]))
  59. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  60. lb = np.array([0, -4])
  61. ub = np.array([1, 0])
  62. for lsq_solver in self.lsq_solvers:
  63. res = lsq_linear(A, b, (lb, ub), method=self.method,
  64. lsq_solver=lsq_solver)
  65. assert_allclose(res.x, np.array([0.005236663400791, -4]))
  66. assert_allclose(res.unbounded_sol[0], unbounded_sol)
  67. def test_bounds_variants(self):
  68. x = np.array([1, 3])
  69. A = self.rnd.uniform(size=(2, 2))
  70. b = A@x
  71. lb = np.array([1, 1])
  72. ub = np.array([2, 2])
  73. bounds_old = (lb, ub)
  74. bounds_new = Bounds(lb, ub)
  75. res_old = lsq_linear(A, b, bounds_old)
  76. res_new = lsq_linear(A, b, bounds_new)
  77. assert not np.allclose(res_new.x, res_new.unbounded_sol[0])
  78. assert_allclose(res_old.x, res_new.x)
  79. def test_np_matrix(self):
  80. # gh-10711
  81. with warnings.catch_warnings():
  82. warnings.simplefilter("ignore", PendingDeprecationWarning)
  83. A = np.matrix([[20, -4, 0, 2, 3], [10, -2, 1, 0, -1]])
  84. k = np.array([20, 15])
  85. lsq_linear(A, k)
  86. def test_dense_rank_deficient(self):
  87. A = np.array([[-0.307, -0.184]])
  88. b = np.array([0.773])
  89. lb = [-0.1, -0.1]
  90. ub = [0.1, 0.1]
  91. for lsq_solver in self.lsq_solvers:
  92. res = lsq_linear(A, b, (lb, ub), method=self.method,
  93. lsq_solver=lsq_solver)
  94. assert_allclose(res.x, [-0.1, -0.1])
  95. assert_allclose(res.unbounded_sol[0], lstsq(A, b, rcond=-1)[0])
  96. A = np.array([
  97. [0.334, 0.668],
  98. [-0.516, -1.032],
  99. [0.192, 0.384],
  100. ])
  101. b = np.array([-1.436, 0.135, 0.909])
  102. lb = [0, -1]
  103. ub = [1, -0.5]
  104. for lsq_solver in self.lsq_solvers:
  105. res = lsq_linear(A, b, (lb, ub), method=self.method,
  106. lsq_solver=lsq_solver)
  107. assert_allclose(res.optimality, 0, atol=1e-11)
  108. assert_allclose(res.unbounded_sol[0], lstsq(A, b, rcond=-1)[0])
  109. def test_full_result(self):
  110. lb = np.array([0, -4])
  111. ub = np.array([1, 0])
  112. res = lsq_linear(A, b, (lb, ub), method=self.method)
  113. assert_allclose(res.x, [0.005236663400791, -4])
  114. assert_allclose(res.unbounded_sol[0], lstsq(A, b, rcond=-1)[0])
  115. r = A.dot(res.x) - b
  116. assert_allclose(res.cost, 0.5 * np.dot(r, r))
  117. assert_allclose(res.fun, r)
  118. assert_allclose(res.optimality, 0.0, atol=1e-12)
  119. assert_equal(res.active_mask, [0, -1])
  120. assert_(res.nit < 15)
  121. assert_(res.status == 1 or res.status == 3)
  122. assert_(isinstance(res.message, str))
  123. assert_(res.success)
  124. # This is a test for issue #9982.
  125. def test_almost_singular(self):
  126. A = np.array(
  127. [[0.8854232310355122, 0.0365312146937765, 0.0365312146836789],
  128. [0.3742460132129041, 0.0130523214078376, 0.0130523214077873],
  129. [0.9680633871281361, 0.0319366128718639, 0.0319366128718388]])
  130. b = np.array(
  131. [0.0055029366538097, 0.0026677442422208, 0.0066612514782381])
  132. result = lsq_linear(A, b, method=self.method)
  133. assert_(result.cost < 1.1e-8)
  134. @pytest.mark.xslow
  135. def test_large_rank_deficient(self):
  136. rng = np.random.default_rng(1290209)
  137. n, m = np.sort(rng.integers(2, 1000, size=2))
  138. m *= 2 # make m >> n
  139. A = 1.0 * rng.integers(-99, 99, size=[m, n])
  140. b = 1.0 * rng.integers(-99, 99, size=[m])
  141. bounds = 1.0 * np.sort(rng.integers(-99, 99, size=(2, n)), axis=0)
  142. bounds[1, :] += 1.0 # ensure up > lb
  143. # Make the A matrix strongly rank deficient by replicating some columns
  144. w = rng.choice(n, n) # Select random columns with duplicates
  145. A = A[:, w]
  146. x_bvls = lsq_linear(A, b, bounds=bounds, method='bvls').x
  147. x_trf = lsq_linear(A, b, bounds=bounds, method='trf').x
  148. cost_bvls = np.sum((A @ x_bvls - b)**2)
  149. cost_trf = np.sum((A @ x_trf - b)**2)
  150. assert_(abs(cost_bvls - cost_trf) < cost_trf*1e-10)
  151. def test_convergence_small_array(self):
  152. A = np.array([[49.0, 41.0, -32.0],
  153. [-19.0, -32.0, -8.0],
  154. [-13.0, 10.0, 69.0]])
  155. b = np.array([-41.0, -90.0, 47.0])
  156. bounds = np.array([[31.0, -44.0, 26.0],
  157. [54.0, -32.0, 28.0]])
  158. x_bvls = lsq_linear(A, b, bounds=bounds, method='bvls').x
  159. x_trf = lsq_linear(A, b, bounds=bounds, method='trf').x
  160. cost_bvls = np.sum((A @ x_bvls - b)**2)
  161. cost_trf = np.sum((A @ x_trf - b)**2)
  162. assert_(abs(cost_bvls - cost_trf) < cost_trf*1e-10)
  163. class SparseMixin:
  164. def test_sparse_and_LinearOperator(self):
  165. m = 5000
  166. n = 1000
  167. rng = np.random.RandomState(0)
  168. A = random_array((m, n), random_state=rng)
  169. b = rng.randn(m)
  170. res = lsq_linear(A, b)
  171. assert_allclose(res.optimality, 0, atol=1e-6)
  172. A = aslinearoperator(A)
  173. res = lsq_linear(A, b)
  174. assert_allclose(res.optimality, 0, atol=1e-6)
  175. @pytest.mark.fail_slow(10)
  176. def test_sparse_bounds(self):
  177. m = 5000
  178. n = 1000
  179. rng = np.random.RandomState(0)
  180. A = random_array((m, n), random_state=rng)
  181. b = rng.randn(m)
  182. lb = rng.randn(n)
  183. ub = lb + 1
  184. res = lsq_linear(A, b, (lb, ub))
  185. assert_allclose(res.optimality, 0.0, atol=1e-6)
  186. res = lsq_linear(A, b, (lb, ub), lsmr_tol=1e-13,
  187. lsmr_maxiter=1500)
  188. assert_allclose(res.optimality, 0.0, atol=1e-6)
  189. res = lsq_linear(A, b, (lb, ub), lsmr_tol='auto')
  190. assert_allclose(res.optimality, 0.0, atol=1e-6)
  191. def test_sparse_ill_conditioned(self):
  192. # Sparse matrix with condition number of ~4 million
  193. data = np.array([1., 1., 1., 1. + 1e-6, 1.])
  194. row = np.array([0, 0, 1, 2, 2])
  195. col = np.array([0, 2, 1, 0, 2])
  196. A = coo_array((data, (row, col)), shape=(3, 3))
  197. # Get the exact solution
  198. exact_sol = lsq_linear(A.toarray(), b, lsq_solver='exact')
  199. # Default lsmr arguments should not fully converge the solution
  200. default_lsmr_sol = lsq_linear(A, b, lsq_solver='lsmr')
  201. with pytest.raises(AssertionError):
  202. assert_allclose(exact_sol.x, default_lsmr_sol.x)
  203. # By increasing the maximum lsmr iters, it will converge
  204. conv_lsmr = lsq_linear(A, b, lsq_solver='lsmr', lsmr_maxiter=10)
  205. assert_allclose(exact_sol.x, conv_lsmr.x)
  206. class TestTRF(BaseMixin, SparseMixin):
  207. method = 'trf'
  208. lsq_solvers = ['exact', 'lsmr']
  209. class TestBVLS(BaseMixin):
  210. method = 'bvls'
  211. lsq_solvers = ['exact']
  212. class TestErrorChecking:
  213. def test_option_lsmr_tol(self):
  214. # Should work with a positive float, string equal to 'auto', or None
  215. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol=1e-2)
  216. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol='auto')
  217. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol=None)
  218. # Should raise error with negative float, strings
  219. # other than 'auto', and integers
  220. err_message = "`lsmr_tol` must be None, 'auto', or positive float."
  221. with pytest.raises(ValueError, match=err_message):
  222. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol=-0.1)
  223. with pytest.raises(ValueError, match=err_message):
  224. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol='foo')
  225. with pytest.raises(ValueError, match=err_message):
  226. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_tol=1)
  227. def test_option_lsmr_maxiter(self):
  228. # Should work with positive integers or None
  229. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_maxiter=1)
  230. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_maxiter=None)
  231. # Should raise error with 0 or negative max iter
  232. err_message = "`lsmr_maxiter` must be None or positive integer."
  233. with pytest.raises(ValueError, match=err_message):
  234. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_maxiter=0)
  235. with pytest.raises(ValueError, match=err_message):
  236. _ = lsq_linear(A, b, lsq_solver='lsmr', lsmr_maxiter=-1)