least_squares.py 42 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044
  1. """Generic interface for least-squares minimization."""
  2. from warnings import warn
  3. import numpy as np
  4. from numpy.linalg import norm
  5. from scipy.sparse.linalg import LinearOperator
  6. from scipy.optimize import _minpack, OptimizeResult
  7. from scipy.optimize._differentiable_functions import VectorFunction
  8. from scipy.optimize._numdiff import group_columns
  9. from scipy.optimize._minimize import Bounds
  10. from scipy._lib._sparse import issparse
  11. from scipy._lib._array_api import array_namespace
  12. from scipy._lib._util import _workers_wrapper
  13. from .trf import trf
  14. from .dogbox import dogbox
  15. from .common import EPS, in_bounds, make_strictly_feasible
  16. from scipy.optimize._optimize import _wrap_callback
  17. TERMINATION_MESSAGES = {
  18. -2: "Stopped because `callback` function raised `StopIteration` or returned `True`",
  19. -1: "Improper input parameters status returned from `leastsq`",
  20. 0: "The maximum number of function evaluations is exceeded.",
  21. 1: "`gtol` termination condition is satisfied.",
  22. 2: "`ftol` termination condition is satisfied.",
  23. 3: "`xtol` termination condition is satisfied.",
  24. 4: "Both `ftol` and `xtol` termination conditions are satisfied."
  25. }
  26. FROM_MINPACK_TO_COMMON = {
  27. 0: -1, # Improper input parameters from MINPACK.
  28. 1: 2,
  29. 2: 3,
  30. 3: 4,
  31. 4: 1,
  32. 5: 0
  33. # There are 6, 7, 8 for too small tolerance parameters,
  34. # but we guard against it by checking ftol, xtol, gtol beforehand.
  35. }
  36. def call_minpack(fun, x0, jac, ftol, xtol, gtol, max_nfev, x_scale, jac_method=None):
  37. n = x0.size
  38. # Compute MINPACK's `diag`, which is inverse of our `x_scale` and
  39. # ``x_scale='jac'`` corresponds to ``diag=None``.
  40. # 1.16.0 - default x_scale changed to 'jac', with diag=None
  41. if x_scale is None or (isinstance(x_scale, str) and x_scale == 'jac'):
  42. diag = None
  43. else:
  44. # x_scale specified, so use that
  45. diag = 1 / x_scale
  46. full_output = True
  47. col_deriv = False
  48. factor = 100.0
  49. if max_nfev is None:
  50. max_nfev = 100 * n
  51. # lmder is typically used for systems with analytic jacobians, with lmdif being
  52. # used if there is only an objective fun (lmdif uses finite differences to estimate
  53. # jacobian). Otherwise they're very similar internally.
  54. # We now do all the finite differencing in VectorFunction, which means we can drop
  55. # lmdif and just use lmder.
  56. # for sending a copy of x0 into _lmder
  57. xp = array_namespace(x0)
  58. x, info, status = _minpack._lmder(
  59. fun, jac, xp.astype(x0, x0.dtype), (), full_output, col_deriv,
  60. ftol, xtol, gtol, max_nfev, factor, diag)
  61. f = info['fvec']
  62. J = jac(x)
  63. cost = 0.5 * np.dot(f, f)
  64. g = J.T.dot(f)
  65. g_norm = norm(g, ord=np.inf)
  66. nfev = info['nfev']
  67. if callable(jac_method):
  68. # user supplied a callable ("analytic") jac
  69. njev = info.get('njev', None)
  70. else:
  71. # If there are no analytic jacobian evaluations we need to set `njev=None`.
  72. njev = None
  73. status = FROM_MINPACK_TO_COMMON[status]
  74. active_mask = np.zeros_like(x0, dtype=int)
  75. return OptimizeResult(
  76. x=x, cost=cost, fun=f, jac=J, grad=g, optimality=g_norm,
  77. active_mask=active_mask, nfev=nfev, njev=njev, status=status)
  78. def prepare_bounds(bounds, n):
  79. lb, ub = (np.asarray(b, dtype=float) for b in bounds)
  80. if lb.ndim == 0:
  81. lb = np.resize(lb, n)
  82. if ub.ndim == 0:
  83. ub = np.resize(ub, n)
  84. return lb, ub
  85. def check_tolerance(ftol, xtol, gtol, method):
  86. def check(tol, name):
  87. if tol is None:
  88. tol = 0
  89. elif tol < EPS:
  90. warn(f"Setting `{name}` below the machine epsilon ({EPS:.2e}) effectively "
  91. f"disables the corresponding termination condition.",
  92. stacklevel=3)
  93. return tol
  94. ftol = check(ftol, "ftol")
  95. xtol = check(xtol, "xtol")
  96. gtol = check(gtol, "gtol")
  97. if method == "lm" and (ftol < EPS or xtol < EPS or gtol < EPS):
  98. raise ValueError("All tolerances must be higher than machine epsilon "
  99. f"({EPS:.2e}) for method 'lm'.")
  100. elif ftol < EPS and xtol < EPS and gtol < EPS:
  101. raise ValueError("At least one of the tolerances must be higher than "
  102. f"machine epsilon ({EPS:.2e}).")
  103. return ftol, xtol, gtol
  104. def check_x_scale(x_scale, x0, method):
  105. # normalise the default scaling
  106. if x_scale is None:
  107. if method == 'lm':
  108. return 'jac'
  109. else: # dogbox, trf
  110. x_scale = 1.0
  111. if isinstance(x_scale, str) and x_scale == 'jac':
  112. return x_scale
  113. try:
  114. x_scale = np.asarray(x_scale, dtype=float)
  115. valid = np.all(np.isfinite(x_scale)) and np.all(x_scale > 0)
  116. except (ValueError, TypeError):
  117. valid = False
  118. if not valid:
  119. raise ValueError("`x_scale` must be 'jac' or array_like with "
  120. "positive numbers.")
  121. if x_scale.ndim == 0:
  122. x_scale = np.resize(x_scale, x0.shape)
  123. if x_scale.shape != x0.shape:
  124. raise ValueError("Inconsistent shapes between `x_scale` and `x0`.")
  125. return x_scale
  126. def check_jac_sparsity(jac_sparsity, m, n):
  127. if jac_sparsity is None:
  128. return None
  129. if not issparse(jac_sparsity):
  130. jac_sparsity = np.atleast_2d(jac_sparsity)
  131. if jac_sparsity.shape != (m, n):
  132. raise ValueError("`jac_sparsity` has wrong shape.")
  133. return jac_sparsity, group_columns(jac_sparsity)
  134. # Loss functions.
  135. def huber(z, rho, cost_only):
  136. mask = z <= 1
  137. rho[0, mask] = z[mask]
  138. rho[0, ~mask] = 2 * z[~mask]**0.5 - 1
  139. if cost_only:
  140. return
  141. rho[1, mask] = 1
  142. rho[1, ~mask] = z[~mask]**-0.5
  143. rho[2, mask] = 0
  144. rho[2, ~mask] = -0.5 * z[~mask]**-1.5
  145. def soft_l1(z, rho, cost_only):
  146. t = 1 + z
  147. rho[0] = 2 * (t**0.5 - 1)
  148. if cost_only:
  149. return
  150. rho[1] = t**-0.5
  151. rho[2] = -0.5 * t**-1.5
  152. def cauchy(z, rho, cost_only):
  153. rho[0] = np.log1p(z)
  154. if cost_only:
  155. return
  156. t = 1 + z
  157. rho[1] = 1 / t
  158. rho[2] = -1 / t**2
  159. def arctan(z, rho, cost_only):
  160. rho[0] = np.arctan(z)
  161. if cost_only:
  162. return
  163. t = 1 + z**2
  164. rho[1] = 1 / t
  165. rho[2] = -2 * z / t**2
  166. IMPLEMENTED_LOSSES = dict(linear=None, huber=huber, soft_l1=soft_l1,
  167. cauchy=cauchy, arctan=arctan)
  168. def construct_loss_function(m, loss, f_scale):
  169. if loss == 'linear':
  170. return None
  171. if not callable(loss):
  172. loss = IMPLEMENTED_LOSSES[loss]
  173. rho = np.empty((3, m))
  174. def loss_function(f, cost_only=False):
  175. z = (f / f_scale) ** 2
  176. loss(z, rho, cost_only=cost_only)
  177. if cost_only:
  178. return 0.5 * f_scale ** 2 * np.sum(rho[0])
  179. rho[0] *= f_scale ** 2
  180. rho[2] /= f_scale ** 2
  181. return rho
  182. else:
  183. def loss_function(f, cost_only=False):
  184. z = (f / f_scale) ** 2
  185. rho = loss(z)
  186. if cost_only:
  187. return 0.5 * f_scale ** 2 * np.sum(rho[0])
  188. rho[0] *= f_scale ** 2
  189. rho[2] /= f_scale ** 2
  190. return rho
  191. return loss_function
  192. class _WrapArgsKwargs:
  193. # Supplies a user function with args and kwargs.
  194. def __init__(self, f, args=(), kwargs=None):
  195. self.f = f
  196. self.args = args
  197. self.kwargs = kwargs or {}
  198. def __call__(self, x):
  199. return self.f(x, *self.args, **self.kwargs)
  200. @_workers_wrapper
  201. def least_squares(
  202. fun, x0, jac='2-point', bounds=(-np.inf, np.inf), method='trf',
  203. ftol=1e-8, xtol=1e-8, gtol=1e-8, x_scale=None, loss='linear',
  204. f_scale=1.0, diff_step=None, tr_solver=None, tr_options=None,
  205. jac_sparsity=None, max_nfev=None, verbose=0, args=(), kwargs=None,
  206. callback=None, workers=None
  207. ):
  208. """Solve a nonlinear least-squares problem with bounds on the variables.
  209. Given the residuals f(x) (an m-D real function of n real
  210. variables) and the loss function rho(s) (a scalar function), `least_squares`
  211. finds a local minimum of the cost function F(x)::
  212. minimize F(x) = 0.5 * sum(rho(f_i(x)**2), i = 0, ..., m - 1)
  213. subject to lb <= x <= ub
  214. The purpose of the loss function rho(s) is to reduce the influence of
  215. outliers on the solution.
  216. Parameters
  217. ----------
  218. fun : callable
  219. Function which computes the vector of residuals, with the signature
  220. ``fun(x, *args, **kwargs)``, i.e., the minimization proceeds with
  221. respect to its first argument. The argument ``x`` passed to this
  222. function is an ndarray of shape (n,) (never a scalar, even for n=1).
  223. It must allocate and return a 1-D array_like of shape (m,) or a scalar.
  224. If the argument ``x`` is complex or the function ``fun`` returns
  225. complex residuals, it must be wrapped in a real function of real
  226. arguments, as shown at the end of the Examples section.
  227. x0 : array_like with shape (n,) or float
  228. Initial guess on independent variables. If float, it will be treated
  229. as a 1-D array with one element. When `method` is 'trf', the initial
  230. guess might be slightly adjusted to lie sufficiently within the given
  231. `bounds`.
  232. jac : {'2-point', '3-point', 'cs', callable}, optional
  233. Method of computing the Jacobian matrix (an m-by-n matrix, where
  234. element (i, j) is the partial derivative of f[i] with respect to
  235. x[j]). The keywords select a finite difference scheme for numerical
  236. estimation. The scheme '3-point' is more accurate, but requires
  237. twice as many operations as '2-point' (default). The scheme 'cs'
  238. uses complex steps, and while potentially the most accurate, it is
  239. applicable only when `fun` correctly handles complex inputs and
  240. can be analytically continued to the complex plane. If callable, it is used as
  241. ``jac(x, *args, **kwargs)`` and should return a good approximation
  242. (or the exact value) for the Jacobian as an array_like (np.atleast_2d
  243. is applied), a sparse array (csr_array preferred for performance) or
  244. a `scipy.sparse.linalg.LinearOperator`.
  245. .. versionchanged:: 1.16.0
  246. An ability to use the '3-point', 'cs' keywords with the 'lm' method.
  247. Previously 'lm' was limited to '2-point' and callable.
  248. bounds : 2-tuple of array_like or `Bounds`, optional
  249. There are two ways to specify bounds:
  250. 1. Instance of `Bounds` class
  251. 2. Lower and upper bounds on independent variables. Defaults to no
  252. bounds. Each array must match the size of `x0` or be a scalar,
  253. in the latter case a bound will be the same for all variables.
  254. Use ``np.inf`` with an appropriate sign to disable bounds on all
  255. or some variables.
  256. method : {'trf', 'dogbox', 'lm'}, optional
  257. Algorithm to perform minimization.
  258. * 'trf' : Trust Region Reflective algorithm, particularly suitable
  259. for large sparse problems with bounds. Generally robust method.
  260. * 'dogbox' : dogleg algorithm with rectangular trust regions,
  261. typical use case is small problems with bounds. Not recommended
  262. for problems with rank-deficient Jacobian.
  263. * 'lm' : Levenberg-Marquardt algorithm as implemented in MINPACK.
  264. Doesn't handle bounds and sparse Jacobians. Usually the most
  265. efficient method for small unconstrained problems.
  266. Default is 'trf'. See Notes for more information.
  267. ftol : float or None, optional
  268. Tolerance for termination by the change of the cost function. Default
  269. is 1e-8. The optimization process is stopped when ``dF < ftol * F``,
  270. and there was an adequate agreement between a local quadratic model and
  271. the true model in the last step.
  272. If None and 'method' is not 'lm', the termination by this condition is
  273. disabled. If 'method' is 'lm', this tolerance must be higher than
  274. machine epsilon.
  275. xtol : float or None, optional
  276. Tolerance for termination by the change of the independent variables.
  277. Default is 1e-8. The exact condition depends on the `method` used:
  278. * For 'trf' and 'dogbox' : ``norm(dx) < xtol * (xtol + norm(x))``.
  279. * For 'lm' : ``Delta < xtol * norm(xs)``, where ``Delta`` is
  280. a trust-region radius and ``xs`` is the value of ``x``
  281. scaled according to `x_scale` parameter (see below).
  282. If None and 'method' is not 'lm', the termination by this condition is
  283. disabled. If 'method' is 'lm', this tolerance must be higher than
  284. machine epsilon.
  285. gtol : float or None, optional
  286. Tolerance for termination by the norm of the gradient. Default is 1e-8.
  287. The exact condition depends on a `method` used:
  288. * For 'trf' : ``norm(g_scaled, ord=np.inf) < gtol``, where
  289. ``g_scaled`` is the value of the gradient scaled to account for
  290. the presence of the bounds [STIR]_.
  291. * For 'dogbox' : ``norm(g_free, ord=np.inf) < gtol``, where
  292. ``g_free`` is the gradient with respect to the variables which
  293. are not in the optimal state on the boundary.
  294. * For 'lm' : the maximum absolute value of the cosine of angles
  295. between columns of the Jacobian and the residual vector is less
  296. than `gtol`, or the residual vector is zero.
  297. If None and 'method' is not 'lm', the termination by this condition is
  298. disabled. If 'method' is 'lm', this tolerance must be higher than
  299. machine epsilon.
  300. x_scale : {None, array_like, 'jac'}, optional
  301. Characteristic scale of each variable. Setting `x_scale` is equivalent
  302. to reformulating the problem in scaled variables ``xs = x / x_scale``.
  303. An alternative view is that the size of a trust region along jth
  304. dimension is proportional to ``x_scale[j]``. Improved convergence may
  305. be achieved by setting `x_scale` such that a step of a given size
  306. along any of the scaled variables has a similar effect on the cost
  307. function. If set to 'jac', the scale is iteratively updated using the
  308. inverse norms of the columns of the Jacobian matrix (as described in
  309. [JJMore]_). The default scaling for each method (i.e.
  310. if ``x_scale is None``) is as follows:
  311. * For 'trf' : ``x_scale == 1``
  312. * For 'dogbox' : ``x_scale == 1``
  313. * For 'lm' : ``x_scale == 'jac'``
  314. .. versionchanged:: 1.16.0
  315. The default keyword value is changed from 1 to None to indicate that
  316. a default approach to scaling is used.
  317. For the 'lm' method the default scaling is changed from 1 to 'jac'.
  318. This has been found to give better performance, and is the same
  319. scaling as performed by ``leastsq``.
  320. loss : str or callable, optional
  321. Determines the loss function. The following keyword values are allowed:
  322. * 'linear' (default) : ``rho(z) = z``. Gives a standard
  323. least-squares problem.
  324. * 'soft_l1' : ``rho(z) = 2 * ((1 + z)**0.5 - 1)``. The smooth
  325. approximation of l1 (absolute value) loss. Usually a good
  326. choice for robust least squares.
  327. * 'huber' : ``rho(z) = z if z <= 1 else 2*z**0.5 - 1``. Works
  328. similarly to 'soft_l1'.
  329. * 'cauchy' : ``rho(z) = ln(1 + z)``. Severely weakens outliers
  330. influence, but may cause difficulties in optimization process.
  331. * 'arctan' : ``rho(z) = arctan(z)``. Limits a maximum loss on
  332. a single residual, has properties similar to 'cauchy'.
  333. If callable, it must take a 1-D ndarray ``z=f**2`` and return an
  334. array_like with shape (3, m) where row 0 contains function values,
  335. row 1 contains first derivatives and row 2 contains second
  336. derivatives. Method 'lm' supports only 'linear' loss.
  337. f_scale : float, optional
  338. Value of soft margin between inlier and outlier residuals, default
  339. is 1.0. The loss function is evaluated as follows
  340. ``rho_(f**2) = C**2 * rho(f**2 / C**2)``, where ``C`` is `f_scale`,
  341. and ``rho`` is determined by `loss` parameter. This parameter has
  342. no effect with ``loss='linear'``, but for other `loss` values it is
  343. of crucial importance.
  344. max_nfev : None or int, optional
  345. For all methods this parameter controls the maximum number of function
  346. evaluations used by each method, separate to those used in numerical
  347. approximation of the jacobian.
  348. If None (default), the value is chosen automatically as 100 * n.
  349. .. versionchanged:: 1.16.0
  350. The default for the 'lm' method is changed to 100 * n, for both a callable
  351. and a numerically estimated jacobian. Previously the default when using an
  352. estimated jacobian was 100 * n * (n + 1), because the method included
  353. evaluations used in the estimation.
  354. diff_step : None or array_like, optional
  355. Determines the relative step size for the finite difference
  356. approximation of the Jacobian. The actual step is computed as
  357. ``x * diff_step``. If None (default), then `diff_step` is taken to be
  358. a conventional "optimal" power of machine epsilon for the finite
  359. difference scheme used [NR]_.
  360. tr_solver : {None, 'exact', 'lsmr'}, optional
  361. Method for solving trust-region subproblems, relevant only for 'trf'
  362. and 'dogbox' methods.
  363. * 'exact' is suitable for not very large problems with dense
  364. Jacobian matrices. The computational complexity per iteration is
  365. comparable to a singular value decomposition of the Jacobian
  366. matrix.
  367. * 'lsmr' is suitable for problems with sparse and large Jacobian
  368. matrices. It uses the iterative procedure
  369. `scipy.sparse.linalg.lsmr` for finding a solution of a linear
  370. least-squares problem and only requires matrix-vector product
  371. evaluations.
  372. If None (default), the solver is chosen based on the type of Jacobian
  373. returned on the first iteration.
  374. tr_options : dict, optional
  375. Keyword options passed to trust-region solver.
  376. * ``tr_solver='exact'``: `tr_options` are ignored.
  377. * ``tr_solver='lsmr'``: options for `scipy.sparse.linalg.lsmr`.
  378. Additionally, ``method='trf'`` supports 'regularize' option
  379. (bool, default is True), which adds a regularization term to the
  380. normal equation, which improves convergence if the Jacobian is
  381. rank-deficient [Byrd]_ (eq. 3.4).
  382. jac_sparsity : {None, array_like, sparse array}, optional
  383. Defines the sparsity structure of the Jacobian matrix for finite
  384. difference estimation, its shape must be (m, n). If the Jacobian has
  385. only few non-zero elements in *each* row, providing the sparsity
  386. structure will greatly speed up the computations [Curtis]_. A zero
  387. entry means that a corresponding element in the Jacobian is identically
  388. zero. If provided, forces the use of 'lsmr' trust-region solver.
  389. If None (default), then dense differencing will be used. Has no effect
  390. for 'lm' method.
  391. verbose : {0, 1, 2}, optional
  392. Level of algorithm's verbosity:
  393. * 0 (default) : work silently.
  394. * 1 : display a termination report.
  395. * 2 : display progress during iterations (not supported by 'lm'
  396. method).
  397. args, kwargs : tuple and dict, optional
  398. Additional arguments passed to `fun` and `jac`. Both empty by default.
  399. The calling signature is ``fun(x, *args, **kwargs)`` and the same for
  400. `jac`.
  401. callback : None or callable, optional
  402. Callback function that is called by the algorithm on each iteration.
  403. This can be used to print or plot the optimization results at each
  404. step, and to stop the optimization algorithm based on some user-defined
  405. condition. Only implemented for the `trf` and `dogbox` methods.
  406. The signature is ``callback(intermediate_result: OptimizeResult)``
  407. `intermediate_result is a `scipy.optimize.OptimizeResult`
  408. which contains the intermediate results of the optimization at the
  409. current iteration.
  410. The callback also supports a signature like: ``callback(x)``
  411. Introspection is used to determine which of the signatures is invoked.
  412. If the `callback` function raises `StopIteration` the optimization algorithm
  413. will stop and return with status code -2.
  414. .. versionadded:: 1.16.0
  415. workers : map-like callable, optional
  416. A map-like callable, such as `multiprocessing.Pool.map` for evaluating
  417. any numerical differentiation in parallel.
  418. This evaluation is carried out as ``workers(fun, iterable)``.
  419. .. versionadded:: 1.16.0
  420. Returns
  421. -------
  422. result : OptimizeResult
  423. `OptimizeResult` with the following fields defined:
  424. x : ndarray, shape (n,)
  425. Solution found.
  426. cost : float
  427. Value of the cost function at the solution.
  428. fun : ndarray, shape (m,)
  429. Vector of residuals at the solution.
  430. jac : ndarray, sparse array or LinearOperator, shape (m, n)
  431. Modified Jacobian matrix at the solution, in the sense that J^T J
  432. is a Gauss-Newton approximation of the Hessian of the cost function.
  433. The type is the same as the one used by the algorithm.
  434. grad : ndarray, shape (m,)
  435. Gradient of the cost function at the solution.
  436. optimality : float
  437. First-order optimality measure. In unconstrained problems, it is
  438. always the uniform norm of the gradient. In constrained problems,
  439. it is the quantity which was compared with `gtol` during iterations.
  440. active_mask : ndarray of int, shape (n,)
  441. Each component shows whether a corresponding constraint is active
  442. (that is, whether a variable is at the bound):
  443. * 0 : a constraint is not active.
  444. * -1 : a lower bound is active.
  445. * 1 : an upper bound is active.
  446. Might be somewhat arbitrary for 'trf' method as it generates a
  447. sequence of strictly feasible iterates and `active_mask` is
  448. determined within a tolerance threshold.
  449. nfev : int
  450. Number of function evaluations done. This number does not include
  451. the function calls used for numerical Jacobian approximation.
  452. .. versionchanged:: 1.16.0
  453. For the 'lm' method the number of function calls used in numerical
  454. Jacobian approximation is no longer included. This is to bring all
  455. methods into line.
  456. njev : int or None
  457. Number of Jacobian evaluations done. If numerical Jacobian
  458. approximation is used in 'lm' method, it is set to None.
  459. status : int
  460. The reason for algorithm termination:
  461. * -2 : terminated because callback raised StopIteration.
  462. * -1 : improper input parameters status returned from MINPACK.
  463. * 0 : the maximum number of function evaluations is exceeded.
  464. * 1 : `gtol` termination condition is satisfied.
  465. * 2 : `ftol` termination condition is satisfied.
  466. * 3 : `xtol` termination condition is satisfied.
  467. * 4 : Both `ftol` and `xtol` termination conditions are satisfied.
  468. message : str
  469. Verbal description of the termination reason.
  470. success : bool
  471. True if one of the convergence criteria is satisfied (`status` > 0).
  472. See Also
  473. --------
  474. leastsq : A legacy wrapper for the MINPACK implementation of the
  475. Levenberg-Marquadt algorithm.
  476. curve_fit : Least-squares minimization applied to a curve-fitting problem.
  477. Notes
  478. -----
  479. Method 'lm' (Levenberg-Marquardt) calls a wrapper over a least-squares
  480. algorithm implemented in MINPACK (lmder). It runs the
  481. Levenberg-Marquardt algorithm formulated as a trust-region type algorithm.
  482. The implementation is based on paper [JJMore]_, it is very robust and
  483. efficient with a lot of smart tricks. It should be your first choice
  484. for unconstrained problems. Note that it doesn't support bounds. Also,
  485. it doesn't work when m < n.
  486. Method 'trf' (Trust Region Reflective) is motivated by the process of
  487. solving a system of equations, which constitute the first-order optimality
  488. condition for a bound-constrained minimization problem as formulated in
  489. [STIR]_. The algorithm iteratively solves trust-region subproblems
  490. augmented by a special diagonal quadratic term and with trust-region shape
  491. determined by the distance from the bounds and the direction of the
  492. gradient. This enhancements help to avoid making steps directly into bounds
  493. and efficiently explore the whole space of variables. To further improve
  494. convergence, the algorithm considers search directions reflected from the
  495. bounds. To obey theoretical requirements, the algorithm keeps iterates
  496. strictly feasible. With dense Jacobians trust-region subproblems are
  497. solved by an exact method very similar to the one described in [JJMore]_
  498. (and implemented in MINPACK). The difference from the MINPACK
  499. implementation is that a singular value decomposition of a Jacobian
  500. matrix is done once per iteration, instead of a QR decomposition and series
  501. of Givens rotation eliminations. For large sparse Jacobians a 2-D subspace
  502. approach of solving trust-region subproblems is used [STIR]_, [Byrd]_.
  503. The subspace is spanned by a scaled gradient and an approximate
  504. Gauss-Newton solution delivered by `scipy.sparse.linalg.lsmr`. When no
  505. constraints are imposed the algorithm is very similar to MINPACK and has
  506. generally comparable performance. The algorithm works quite robust in
  507. unbounded and bounded problems, thus it is chosen as a default algorithm.
  508. Method 'dogbox' operates in a trust-region framework, but considers
  509. rectangular trust regions as opposed to conventional ellipsoids [Voglis]_.
  510. The intersection of a current trust region and initial bounds is again
  511. rectangular, so on each iteration a quadratic minimization problem subject
  512. to bound constraints is solved approximately by Powell's dogleg method
  513. [NumOpt]_. The required Gauss-Newton step can be computed exactly for
  514. dense Jacobians or approximately by `scipy.sparse.linalg.lsmr` for large
  515. sparse Jacobians. The algorithm is likely to exhibit slow convergence when
  516. the rank of Jacobian is less than the number of variables. The algorithm
  517. often outperforms 'trf' in bounded problems with a small number of
  518. variables.
  519. Robust loss functions are implemented as described in [BA]_. The idea
  520. is to modify a residual vector and a Jacobian matrix on each iteration
  521. such that computed gradient and Gauss-Newton Hessian approximation match
  522. the true gradient and Hessian approximation of the cost function. Then
  523. the algorithm proceeds in a normal way, i.e., robust loss functions are
  524. implemented as a simple wrapper over standard least-squares algorithms.
  525. .. versionadded:: 0.17.0
  526. References
  527. ----------
  528. .. [STIR] M. A. Branch, T. F. Coleman, and Y. Li, "A Subspace, Interior,
  529. and Conjugate Gradient Method for Large-Scale Bound-Constrained
  530. Minimization Problems," SIAM Journal on Scientific Computing,
  531. Vol. 21, Number 1, pp 1-23, 1999.
  532. .. [NR] William H. Press et. al., "Numerical Recipes. The Art of Scientific
  533. Computing. 3rd edition", Sec. 5.7.
  534. .. [Byrd] R. H. Byrd, R. B. Schnabel and G. A. Shultz, "Approximate
  535. solution of the trust region problem by minimization over
  536. two-dimensional subspaces", Math. Programming, 40, pp. 247-263,
  537. 1988.
  538. .. [Curtis] A. Curtis, M. J. D. Powell, and J. Reid, "On the estimation of
  539. sparse Jacobian matrices", Journal of the Institute of
  540. Mathematics and its Applications, 13, pp. 117-120, 1974.
  541. .. [JJMore] J. J. More, "The Levenberg-Marquardt Algorithm: Implementation
  542. and Theory," Numerical Analysis, ed. G. A. Watson, Lecture
  543. Notes in Mathematics 630, Springer Verlag, pp. 105-116, 1977.
  544. .. [Voglis] C. Voglis and I. E. Lagaris, "A Rectangular Trust Region
  545. Dogleg Approach for Unconstrained and Bound Constrained
  546. Nonlinear Optimization", WSEAS International Conference on
  547. Applied Mathematics, Corfu, Greece, 2004.
  548. .. [NumOpt] J. Nocedal and S. J. Wright, "Numerical optimization,
  549. 2nd edition", Chapter 4.
  550. .. [BA] B. Triggs et. al., "Bundle Adjustment - A Modern Synthesis",
  551. Proceedings of the International Workshop on Vision Algorithms:
  552. Theory and Practice, pp. 298-372, 1999.
  553. Examples
  554. --------
  555. In this example we find a minimum of the Rosenbrock function without bounds
  556. on independent variables.
  557. >>> import numpy as np
  558. >>> def fun_rosenbrock(x):
  559. ... return np.array([10 * (x[1] - x[0]**2), (1 - x[0])])
  560. Notice that we only provide the vector of the residuals. The algorithm
  561. constructs the cost function as a sum of squares of the residuals, which
  562. gives the Rosenbrock function. The exact minimum is at ``x = [1.0, 1.0]``.
  563. >>> from scipy.optimize import least_squares
  564. >>> x0_rosenbrock = np.array([2, 2])
  565. >>> res_1 = least_squares(fun_rosenbrock, x0_rosenbrock)
  566. >>> res_1.x
  567. array([ 1., 1.])
  568. >>> res_1.cost
  569. 9.8669242910846867e-30
  570. >>> res_1.optimality
  571. 8.8928864934219529e-14
  572. We now constrain the variables, in such a way that the previous solution
  573. becomes infeasible. Specifically, we require that ``x[1] >= 1.5``, and
  574. ``x[0]`` left unconstrained. To this end, we specify the `bounds` parameter
  575. to `least_squares` in the form ``bounds=([-np.inf, 1.5], np.inf)``.
  576. We also provide the analytic Jacobian:
  577. >>> def jac_rosenbrock(x):
  578. ... return np.array([
  579. ... [-20 * x[0], 10],
  580. ... [-1, 0]])
  581. Putting this all together, we see that the new solution lies on the bound:
  582. >>> res_2 = least_squares(fun_rosenbrock, x0_rosenbrock, jac_rosenbrock,
  583. ... bounds=([-np.inf, 1.5], np.inf))
  584. >>> res_2.x
  585. array([ 1.22437075, 1.5 ])
  586. >>> res_2.cost
  587. 0.025213093946805685
  588. >>> res_2.optimality
  589. 1.5885401433157753e-07
  590. Now we solve a system of equations (i.e., the cost function should be zero
  591. at a minimum) for a Broyden tridiagonal vector-valued function of 100000
  592. variables:
  593. >>> def fun_broyden(x):
  594. ... f = (3 - x) * x + 1
  595. ... f[1:] -= x[:-1]
  596. ... f[:-1] -= 2 * x[1:]
  597. ... return f
  598. The corresponding Jacobian matrix is sparse. We tell the algorithm to
  599. estimate it by finite differences and provide the sparsity structure of
  600. Jacobian to significantly speed up this process.
  601. >>> from scipy.sparse import lil_array
  602. >>> def sparsity_broyden(n):
  603. ... sparsity = lil_array((n, n), dtype=int)
  604. ... i = np.arange(n)
  605. ... sparsity[i, i] = 1
  606. ... i = np.arange(1, n)
  607. ... sparsity[i, i - 1] = 1
  608. ... i = np.arange(n - 1)
  609. ... sparsity[i, i + 1] = 1
  610. ... return sparsity
  611. ...
  612. >>> n = 100000
  613. >>> x0_broyden = -np.ones(n)
  614. ...
  615. >>> res_3 = least_squares(fun_broyden, x0_broyden,
  616. ... jac_sparsity=sparsity_broyden(n))
  617. >>> res_3.cost
  618. 4.5687069299604613e-23
  619. >>> res_3.optimality
  620. 1.1650454296851518e-11
  621. Let's also solve a curve fitting problem using robust loss function to
  622. take care of outliers in the data. Define the model function as
  623. ``y = a + b * exp(c * t)``, where t is a predictor variable, y is an
  624. observation and a, b, c are parameters to estimate.
  625. First, define the function which generates the data with noise and
  626. outliers, define the model parameters, and generate data:
  627. >>> from numpy.random import default_rng
  628. >>> rng = default_rng()
  629. >>> def gen_data(t, a, b, c, noise=0., n_outliers=0, seed=None):
  630. ... rng = default_rng(seed)
  631. ...
  632. ... y = a + b * np.exp(t * c)
  633. ...
  634. ... error = noise * rng.standard_normal(t.size)
  635. ... outliers = rng.integers(0, t.size, n_outliers)
  636. ... error[outliers] *= 10
  637. ...
  638. ... return y + error
  639. ...
  640. >>> a = 0.5
  641. >>> b = 2.0
  642. >>> c = -1
  643. >>> t_min = 0
  644. >>> t_max = 10
  645. >>> n_points = 15
  646. ...
  647. >>> t_train = np.linspace(t_min, t_max, n_points)
  648. >>> y_train = gen_data(t_train, a, b, c, noise=0.1, n_outliers=3)
  649. Define function for computing residuals and initial estimate of
  650. parameters.
  651. >>> def fun(x, t, y):
  652. ... return x[0] + x[1] * np.exp(x[2] * t) - y
  653. ...
  654. >>> x0 = np.array([1.0, 1.0, 0.0])
  655. Compute a standard least-squares solution:
  656. >>> res_lsq = least_squares(fun, x0, args=(t_train, y_train))
  657. Now compute two solutions with two different robust loss functions. The
  658. parameter `f_scale` is set to 0.1, meaning that inlier residuals should
  659. not significantly exceed 0.1 (the noise level used).
  660. >>> res_soft_l1 = least_squares(fun, x0, loss='soft_l1', f_scale=0.1,
  661. ... args=(t_train, y_train))
  662. >>> res_log = least_squares(fun, x0, loss='cauchy', f_scale=0.1,
  663. ... args=(t_train, y_train))
  664. And, finally, plot all the curves. We see that by selecting an appropriate
  665. `loss` we can get estimates close to optimal even in the presence of
  666. strong outliers. But keep in mind that generally it is recommended to try
  667. 'soft_l1' or 'huber' losses first (if at all necessary) as the other two
  668. options may cause difficulties in optimization process.
  669. >>> t_test = np.linspace(t_min, t_max, n_points * 10)
  670. >>> y_true = gen_data(t_test, a, b, c)
  671. >>> y_lsq = gen_data(t_test, *res_lsq.x)
  672. >>> y_soft_l1 = gen_data(t_test, *res_soft_l1.x)
  673. >>> y_log = gen_data(t_test, *res_log.x)
  674. ...
  675. >>> import matplotlib.pyplot as plt
  676. >>> plt.plot(t_train, y_train, 'o')
  677. >>> plt.plot(t_test, y_true, 'k', linewidth=2, label='true')
  678. >>> plt.plot(t_test, y_lsq, label='linear loss')
  679. >>> plt.plot(t_test, y_soft_l1, label='soft_l1 loss')
  680. >>> plt.plot(t_test, y_log, label='cauchy loss')
  681. >>> plt.xlabel("t")
  682. >>> plt.ylabel("y")
  683. >>> plt.legend()
  684. >>> plt.show()
  685. In the next example, we show how complex-valued residual functions of
  686. complex variables can be optimized with ``least_squares()``. Consider the
  687. following function:
  688. >>> def f(z):
  689. ... return z - (0.5 + 0.5j)
  690. We wrap it into a function of real variables that returns real residuals
  691. by simply handling the real and imaginary parts as independent variables:
  692. >>> def f_wrap(x):
  693. ... fx = f(x[0] + 1j*x[1])
  694. ... return np.array([fx.real, fx.imag])
  695. Thus, instead of the original m-D complex function of n complex
  696. variables we optimize a 2m-D real function of 2n real variables:
  697. >>> from scipy.optimize import least_squares
  698. >>> res_wrapped = least_squares(f_wrap, (0.1, 0.1), bounds=([0, 0], [1, 1]))
  699. >>> z = res_wrapped.x[0] + res_wrapped.x[1]*1j
  700. >>> z
  701. (0.49999999999925893+0.49999999999925893j)
  702. """
  703. if method not in ['trf', 'dogbox', 'lm']:
  704. raise ValueError("`method` must be 'trf', 'dogbox' or 'lm'.")
  705. if jac not in ['2-point', '3-point', 'cs'] and not callable(jac):
  706. raise ValueError("`jac` must be '2-point', '3-point', 'cs' or "
  707. "callable.")
  708. if tr_solver not in [None, 'exact', 'lsmr']:
  709. raise ValueError("`tr_solver` must be None, 'exact' or 'lsmr'.")
  710. if loss not in IMPLEMENTED_LOSSES and not callable(loss):
  711. raise ValueError(f"`loss` must be one of {IMPLEMENTED_LOSSES.keys()}"
  712. " or a callable.")
  713. if method == 'lm' and loss != 'linear':
  714. raise ValueError("method='lm' supports only 'linear' loss function.")
  715. if verbose not in [0, 1, 2]:
  716. raise ValueError("`verbose` must be in [0, 1, 2].")
  717. if max_nfev is not None and max_nfev <= 0:
  718. raise ValueError("`max_nfev` must be None or positive integer.")
  719. if np.iscomplexobj(x0):
  720. raise ValueError("`x0` must be real.")
  721. x0 = np.atleast_1d(x0).astype(float)
  722. if x0.ndim > 1:
  723. raise ValueError("`x0` must have at most 1 dimension.")
  724. if isinstance(bounds, Bounds):
  725. lb, ub = bounds.lb, bounds.ub
  726. bounds = (lb, ub)
  727. else:
  728. if len(bounds) == 2:
  729. lb, ub = prepare_bounds(bounds, x0.shape[0])
  730. else:
  731. raise ValueError("`bounds` must contain 2 elements.")
  732. if method == 'lm' and not np.all((lb == -np.inf) & (ub == np.inf)):
  733. raise ValueError("Method 'lm' doesn't support bounds.")
  734. if lb.shape != x0.shape or ub.shape != x0.shape:
  735. raise ValueError("Inconsistent shapes between bounds and `x0`.")
  736. if np.any(lb >= ub):
  737. raise ValueError("Each lower bound must be strictly less than each "
  738. "upper bound.")
  739. if not in_bounds(x0, lb, ub):
  740. raise ValueError("Initial guess is outside of provided bounds")
  741. x_scale = check_x_scale(x_scale, x0, method)
  742. ftol, xtol, gtol = check_tolerance(ftol, xtol, gtol, method)
  743. if method == 'trf':
  744. x0 = make_strictly_feasible(x0, lb, ub)
  745. if tr_options is None:
  746. tr_options = {}
  747. ###########################################################################
  748. # assemble VectorFunction
  749. ###########################################################################
  750. # first wrap the args/kwargs
  751. fun_wrapped = _WrapArgsKwargs(fun, args=args, kwargs=kwargs)
  752. jac_wrapped = jac
  753. if callable(jac):
  754. jac_wrapped = _WrapArgsKwargs(jac, args=args, kwargs=kwargs)
  755. def _dummy_hess(x, *args):
  756. # we don't care about Hessian evaluations
  757. return x
  758. vector_fun = VectorFunction(
  759. fun_wrapped,
  760. x0,
  761. jac_wrapped,
  762. _dummy_hess,
  763. finite_diff_rel_step=diff_step,
  764. finite_diff_jac_sparsity=jac_sparsity,
  765. finite_diff_bounds=bounds,
  766. workers=workers
  767. )
  768. ###########################################################################
  769. f0 = vector_fun.fun(x0)
  770. J0 = vector_fun.jac(x0)
  771. if f0.ndim != 1:
  772. raise ValueError("`fun` must return at most 1-d array_like. "
  773. f"f0.shape: {f0.shape}")
  774. if not np.all(np.isfinite(f0)):
  775. raise ValueError("Residuals are not finite in the initial point.")
  776. n = x0.size
  777. m = f0.size
  778. if method == 'lm' and m < n:
  779. raise ValueError("Method 'lm' doesn't work when the number of "
  780. "residuals is less than the number of variables.")
  781. loss_function = construct_loss_function(m, loss, f_scale)
  782. if callable(loss):
  783. rho = loss_function(f0)
  784. if rho.shape != (3, m):
  785. raise ValueError("The return value of `loss` callable has wrong "
  786. "shape.")
  787. initial_cost = 0.5 * np.sum(rho[0])
  788. elif loss_function is not None:
  789. initial_cost = loss_function(f0, cost_only=True)
  790. else:
  791. initial_cost = 0.5 * np.dot(f0, f0)
  792. if not callable(jac):
  793. # Estimate Jacobian by finite differences.
  794. if method == 'lm':
  795. if jac_sparsity is not None:
  796. raise ValueError("method='lm' does not support "
  797. "`jac_sparsity`.")
  798. else:
  799. # this will raise a ValueError if the jac_sparsity isn't correct
  800. _ = check_jac_sparsity(jac_sparsity, m, n)
  801. if jac_sparsity is not None and tr_solver == 'exact':
  802. raise ValueError("tr_solver='exact' is incompatible "
  803. "with `jac_sparsity`.")
  804. if J0.shape != (m, n):
  805. raise ValueError(
  806. f"The return value of `jac` has wrong shape: expected {(m, n)}, "
  807. f"actual {J0.shape}."
  808. )
  809. if not isinstance(J0, np.ndarray):
  810. if method == 'lm':
  811. raise ValueError("method='lm' works only with dense "
  812. "Jacobian matrices.")
  813. if tr_solver == 'exact':
  814. raise ValueError(
  815. "tr_solver='exact' works only with dense "
  816. "Jacobian matrices.")
  817. jac_scale = isinstance(x_scale, str) and x_scale == 'jac'
  818. if isinstance(J0, LinearOperator) and jac_scale:
  819. raise ValueError("x_scale='jac' can't be used when `jac` "
  820. "returns LinearOperator.")
  821. if tr_solver is None:
  822. if isinstance(J0, np.ndarray):
  823. tr_solver = 'exact'
  824. else:
  825. tr_solver = 'lsmr'
  826. # Wrap callback function. If callback is None, callback_wrapped also is None
  827. callback_wrapped = _wrap_callback(callback)
  828. if method == 'lm':
  829. if callback is not None:
  830. warn("Callback function specified, but not supported with `lm` method.",
  831. stacklevel=2)
  832. result = call_minpack(vector_fun.fun, x0, vector_fun.jac, ftol, xtol, gtol,
  833. max_nfev, x_scale, jac_method=jac)
  834. elif method == 'trf':
  835. result = trf(vector_fun.fun, vector_fun.jac, x0, f0, J0, lb, ub, ftol, xtol,
  836. gtol, max_nfev, x_scale, loss_function, tr_solver,
  837. tr_options.copy(), verbose, callback=callback_wrapped)
  838. elif method == 'dogbox':
  839. if tr_solver == 'lsmr' and 'regularize' in tr_options:
  840. warn("The keyword 'regularize' in `tr_options` is not relevant "
  841. "for 'dogbox' method.",
  842. stacklevel=2)
  843. tr_options = tr_options.copy()
  844. del tr_options['regularize']
  845. result = dogbox(vector_fun.fun, vector_fun.jac, x0, f0, J0, lb, ub, ftol,
  846. xtol, gtol, max_nfev, x_scale, loss_function,
  847. tr_solver, tr_options, verbose, callback=callback_wrapped)
  848. result.message = TERMINATION_MESSAGES[result.status]
  849. result.success = result.status > 0
  850. if verbose >= 1:
  851. print(result.message)
  852. print(f"Function evaluations {result.nfev}, initial cost {initial_cost:.4e}, "
  853. f"final cost {result.cost:.4e}, "
  854. f"first-order optimality {result.optimality:.2e}.")
  855. return result