test_least_squares.py 37 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999
  1. import warnings
  2. from itertools import product
  3. from multiprocessing import Pool
  4. import numpy as np
  5. from numpy.linalg import norm
  6. from numpy.testing import (assert_, assert_allclose,
  7. assert_equal)
  8. import pytest
  9. from pytest import raises as assert_raises
  10. from scipy.sparse import issparse, lil_array
  11. from scipy.sparse.linalg import aslinearoperator
  12. from scipy.optimize import least_squares, Bounds, minimize
  13. from scipy.optimize._lsq.least_squares import IMPLEMENTED_LOSSES
  14. from scipy.optimize._lsq.common import EPS, make_strictly_feasible, CL_scaling_vector
  15. from scipy.optimize import OptimizeResult
  16. def fun_trivial(x, a=0):
  17. return (x - a)**2 + 5.0
  18. def jac_trivial(x, a=0.0):
  19. return 2 * (x - a)
  20. def fun_2d_trivial(x):
  21. return np.array([x[0], x[1]])
  22. def jac_2d_trivial(x):
  23. return np.identity(2)
  24. def fun_rosenbrock(x):
  25. return np.array([10 * (x[1] - x[0]**2), (1 - x[0])])
  26. class Fun_Rosenbrock:
  27. def __init__(self):
  28. self.nfev = 0
  29. def __call__(self, x, a=0):
  30. self.nfev += 1
  31. return fun_rosenbrock(x)
  32. def jac_rosenbrock(x):
  33. return np.array([
  34. [-20 * x[0], 10],
  35. [-1, 0]
  36. ])
  37. def jac_rosenbrock_bad_dim(x):
  38. return np.array([
  39. [-20 * x[0], 10],
  40. [-1, 0],
  41. [0.0, 0.0]
  42. ])
  43. def fun_rosenbrock_cropped(x):
  44. return fun_rosenbrock(x)[0]
  45. def jac_rosenbrock_cropped(x):
  46. return jac_rosenbrock(x)[0]
  47. # When x is 1-D array, return is 2-D array.
  48. def fun_wrong_dimensions(x):
  49. return np.array([x, x**2, x**3])
  50. def jac_wrong_dimensions(x, a=0.0):
  51. return np.atleast_3d(jac_trivial(x, a=a))
  52. def fun_bvp(x):
  53. n = int(np.sqrt(x.shape[0]))
  54. u = np.zeros((n + 2, n + 2))
  55. x = x.reshape((n, n))
  56. u[1:-1, 1:-1] = x
  57. y = u[:-2, 1:-1] + u[2:, 1:-1] + u[1:-1, :-2] + u[1:-1, 2:] - 4 * x + x**3
  58. return y.ravel()
  59. class BroydenTridiagonal:
  60. def __init__(self, n=100, mode='sparse'):
  61. rng = np.random.default_rng(123440)
  62. self.n = n
  63. self.x0 = -np.ones(n)
  64. self.lb = np.linspace(-2, -1.5, n)
  65. self.ub = np.linspace(-0.8, 0.0, n)
  66. self.lb += 0.1 * rng.standard_normal(n)
  67. self.ub += 0.1 * rng.standard_normal(n)
  68. self.x0 += 0.1 * rng.standard_normal(n)
  69. self.x0 = make_strictly_feasible(self.x0, self.lb, self.ub)
  70. if mode == 'sparse':
  71. self.sparsity = lil_array((n, n), dtype=int)
  72. i = np.arange(n)
  73. self.sparsity[i, i] = 1
  74. i = np.arange(1, n)
  75. self.sparsity[i, i - 1] = 1
  76. i = np.arange(n - 1)
  77. self.sparsity[i, i + 1] = 1
  78. self.jac = self._jac
  79. elif mode == 'operator':
  80. self.jac = lambda x: aslinearoperator(self._jac(x))
  81. elif mode == 'dense':
  82. self.sparsity = None
  83. self.jac = lambda x: self._jac(x).toarray()
  84. else:
  85. assert_(False)
  86. def fun(self, x):
  87. f = (3 - x) * x + 1
  88. f[1:] -= x[:-1]
  89. f[:-1] -= 2 * x[1:]
  90. return f
  91. def _jac(self, x):
  92. J = lil_array((self.n, self.n))
  93. i = np.arange(self.n)
  94. J[i, i] = 3 - 2 * x
  95. i = np.arange(1, self.n)
  96. J[i, i - 1] = -1
  97. i = np.arange(self.n - 1)
  98. J[i, i + 1] = -2
  99. return J
  100. class ExponentialFittingProblem:
  101. """Provide data and function for exponential fitting in the form
  102. y = a + exp(b * x) + noise."""
  103. def __init__(self, a, b, noise, n_outliers=1, x_range=(-1, 1),
  104. n_points=11, rng=None):
  105. rng = np.random.default_rng(rng)
  106. self.m = n_points
  107. self.n = 2
  108. self.p0 = np.zeros(2)
  109. self.x = np.linspace(x_range[0], x_range[1], n_points)
  110. self.y = a + np.exp(b * self.x)
  111. self.y += noise * rng.standard_normal(self.m)
  112. outliers = rng.integers(0, self.m, n_outliers)
  113. self.y[outliers] += 50 * noise * rng.random(n_outliers)
  114. self.p_opt = np.array([a, b])
  115. def fun(self, p):
  116. return p[0] + np.exp(p[1] * self.x) - self.y
  117. def jac(self, p):
  118. J = np.empty((self.m, self.n))
  119. J[:, 0] = 1
  120. J[:, 1] = self.x * np.exp(p[1] * self.x)
  121. return J
  122. def cubic_soft_l1(z):
  123. rho = np.empty((3, z.size))
  124. t = 1 + z
  125. rho[0] = 3 * (t**(1/3) - 1)
  126. rho[1] = t ** (-2/3)
  127. rho[2] = -2/3 * t**(-5/3)
  128. return rho
  129. LOSSES = list(IMPLEMENTED_LOSSES.keys()) + [cubic_soft_l1]
  130. class BaseMixin:
  131. msg = "jac='(3-point|cs)' works equivalently to '2-point' for method='lm'"
  132. def test_basic(self):
  133. # Test that the basic calling sequence works.
  134. res = least_squares(fun_trivial, 2., method=self.method)
  135. assert_allclose(res.x, 0, atol=1e-4)
  136. assert_allclose(res.fun, fun_trivial(res.x))
  137. def test_args_kwargs(self):
  138. # Test that args and kwargs are passed correctly to the functions.
  139. a = 3.0
  140. for jac in ['2-point', '3-point', 'cs', jac_trivial]:
  141. with warnings.catch_warnings():
  142. warnings.filterwarnings("ignore", self.msg, UserWarning)
  143. res = least_squares(fun_trivial, 2.0, jac, args=(a,),
  144. method=self.method)
  145. res1 = least_squares(fun_trivial, 2.0, jac, kwargs={'a': a},
  146. method=self.method)
  147. assert_allclose(res.x, a, rtol=1e-4)
  148. assert_allclose(res1.x, a, rtol=1e-4)
  149. assert_raises(TypeError, least_squares, fun_trivial, 2.0,
  150. args=(3, 4,), method=self.method)
  151. assert_raises(TypeError, least_squares, fun_trivial, 2.0,
  152. kwargs={'kaboom': 3}, method=self.method)
  153. def test_jac_options(self):
  154. for jac in ['2-point', '3-point', 'cs', jac_trivial]:
  155. with warnings.catch_warnings():
  156. warnings.filterwarnings("ignore", self.msg, UserWarning)
  157. res = least_squares(fun_trivial, 2.0, jac, method=self.method)
  158. assert_allclose(res.x, 0, atol=1e-4)
  159. assert_raises(ValueError, least_squares, fun_trivial, 2.0, jac='oops',
  160. method=self.method)
  161. def test_nfev_options(self):
  162. for max_nfev in [None, 20]:
  163. res = least_squares(fun_trivial, 2.0, max_nfev=max_nfev,
  164. method=self.method)
  165. assert_allclose(res.x, 0, atol=1e-4)
  166. def test_x_scale_options(self):
  167. for x_scale in [1.0, np.array([0.5]), 'jac']:
  168. res = least_squares(fun_trivial, 2.0, x_scale=x_scale)
  169. assert_allclose(res.x, 0)
  170. assert_raises(ValueError, least_squares, fun_trivial,
  171. 2.0, x_scale='auto', method=self.method)
  172. assert_raises(ValueError, least_squares, fun_trivial,
  173. 2.0, x_scale=-1.0, method=self.method)
  174. assert_raises(ValueError, least_squares, fun_trivial,
  175. 2.0, x_scale=1.0+2.0j, method=self.method)
  176. def test_diff_step(self):
  177. res1 = least_squares(fun_trivial, 2.0, diff_step=1e-1,
  178. method=self.method)
  179. res3 = least_squares(fun_trivial, 2.0,
  180. diff_step=None, method=self.method)
  181. assert_allclose(res1.x, 0, atol=1e-4)
  182. assert_allclose(res3.x, 0, atol=1e-4)
  183. def test_incorrect_options_usage(self):
  184. assert_raises(TypeError, least_squares, fun_trivial, 2.0,
  185. method=self.method, options={'no_such_option': 100})
  186. assert_raises(TypeError, least_squares, fun_trivial, 2.0,
  187. method=self.method, options={'max_nfev': 100})
  188. def test_full_result(self):
  189. # MINPACK doesn't work very well with factor=100 on this problem,
  190. # thus using low 'atol'.
  191. res = least_squares(fun_trivial, 2.0, method=self.method)
  192. assert_allclose(res.x, 0, atol=1e-4)
  193. assert_allclose(res.cost, 12.5)
  194. assert_allclose(res.fun, 5)
  195. assert_allclose(res.jac, 0, atol=1e-4)
  196. assert_allclose(res.grad, 0, atol=1e-2)
  197. assert_allclose(res.optimality, 0, atol=1e-2)
  198. assert_equal(res.active_mask, 0)
  199. if self.method == 'lm':
  200. assert_(res.njev is None)
  201. else:
  202. assert_(res.nfev < 10)
  203. assert_(res.njev < 10)
  204. assert_(res.status > 0)
  205. assert_(res.success)
  206. def test_full_result_single_fev(self):
  207. # MINPACK checks the number of nfev after the iteration,
  208. # so it's hard to tell what he is going to compute.
  209. if self.method == 'lm':
  210. return
  211. res = least_squares(fun_trivial, 2.0, method=self.method,
  212. max_nfev=1)
  213. assert_equal(res.x, np.array([2]))
  214. assert_equal(res.cost, 40.5)
  215. assert_equal(res.fun, np.array([9]))
  216. assert_equal(res.jac, np.array([[4]]))
  217. assert_equal(res.grad, np.array([36]))
  218. assert_equal(res.optimality, 36)
  219. assert_equal(res.active_mask, np.array([0]))
  220. assert_equal(res.nfev, 1)
  221. assert_equal(res.njev, 1)
  222. assert_equal(res.status, 0)
  223. assert_equal(res.success, 0)
  224. def test_nfev(self):
  225. # checks that the true number of nfev are being consumed
  226. for i in range(1, 3):
  227. rng = np.random.default_rng(128908)
  228. x0 = rng.uniform(size=2) * 10
  229. ftrivial = Fun_Rosenbrock()
  230. res = least_squares(
  231. ftrivial, x0, jac=jac_rosenbrock, method=self.method, max_nfev=i
  232. )
  233. assert res.nfev == ftrivial.nfev
  234. def test_rosenbrock(self):
  235. x0 = [-2, 1]
  236. x_opt = [1, 1]
  237. for jac, x_scale, tr_solver in product(
  238. ['2-point', '3-point', 'cs', jac_rosenbrock],
  239. [1.0, np.array([1.0, 0.2]), 'jac'],
  240. ['exact', 'lsmr']):
  241. with warnings.catch_warnings():
  242. warnings.filterwarnings("ignore", self.msg, UserWarning)
  243. res = least_squares(fun_rosenbrock, x0, jac, x_scale=x_scale,
  244. tr_solver=tr_solver, method=self.method)
  245. assert_allclose(res.x, x_opt)
  246. def test_rosenbrock_cropped(self):
  247. x0 = [-2, 1]
  248. if self.method == 'lm':
  249. assert_raises(ValueError, least_squares, fun_rosenbrock_cropped,
  250. x0, method='lm')
  251. else:
  252. for jac, x_scale, tr_solver in product(
  253. ['2-point', '3-point', 'cs', jac_rosenbrock_cropped],
  254. [1.0, np.array([1.0, 0.2]), 'jac'],
  255. ['exact', 'lsmr']):
  256. res = least_squares(
  257. fun_rosenbrock_cropped, x0, jac, x_scale=x_scale,
  258. tr_solver=tr_solver, method=self.method)
  259. assert_allclose(res.cost, 0, atol=1e-14)
  260. def test_fun_wrong_dimensions(self):
  261. assert_raises(ValueError, least_squares, fun_wrong_dimensions,
  262. 2.0, method=self.method)
  263. def test_jac_wrong_dimensions(self):
  264. assert_raises(ValueError, least_squares, fun_trivial,
  265. 2.0, jac_wrong_dimensions, method=self.method)
  266. def test_fun_and_jac_inconsistent_dimensions(self):
  267. x0 = [1, 2]
  268. assert_raises(ValueError, least_squares, fun_rosenbrock, x0,
  269. jac_rosenbrock_bad_dim, method=self.method)
  270. def test_x0_multidimensional(self):
  271. x0 = np.ones(4).reshape(2, 2)
  272. assert_raises(ValueError, least_squares, fun_trivial, x0,
  273. method=self.method)
  274. def test_x0_complex_scalar(self):
  275. x0 = 2.0 + 0.0*1j
  276. assert_raises(ValueError, least_squares, fun_trivial, x0,
  277. method=self.method)
  278. def test_x0_complex_array(self):
  279. x0 = [1.0, 2.0 + 0.0*1j]
  280. assert_raises(ValueError, least_squares, fun_trivial, x0,
  281. method=self.method)
  282. def test_bvp(self):
  283. # This test was introduced with fix #5556. It turned out that
  284. # dogbox solver had a bug with trust-region radius update, which
  285. # could block its progress and create an infinite loop. And this
  286. # discrete boundary value problem is the one which triggers it.
  287. n = 10
  288. x0 = np.ones(n**2)
  289. if self.method == 'lm':
  290. max_nfev = 5000 # To account for Jacobian estimation.
  291. else:
  292. max_nfev = 100
  293. res = least_squares(fun_bvp, x0, ftol=1e-2, method=self.method,
  294. max_nfev=max_nfev)
  295. assert_(res.nfev < max_nfev)
  296. assert_(res.cost < 0.5)
  297. def test_error_raised_when_all_tolerances_below_eps(self):
  298. # Test that all 0 tolerances are not allowed.
  299. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  300. method=self.method, ftol=None, xtol=None, gtol=None)
  301. def test_convergence_with_only_one_tolerance_enabled(self):
  302. if self.method == 'lm':
  303. return # should not do test
  304. x0 = [-2, 1]
  305. x_opt = [1, 1]
  306. for ftol, xtol, gtol in [(1e-8, None, None),
  307. (None, 1e-8, None),
  308. (None, None, 1e-8)]:
  309. res = least_squares(fun_rosenbrock, x0, jac=jac_rosenbrock,
  310. ftol=ftol, gtol=gtol, xtol=xtol,
  311. method=self.method)
  312. assert_allclose(res.x, x_opt)
  313. # This test is thread safe, but it is too slow and opens
  314. # too many file descriptors to run it in parallel.
  315. @pytest.mark.parallel_threads_limit(1)
  316. @pytest.mark.fail_slow(5.0)
  317. def test_workers(self):
  318. serial = least_squares(fun_trivial, 2.0, method=self.method)
  319. reses = []
  320. for workers in [None, 2]:
  321. res = least_squares(
  322. fun_trivial, 2.0, method=self.method, workers=workers
  323. )
  324. reses.append(res)
  325. with Pool(2) as workers:
  326. res = least_squares(
  327. fun_trivial, 2.0, method=self.method, workers=workers.map
  328. )
  329. reses.append(res)
  330. for res in reses:
  331. assert res.success
  332. assert_equal(res.x, serial.x)
  333. assert_equal(res.nfev, serial.nfev)
  334. assert_equal(res.njev, serial.njev)
  335. class BoundsMixin:
  336. def test_inconsistent(self):
  337. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  338. bounds=(10.0, 0.0), method=self.method)
  339. def test_infeasible(self):
  340. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  341. bounds=(3., 4), method=self.method)
  342. def test_wrong_number(self):
  343. assert_raises(ValueError, least_squares, fun_trivial, 2.,
  344. bounds=(1., 2, 3), method=self.method)
  345. def test_inconsistent_shape(self):
  346. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  347. bounds=(1.0, [2.0, 3.0]), method=self.method)
  348. # 1-D array won't be broadcast
  349. assert_raises(ValueError, least_squares, fun_rosenbrock, [1.0, 2.0],
  350. bounds=([0.0], [3.0, 4.0]), method=self.method)
  351. def test_in_bounds(self):
  352. for jac in ['2-point', '3-point', 'cs', jac_trivial]:
  353. res = least_squares(fun_trivial, 2.0, jac=jac,
  354. bounds=(-1.0, 3.0), method=self.method)
  355. assert_allclose(res.x, 0.0, atol=1e-4)
  356. assert_equal(res.active_mask, [0])
  357. assert_(-1 <= res.x <= 3)
  358. res = least_squares(fun_trivial, 2.0, jac=jac,
  359. bounds=(0.5, 3.0), method=self.method)
  360. assert_allclose(res.x, 0.5, atol=1e-4)
  361. assert_equal(res.active_mask, [-1])
  362. assert_(0.5 <= res.x <= 3)
  363. def test_bounds_shape(self):
  364. def get_bounds_direct(lb, ub):
  365. return lb, ub
  366. def get_bounds_instances(lb, ub):
  367. return Bounds(lb, ub)
  368. for jac in ['2-point', '3-point', 'cs', jac_2d_trivial]:
  369. for bounds_func in [get_bounds_direct, get_bounds_instances]:
  370. x0 = [1.0, 1.0]
  371. res = least_squares(fun_2d_trivial, x0, jac=jac)
  372. assert_allclose(res.x, [0.0, 0.0])
  373. res = least_squares(fun_2d_trivial, x0, jac=jac,
  374. bounds=bounds_func(0.5, [2.0, 2.0]),
  375. method=self.method)
  376. assert_allclose(res.x, [0.5, 0.5])
  377. res = least_squares(fun_2d_trivial, x0, jac=jac,
  378. bounds=bounds_func([0.3, 0.2], 3.0),
  379. method=self.method)
  380. assert_allclose(res.x, [0.3, 0.2])
  381. res = least_squares(
  382. fun_2d_trivial, x0, jac=jac,
  383. bounds=bounds_func([-1, 0.5], [1.0, 3.0]),
  384. method=self.method)
  385. assert_allclose(res.x, [0.0, 0.5], atol=1e-5)
  386. def test_bounds_instances(self):
  387. res = least_squares(fun_trivial, 0.5, bounds=Bounds())
  388. assert_allclose(res.x, 0.0, atol=1e-4)
  389. res = least_squares(fun_trivial, 3.0, bounds=Bounds(lb=1.0))
  390. assert_allclose(res.x, 1.0, atol=1e-4)
  391. res = least_squares(fun_trivial, 0.5, bounds=Bounds(lb=-1.0, ub=1.0))
  392. assert_allclose(res.x, 0.0, atol=1e-4)
  393. res = least_squares(fun_trivial, -3.0, bounds=Bounds(ub=-1.0))
  394. assert_allclose(res.x, -1.0, atol=1e-4)
  395. res = least_squares(fun_2d_trivial, [0.5, 0.5],
  396. bounds=Bounds(lb=[-1.0, -1.0], ub=1.0))
  397. assert_allclose(res.x, [0.0, 0.0], atol=1e-5)
  398. res = least_squares(fun_2d_trivial, [0.5, 0.5],
  399. bounds=Bounds(lb=[0.1, 0.1]))
  400. assert_allclose(res.x, [0.1, 0.1], atol=1e-5)
  401. @pytest.mark.fail_slow(10)
  402. def test_rosenbrock_bounds(self):
  403. x0_1 = np.array([-2.0, 1.0])
  404. x0_2 = np.array([2.0, 2.0])
  405. x0_3 = np.array([-2.0, 2.0])
  406. x0_4 = np.array([0.0, 2.0])
  407. x0_5 = np.array([-1.2, 1.0])
  408. problems = [
  409. (x0_1, ([-np.inf, -1.5], np.inf)),
  410. (x0_2, ([-np.inf, 1.5], np.inf)),
  411. (x0_3, ([-np.inf, 1.5], np.inf)),
  412. (x0_4, ([-np.inf, 1.5], [1.0, np.inf])),
  413. (x0_2, ([1.0, 1.5], [3.0, 3.0])),
  414. (x0_5, ([-50.0, 0.0], [0.5, 100]))
  415. ]
  416. for x0, bounds in problems:
  417. for jac, x_scale, tr_solver in product(
  418. ['2-point', '3-point', 'cs', jac_rosenbrock],
  419. [1.0, [1.0, 0.5], 'jac'],
  420. ['exact', 'lsmr']):
  421. res = least_squares(fun_rosenbrock, x0, jac, bounds,
  422. x_scale=x_scale, tr_solver=tr_solver,
  423. method=self.method)
  424. assert_allclose(res.optimality, 0.0, atol=1e-5)
  425. class SparseMixin:
  426. def test_exact_tr_solver(self):
  427. p = BroydenTridiagonal()
  428. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  429. tr_solver='exact', method=self.method)
  430. assert_raises(ValueError, least_squares, p.fun, p.x0,
  431. tr_solver='exact', jac_sparsity=p.sparsity,
  432. method=self.method)
  433. def test_equivalence(self):
  434. sparse = BroydenTridiagonal(mode='sparse')
  435. dense = BroydenTridiagonal(mode='dense')
  436. res_sparse = least_squares(
  437. sparse.fun, sparse.x0, jac=sparse.jac,
  438. method=self.method)
  439. res_dense = least_squares(
  440. dense.fun, dense.x0, jac=sparse.jac,
  441. method=self.method)
  442. assert_equal(res_sparse.nfev, res_dense.nfev)
  443. assert_allclose(res_sparse.x, res_dense.x, atol=1e-20)
  444. assert_allclose(res_sparse.cost, 0, atol=1e-20)
  445. assert_allclose(res_dense.cost, 0, atol=1e-20)
  446. def test_tr_options(self):
  447. p = BroydenTridiagonal()
  448. res = least_squares(p.fun, p.x0, p.jac, method=self.method,
  449. tr_options={'btol': 1e-10})
  450. assert_allclose(res.cost, 0, atol=1e-20)
  451. def test_wrong_parameters(self):
  452. p = BroydenTridiagonal()
  453. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  454. tr_solver='best', method=self.method)
  455. assert_raises(TypeError, least_squares, p.fun, p.x0, p.jac,
  456. tr_solver='lsmr', tr_options={'tol': 1e-10})
  457. def test_solver_selection(self):
  458. sparse = BroydenTridiagonal(mode='sparse')
  459. dense = BroydenTridiagonal(mode='dense')
  460. res_sparse = least_squares(sparse.fun, sparse.x0, jac=sparse.jac,
  461. method=self.method)
  462. res_dense = least_squares(dense.fun, dense.x0, jac=dense.jac,
  463. method=self.method)
  464. assert_allclose(res_sparse.cost, 0, atol=1e-20)
  465. assert_allclose(res_dense.cost, 0, atol=1e-20)
  466. assert_(issparse(res_sparse.jac))
  467. assert_(isinstance(res_dense.jac, np.ndarray))
  468. def test_numerical_jac(self):
  469. p = BroydenTridiagonal()
  470. for jac in ['2-point', '3-point', 'cs']:
  471. res_dense = least_squares(p.fun, p.x0, jac, method=self.method)
  472. res_sparse = least_squares(
  473. p.fun, p.x0, jac,method=self.method,
  474. jac_sparsity=p.sparsity)
  475. assert_equal(res_dense.nfev, res_sparse.nfev)
  476. assert_allclose(res_dense.x, res_sparse.x, atol=1e-20)
  477. assert_allclose(res_dense.cost, 0, atol=1e-20)
  478. assert_allclose(res_sparse.cost, 0, atol=1e-20)
  479. @pytest.mark.fail_slow(10)
  480. def test_with_bounds(self):
  481. p = BroydenTridiagonal()
  482. for jac, jac_sparsity in product(
  483. [p.jac, '2-point', '3-point', 'cs'], [None, p.sparsity]):
  484. res_1 = least_squares(
  485. p.fun, p.x0, jac, bounds=(p.lb, np.inf),
  486. method=self.method,jac_sparsity=jac_sparsity)
  487. res_2 = least_squares(
  488. p.fun, p.x0, jac, bounds=(-np.inf, p.ub),
  489. method=self.method, jac_sparsity=jac_sparsity)
  490. res_3 = least_squares(
  491. p.fun, p.x0, jac, bounds=(p.lb, p.ub),
  492. method=self.method, jac_sparsity=jac_sparsity)
  493. assert_allclose(res_1.optimality, 0, atol=1e-10)
  494. assert_allclose(res_2.optimality, 0, atol=1e-10)
  495. assert_allclose(res_3.optimality, 0, atol=1e-10)
  496. def test_wrong_jac_sparsity(self):
  497. p = BroydenTridiagonal()
  498. sparsity = p.sparsity[:-1]
  499. assert_raises(ValueError, least_squares, p.fun, p.x0,
  500. jac_sparsity=sparsity, method=self.method)
  501. def test_linear_operator(self):
  502. p = BroydenTridiagonal(mode='operator')
  503. res = least_squares(p.fun, p.x0, p.jac, method=self.method)
  504. assert_allclose(res.cost, 0.0, atol=1e-20)
  505. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  506. method=self.method, tr_solver='exact')
  507. def test_x_scale_jac_scale(self):
  508. p = BroydenTridiagonal()
  509. res = least_squares(p.fun, p.x0, p.jac, method=self.method,
  510. x_scale='jac')
  511. assert_allclose(res.cost, 0.0, atol=1e-20)
  512. p = BroydenTridiagonal(mode='operator')
  513. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  514. method=self.method, x_scale='jac')
  515. class LossFunctionMixin:
  516. def test_options(self):
  517. for loss in LOSSES:
  518. res = least_squares(fun_trivial, 2.0, loss=loss,
  519. method=self.method)
  520. assert_allclose(res.x, 0, atol=1e-15)
  521. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  522. loss='hinge', method=self.method)
  523. def test_fun(self):
  524. # Test that res.fun is actual residuals, and not modified by loss
  525. # function stuff.
  526. for loss in LOSSES:
  527. res = least_squares(fun_trivial, 2.0, loss=loss,
  528. method=self.method)
  529. assert_equal(res.fun, fun_trivial(res.x))
  530. def test_grad(self):
  531. # Test that res.grad is true gradient of loss function at the
  532. # solution. Use max_nfev = 1, to avoid reaching minimum.
  533. x = np.array([2.0]) # res.x will be this.
  534. res = least_squares(fun_trivial, x, jac_trivial, loss='linear',
  535. max_nfev=1, method=self.method)
  536. assert_equal(res.grad, 2 * x * (x**2 + 5))
  537. res = least_squares(fun_trivial, x, jac_trivial, loss='huber',
  538. max_nfev=1, method=self.method)
  539. assert_equal(res.grad, 2 * x)
  540. res = least_squares(fun_trivial, x, jac_trivial, loss='soft_l1',
  541. max_nfev=1, method=self.method)
  542. assert_allclose(res.grad,
  543. 2 * x * (x**2 + 5) / (1 + (x**2 + 5)**2)**0.5)
  544. res = least_squares(fun_trivial, x, jac_trivial, loss='cauchy',
  545. max_nfev=1, method=self.method)
  546. assert_allclose(res.grad, 2 * x * (x**2 + 5) / (1 + (x**2 + 5)**2))
  547. res = least_squares(fun_trivial, x, jac_trivial, loss='arctan',
  548. max_nfev=1, method=self.method)
  549. assert_allclose(res.grad, 2 * x * (x**2 + 5) / (1 + (x**2 + 5)**4))
  550. res = least_squares(fun_trivial, x, jac_trivial, loss=cubic_soft_l1,
  551. max_nfev=1, method=self.method)
  552. assert_allclose(res.grad,
  553. 2 * x * (x**2 + 5) / (1 + (x**2 + 5)**2)**(2/3))
  554. def test_jac(self):
  555. # Test that res.jac.T.dot(res.jac) gives Gauss-Newton approximation
  556. # of Hessian. This approximation is computed by doubly differentiating
  557. # the cost function and dropping the part containing second derivative
  558. # of f. For a scalar function it is computed as
  559. # H = (rho' + 2 * rho'' * f**2) * f'**2, if the expression inside the
  560. # brackets is less than EPS it is replaced by EPS. Here, we check
  561. # against the root of H.
  562. x = 2.0 # res.x will be this.
  563. f = x**2 + 5 # res.fun will be this.
  564. res = least_squares(fun_trivial, x, jac_trivial, loss='linear',
  565. max_nfev=1, method=self.method)
  566. assert_equal(res.jac, 2 * x)
  567. # For `huber` loss the Jacobian correction is identically zero
  568. # in outlier region, in such cases it is modified to be equal EPS**0.5.
  569. res = least_squares(fun_trivial, x, jac_trivial, loss='huber',
  570. max_nfev=1, method=self.method)
  571. assert_equal(res.jac, 2 * x * EPS**0.5)
  572. # Now, let's apply `loss_scale` to turn the residual into an inlier.
  573. # The loss function becomes linear.
  574. res = least_squares(fun_trivial, x, jac_trivial, loss='huber',
  575. f_scale=10, max_nfev=1)
  576. assert_equal(res.jac, 2 * x)
  577. # 'soft_l1' always gives a positive scaling.
  578. res = least_squares(fun_trivial, x, jac_trivial, loss='soft_l1',
  579. max_nfev=1, method=self.method)
  580. assert_allclose(res.jac, 2 * x * (1 + f**2)**-0.75)
  581. # For 'cauchy' the correction term turns out to be negative, and it
  582. # replaced by EPS**0.5.
  583. res = least_squares(fun_trivial, x, jac_trivial, loss='cauchy',
  584. max_nfev=1, method=self.method)
  585. assert_allclose(res.jac, 2 * x * EPS**0.5)
  586. # Now use scaling to turn the residual to inlier.
  587. res = least_squares(fun_trivial, x, jac_trivial, loss='cauchy',
  588. f_scale=10, max_nfev=1, method=self.method)
  589. fs = f / 10
  590. assert_allclose(res.jac, 2 * x * (1 - fs**2)**0.5 / (1 + fs**2))
  591. # 'arctan' gives an outlier.
  592. res = least_squares(fun_trivial, x, jac_trivial, loss='arctan',
  593. max_nfev=1, method=self.method)
  594. assert_allclose(res.jac, 2 * x * EPS**0.5)
  595. # Turn to inlier.
  596. res = least_squares(fun_trivial, x, jac_trivial, loss='arctan',
  597. f_scale=20.0, max_nfev=1, method=self.method)
  598. fs = f / 20
  599. assert_allclose(res.jac, 2 * x * (1 - 3 * fs**4)**0.5 / (1 + fs**4))
  600. # cubic_soft_l1 will give an outlier.
  601. res = least_squares(fun_trivial, x, jac_trivial, loss=cubic_soft_l1,
  602. max_nfev=1)
  603. assert_allclose(res.jac, 2 * x * EPS**0.5)
  604. # Turn to inlier.
  605. res = least_squares(fun_trivial, x, jac_trivial,
  606. loss=cubic_soft_l1, f_scale=6, max_nfev=1)
  607. fs = f / 6
  608. assert_allclose(res.jac,
  609. 2 * x * (1 - fs**2 / 3)**0.5 * (1 + fs**2)**(-5/6))
  610. def test_robustness(self):
  611. for noise in [0.1, 1.0]:
  612. p = ExponentialFittingProblem(1, 0.1, noise, rng=12220903)
  613. for jac in ['2-point', '3-point', 'cs', p.jac]:
  614. res_lsq = least_squares(p.fun, p.p0, jac=jac,
  615. method=self.method)
  616. assert_allclose(res_lsq.optimality, 0, atol=1e-2)
  617. for loss in LOSSES:
  618. if loss == 'linear':
  619. continue
  620. res_robust = least_squares(
  621. p.fun, p.p0, jac=jac, loss=loss, f_scale=noise,
  622. method=self.method)
  623. assert_allclose(res_robust.optimality, 0, atol=1e-2)
  624. assert_(norm(res_robust.x - p.p_opt) <
  625. norm(res_lsq.x - p.p_opt))
  626. class TestDogbox(BaseMixin, BoundsMixin, SparseMixin, LossFunctionMixin):
  627. method = 'dogbox'
  628. class TestTRF(BaseMixin, BoundsMixin, SparseMixin, LossFunctionMixin):
  629. method = 'trf'
  630. def test_lsmr_regularization(self):
  631. p = BroydenTridiagonal()
  632. for regularize in [True, False]:
  633. res = least_squares(p.fun, p.x0, p.jac, method='trf',
  634. tr_options={'regularize': regularize})
  635. assert_allclose(res.cost, 0, atol=1e-20)
  636. class TestLM(BaseMixin):
  637. method = 'lm'
  638. def test_bounds_not_supported(self):
  639. assert_raises(ValueError, least_squares, fun_trivial,
  640. 2.0, bounds=(-3.0, 3.0), method='lm')
  641. def test_m_less_n_not_supported(self):
  642. x0 = [-2, 1]
  643. assert_raises(ValueError, least_squares, fun_rosenbrock_cropped, x0,
  644. method='lm')
  645. def test_sparse_not_supported(self):
  646. p = BroydenTridiagonal()
  647. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  648. method='lm')
  649. def test_jac_sparsity_not_supported(self):
  650. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  651. jac_sparsity=[1], method='lm')
  652. def test_LinearOperator_not_supported(self):
  653. p = BroydenTridiagonal(mode="operator")
  654. assert_raises(ValueError, least_squares, p.fun, p.x0, p.jac,
  655. method='lm')
  656. def test_loss(self):
  657. res = least_squares(fun_trivial, 2.0, loss='linear', method='lm')
  658. assert_allclose(res.x, 0.0, atol=1e-4)
  659. assert_raises(ValueError, least_squares, fun_trivial, 2.0,
  660. method='lm', loss='huber')
  661. def test_callback_with_lm_method(self):
  662. def callback(x):
  663. assert(False) # Dummy callback function
  664. with warnings.catch_warnings():
  665. warnings.filterwarnings(
  666. "ignore",
  667. "Callback function specified, but not supported with `lm` method.",
  668. UserWarning,
  669. )
  670. least_squares(fun_trivial, x0=[0], method='lm', callback=callback)
  671. def test_basic():
  672. # test that 'method' arg is really optional
  673. res = least_squares(fun_trivial, 2.0)
  674. assert_allclose(res.x, 0, atol=1e-10)
  675. def test_callback():
  676. # test that callback function works as expected
  677. results = []
  678. def my_callback_optimresult(intermediate_result: OptimizeResult):
  679. results.append(intermediate_result)
  680. def my_callback_x(x):
  681. r = OptimizeResult()
  682. r.nit = 1
  683. r.x = x
  684. results.append(r)
  685. return False
  686. def my_callback_optimresult_stop_exception(
  687. intermediate_result: OptimizeResult):
  688. results.append(intermediate_result)
  689. raise StopIteration
  690. def my_callback_x_stop_exception(x):
  691. r = OptimizeResult()
  692. r.nit = 1
  693. r.x = x
  694. results.append(r)
  695. raise StopIteration
  696. # Try for different function signatures and stop methods
  697. callbacks_nostop = [my_callback_optimresult, my_callback_x]
  698. callbacks_stop = [my_callback_optimresult_stop_exception,
  699. my_callback_x_stop_exception]
  700. # Try for all the implemented methods: trf, trf_bounds and dogbox
  701. calls = [
  702. lambda callback: least_squares(fun_trivial, 5.0, method='trf',
  703. callback=callback),
  704. lambda callback: least_squares(fun_trivial, 5.0, method='trf',
  705. bounds=(-8.0, 8.0), callback=callback),
  706. lambda callback: least_squares(fun_trivial, 5.0, method='dogbox',
  707. callback=callback)
  708. ]
  709. for mycallback, call in product(callbacks_nostop, calls):
  710. results.clear()
  711. # Call the different implemented methods
  712. res = call(mycallback)
  713. # Check that callback was called
  714. assert len(results) > 0
  715. # Check that results data makes sense
  716. assert results[-1].nit > 0
  717. # Check that it didn't stop because of the callback
  718. assert res.status != -2
  719. # final callback x should be same as final result
  720. assert_allclose(results[-1].x, res.x)
  721. for mycallback, call in product(callbacks_stop, calls):
  722. results.clear()
  723. # Call the different implemented methods
  724. res = call(mycallback)
  725. # Check that callback was called
  726. assert len(results) > 0
  727. # Check that only one iteration was run
  728. assert results[-1].nit == 1
  729. # Check that it stopped because of the callback
  730. assert res.status == -2
  731. def test_small_tolerances_for_lm():
  732. for ftol, xtol, gtol in [(None, 1e-13, 1e-13),
  733. (1e-13, None, 1e-13),
  734. (1e-13, 1e-13, None)]:
  735. assert_raises(ValueError, least_squares, fun_trivial, 2.0, xtol=xtol,
  736. ftol=ftol, gtol=gtol, method='lm')
  737. def test_fp32_gh12991():
  738. # checks that smaller FP sizes can be used in least_squares
  739. # this is the minimum working example reported for gh12991
  740. rng = np.random.default_rng(1978)
  741. x = np.linspace(0, 1, 100, dtype=np.float32)
  742. y = rng.random(size=100, dtype=np.float32)
  743. # changed in gh21872. These functions should've been working in fp32 to force
  744. # approx_derivative to work in fp32. One of the initial steps in least_squares
  745. # is to force x0 (p) to be a float, meaning that the output of func and err would
  746. # be in float64, unless forced to be in float32
  747. def func(p, x):
  748. return (p[0] + p[1] * x).astype(np.float32)
  749. def err(p, x, y):
  750. return (func(p, x) - y).astype(np.float32)
  751. def mse(p, x, y):
  752. return np.sum(err(p, x, y)**2)
  753. res = least_squares(err, [-1.0, -1.0], args=(x, y))
  754. # previously the initial jacobian calculated for this would be all 0
  755. # and the minimize would terminate immediately, with nfev=1, would
  756. # report a successful minimization (it shouldn't have done), but be
  757. # unchanged from the initial solution.
  758. # It was terminating early because the underlying approx_derivative
  759. # used a step size for FP64 when the working space was FP32.
  760. assert res.nfev > 2
  761. # compare output to solver that doesn't use derivatives
  762. res2 = minimize(
  763. mse,
  764. [-1.0, 1.0],
  765. method='cobyqa',
  766. args=(x, y),
  767. options={'final_tr_radius': 1e-6}
  768. )
  769. assert_allclose(res.x, res2.x, atol=9e-5)
  770. def test_gh_18793_and_19351():
  771. answer = 1e-12
  772. initial_guess = 1.1e-12
  773. def chi2(x):
  774. return (x-answer)**2
  775. gtol = 1e-15
  776. res = least_squares(chi2, x0=initial_guess, gtol=1e-15, bounds=(0, np.inf))
  777. # Original motivation: gh-18793
  778. # if we choose an initial condition that is close to the solution
  779. # we shouldn't return an answer that is further away from the solution
  780. # Update: gh-19351
  781. # However this requirement does not go well with 'trf' algorithm logic.
  782. # Some regressions were reported after the presumed fix.
  783. # The returned solution is good as long as it satisfies the convergence
  784. # conditions.
  785. # Specifically in this case the scaled gradient will be sufficiently low.
  786. scaling, _ = CL_scaling_vector(res.x, res.grad,
  787. np.atleast_1d(0), np.atleast_1d(np.inf))
  788. assert res.status == 1 # Converged by gradient
  789. assert np.linalg.norm(res.grad * scaling, ord=np.inf) < gtol
  790. def test_gh_19103():
  791. # Checks that least_squares trf method selects a strictly feasible point,
  792. # and thus succeeds instead of failing,
  793. # when the initial guess is reported exactly at a boundary point.
  794. # This is a reduced example from gh191303
  795. ydata = np.array([0.] * 66 + [
  796. 1., 0., 0., 0., 0., 0., 1., 1., 0., 0., 1.,
  797. 1., 1., 1., 0., 0., 0., 1., 0., 0., 2., 1.,
  798. 0., 3., 1., 6., 5., 0., 0., 2., 8., 4., 4.,
  799. 6., 9., 7., 2., 7., 8., 2., 13., 9., 8., 11.,
  800. 10., 13., 14., 19., 11., 15., 18., 26., 19., 32., 29.,
  801. 28., 36., 32., 35., 36., 43., 52., 32., 58., 56., 52.,
  802. 67., 53., 72., 88., 77., 95., 94., 84., 86., 101., 107.,
  803. 108., 118., 96., 115., 138., 137.,
  804. ])
  805. xdata = np.arange(0, ydata.size) * 0.1
  806. def exponential_wrapped(params):
  807. A, B, x0 = params
  808. return A * np.exp(B * (xdata - x0)) - ydata
  809. x0 = [0.01, 1., 5.]
  810. bounds = ((0.01, 0, 0), (np.inf, 10, 20.9))
  811. res = least_squares(exponential_wrapped, x0, method='trf', bounds=bounds)
  812. assert res.success