least_squares.py 39 KB

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