_minpack_py.py 44 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178
  1. import warnings
  2. from . import _minpack
  3. import numpy as np
  4. from numpy import (atleast_1d, triu, shape, transpose, zeros, prod, greater,
  5. asarray, inf,
  6. finfo, inexact, issubdtype, dtype)
  7. from scipy import linalg
  8. from scipy.linalg import svd, cholesky, solve_triangular, LinAlgError
  9. from scipy._lib._util import _asarray_validated, _contains_nan
  10. from scipy._lib._util import getfullargspec_no_self as _getfullargspec
  11. import scipy._lib.array_api_extra as xpx
  12. from ._optimize import OptimizeResult, _check_unknown_options, OptimizeWarning
  13. from ._lsq import least_squares
  14. # from ._lsq.common import make_strictly_feasible
  15. from ._lsq.least_squares import prepare_bounds
  16. from scipy.optimize._minimize import Bounds
  17. __all__ = ['fsolve', 'leastsq', 'fixed_point', 'curve_fit']
  18. def _check_func(checker, argname, thefunc, x0, args, numinputs,
  19. output_shape=None):
  20. res = atleast_1d(thefunc(*((x0[:numinputs],) + args)))
  21. if (output_shape is not None) and (shape(res) != output_shape):
  22. if (output_shape[0] != 1):
  23. if len(output_shape) > 1:
  24. if output_shape[1] == 1:
  25. return shape(res)
  26. msg = f"{checker}: there is a mismatch between the input and output " \
  27. f"shape of the '{argname}' argument"
  28. func_name = getattr(thefunc, '__name__', None)
  29. if func_name:
  30. msg += f" '{func_name}'."
  31. else:
  32. msg += "."
  33. msg += f'Shape should be {output_shape} but it is {shape(res)}.'
  34. raise TypeError(msg)
  35. if issubdtype(res.dtype, inexact):
  36. dt = res.dtype
  37. else:
  38. dt = dtype(float)
  39. return shape(res), dt
  40. def fsolve(func, x0, args=(), fprime=None, full_output=0,
  41. col_deriv=0, xtol=1.49012e-8, maxfev=0, band=None,
  42. epsfcn=None, factor=100, diag=None):
  43. """
  44. Find the roots of a function.
  45. Return the roots of the (non-linear) equations defined by
  46. ``func(x) = 0`` given a starting estimate.
  47. Parameters
  48. ----------
  49. func : callable ``f(x, *args)``
  50. A function that takes at least one (possibly vector) argument,
  51. and returns a value of the same length.
  52. x0 : ndarray
  53. The starting estimate for the roots of ``func(x) = 0``.
  54. args : tuple, optional
  55. Any extra arguments to `func`.
  56. fprime : callable ``f(x, *args)``, optional
  57. A function to compute the Jacobian of `func` with derivatives
  58. across the rows. By default, the Jacobian will be estimated.
  59. full_output : bool, optional
  60. If True, return optional outputs.
  61. col_deriv : bool, optional
  62. Specify whether the Jacobian function computes derivatives down
  63. the columns (faster, because there is no transpose operation).
  64. xtol : float, optional
  65. The calculation will terminate if the relative error between two
  66. consecutive iterates is at most `xtol`.
  67. maxfev : int, optional
  68. The maximum number of calls to the function. If zero, then
  69. ``100*(N+1)`` is the maximum where N is the number of elements
  70. in `x0`.
  71. band : tuple, optional
  72. If set to a two-sequence containing the number of sub- and
  73. super-diagonals within the band of the Jacobi matrix, the
  74. Jacobi matrix is considered banded (only for ``fprime=None``).
  75. epsfcn : float, optional
  76. A suitable step length for the forward-difference
  77. approximation of the Jacobian (for ``fprime=None``). If
  78. `epsfcn` is less than the machine precision, it is assumed
  79. that the relative errors in the functions are of the order of
  80. the machine precision.
  81. factor : float, optional
  82. A parameter determining the initial step bound
  83. (``factor * || diag * x||``). Should be in the interval
  84. ``(0.1, 100)``.
  85. diag : sequence, optional
  86. N positive entries that serve as a scale factors for the
  87. variables.
  88. Returns
  89. -------
  90. x : ndarray
  91. The solution (or the result of the last iteration for
  92. an unsuccessful call).
  93. infodict : dict
  94. A dictionary of optional outputs with the keys:
  95. ``nfev``
  96. number of function calls
  97. ``njev``
  98. number of Jacobian calls
  99. ``fvec``
  100. function evaluated at the output
  101. ``fjac``
  102. the orthogonal matrix, q, produced by the QR
  103. factorization of the final approximate Jacobian
  104. matrix, stored column wise
  105. ``r``
  106. upper triangular matrix produced by QR factorization
  107. of the same matrix
  108. ``qtf``
  109. the vector ``(transpose(q) * fvec)``
  110. ier : int
  111. An integer flag. Set to 1 if a solution was found, otherwise refer
  112. to `mesg` for more information.
  113. mesg : str
  114. If no solution is found, `mesg` details the cause of failure.
  115. See Also
  116. --------
  117. root : Interface to root finding algorithms for multivariate
  118. functions. See the ``method='hybr'`` in particular.
  119. Notes
  120. -----
  121. ``fsolve`` is a wrapper around MINPACK's hybrd and hybrj algorithms.
  122. Examples
  123. --------
  124. Find a solution to the system of equations:
  125. ``x0*cos(x1) = 4, x1*x0 - x1 = 5``.
  126. >>> import numpy as np
  127. >>> from scipy.optimize import fsolve
  128. >>> def func(x):
  129. ... return [x[0] * np.cos(x[1]) - 4,
  130. ... x[1] * x[0] - x[1] - 5]
  131. >>> root = fsolve(func, [1, 1])
  132. >>> root
  133. array([6.50409711, 0.90841421])
  134. >>> np.isclose(func(root), [0.0, 0.0]) # func(root) should be almost 0.0.
  135. array([ True, True])
  136. """
  137. def _wrapped_func(*fargs):
  138. """
  139. Wrapped `func` to track the number of times
  140. the function has been called.
  141. """
  142. _wrapped_func.nfev += 1
  143. return func(*fargs)
  144. _wrapped_func.nfev = 0
  145. options = {'col_deriv': col_deriv,
  146. 'xtol': xtol,
  147. 'maxfev': maxfev,
  148. 'band': band,
  149. 'eps': epsfcn,
  150. 'factor': factor,
  151. 'diag': diag}
  152. res = _root_hybr(_wrapped_func, x0, args, jac=fprime, **options)
  153. res.nfev = _wrapped_func.nfev
  154. if full_output:
  155. x = res['x']
  156. info = {k: res.get(k)
  157. for k in ('nfev', 'njev', 'fjac', 'r', 'qtf') if k in res}
  158. info['fvec'] = res['fun']
  159. return x, info, res['status'], res['message']
  160. else:
  161. status = res['status']
  162. msg = res['message']
  163. if status == 0:
  164. raise TypeError(msg)
  165. elif status == 1:
  166. pass
  167. elif status in [2, 3, 4, 5]:
  168. warnings.warn(msg, RuntimeWarning, stacklevel=2)
  169. else:
  170. raise TypeError(msg)
  171. return res['x']
  172. def _root_hybr(func, x0, args=(), jac=None,
  173. col_deriv=0, xtol=1.49012e-08, maxfev=0, band=None, eps=None,
  174. factor=100, diag=None, **unknown_options):
  175. """
  176. Find the roots of a multivariate function using MINPACK's hybrd and
  177. hybrj routines (modified Powell method).
  178. Options
  179. -------
  180. col_deriv : bool
  181. Specify whether the Jacobian function computes derivatives down
  182. the columns (faster, because there is no transpose operation).
  183. xtol : float
  184. The calculation will terminate if the relative error between two
  185. consecutive iterates is at most `xtol`.
  186. maxfev : int
  187. The maximum number of calls to the function. If zero, then
  188. ``100*(N+1)`` is the maximum where N is the number of elements
  189. in `x0`.
  190. band : tuple
  191. If set to a two-sequence containing the number of sub- and
  192. super-diagonals within the band of the Jacobi matrix, the
  193. Jacobi matrix is considered banded (only for ``jac=None``).
  194. eps : float
  195. A suitable step length for the forward-difference
  196. approximation of the Jacobian (for ``jac=None``). If
  197. `eps` is less than the machine precision, it is assumed
  198. that the relative errors in the functions are of the order of
  199. the machine precision.
  200. factor : float
  201. A parameter determining the initial step bound
  202. (``factor * || diag * x||``). Should be in the interval
  203. ``(0.1, 100)``.
  204. diag : sequence
  205. N positive entries that serve as a scale factors for the
  206. variables.
  207. """
  208. _check_unknown_options(unknown_options)
  209. epsfcn = eps
  210. x0 = asarray(x0).flatten()
  211. n = len(x0)
  212. if not isinstance(args, tuple):
  213. args = (args,)
  214. shape, dtype = _check_func('fsolve', 'func', func, x0, args, n, (n,))
  215. if epsfcn is None:
  216. epsfcn = finfo(dtype).eps
  217. Dfun = jac
  218. if Dfun is None:
  219. if band is None:
  220. ml, mu = -10, -10
  221. else:
  222. ml, mu = band[:2]
  223. if maxfev == 0:
  224. maxfev = 200 * (n + 1)
  225. retval = _minpack._hybrd(func, x0, args, 1, xtol, maxfev,
  226. ml, mu, epsfcn, factor, diag)
  227. else:
  228. _check_func('fsolve', 'fprime', Dfun, x0, args, n, (n, n))
  229. if (maxfev == 0):
  230. maxfev = 100 * (n + 1)
  231. retval = _minpack._hybrj(func, Dfun, x0, args, 1,
  232. col_deriv, xtol, maxfev, factor, diag)
  233. x, status = retval[0], retval[-1]
  234. errors = {0: "Improper input parameters were entered.",
  235. 1: "The solution converged.",
  236. 2: "The number of calls to function has "
  237. f"reached maxfev = {maxfev}.",
  238. 3: f"xtol={xtol:f} is too small, no further improvement "
  239. "in the approximate\n solution is possible.",
  240. 4: "The iteration is not making good progress, as measured "
  241. "by the \n improvement from the last five "
  242. "Jacobian evaluations.",
  243. 5: "The iteration is not making good progress, "
  244. "as measured by the \n improvement from the last "
  245. "ten iterations.",
  246. 'unknown': "An error occurred."}
  247. info = retval[1]
  248. info['fun'] = info.pop('fvec')
  249. sol = OptimizeResult(x=x, success=(status == 1), status=status,
  250. method="hybr")
  251. sol.update(info)
  252. try:
  253. sol['message'] = errors[status]
  254. except KeyError:
  255. sol['message'] = errors['unknown']
  256. return sol
  257. LEASTSQ_SUCCESS = [1, 2, 3, 4]
  258. LEASTSQ_FAILURE = [5, 6, 7, 8]
  259. def leastsq(func, x0, args=(), Dfun=None, full_output=False,
  260. col_deriv=False, ftol=1.49012e-8, xtol=1.49012e-8,
  261. gtol=0.0, maxfev=0, epsfcn=None, factor=100, diag=None):
  262. """
  263. Minimize the sum of squares of a set of equations.
  264. ::
  265. x = arg min(sum(func(y)**2,axis=0))
  266. y
  267. Parameters
  268. ----------
  269. func : callable
  270. Should take at least one (possibly length ``N`` vector) argument and
  271. returns ``M`` floating point numbers. It must not return NaNs or
  272. fitting might fail. ``M`` must be greater than or equal to ``N``.
  273. x0 : ndarray
  274. The starting estimate for the minimization.
  275. args : tuple, optional
  276. Any extra arguments to func are placed in this tuple.
  277. Dfun : callable, optional
  278. A function or method to compute the Jacobian of func with derivatives
  279. across the rows. If this is None, the Jacobian will be estimated.
  280. full_output : bool, optional
  281. If ``True``, return all optional outputs (not just `x` and `ier`).
  282. col_deriv : bool, optional
  283. If ``True``, specify that the Jacobian function computes derivatives
  284. down the columns (faster, because there is no transpose operation).
  285. ftol : float, optional
  286. Relative error desired in the sum of squares.
  287. xtol : float, optional
  288. Relative error desired in the approximate solution.
  289. gtol : float, optional
  290. Orthogonality desired between the function vector and the columns of
  291. the Jacobian.
  292. maxfev : int, optional
  293. The maximum number of calls to the function. If `Dfun` is provided,
  294. then the default `maxfev` is 100*(N+1) where N is the number of elements
  295. in x0, otherwise the default `maxfev` is 200*(N+1).
  296. epsfcn : float, optional
  297. A variable used in determining a suitable step length for the forward-
  298. difference approximation of the Jacobian (for Dfun=None).
  299. Normally the actual step length will be sqrt(epsfcn)*x
  300. If epsfcn is less than the machine precision, it is assumed that the
  301. relative errors are of the order of the machine precision.
  302. factor : float, optional
  303. A parameter determining the initial step bound
  304. (``factor * || diag * x||``). Should be in interval ``(0.1, 100)``.
  305. diag : sequence, optional
  306. N positive entries that serve as a scale factors for the variables.
  307. Returns
  308. -------
  309. x : ndarray
  310. The solution (or the result of the last iteration for an unsuccessful
  311. call).
  312. cov_x : ndarray
  313. The inverse of the Hessian. `fjac` and `ipvt` are used to construct an
  314. estimate of the Hessian. A value of None indicates a singular matrix,
  315. which means the curvature in parameters `x` is numerically flat. To
  316. obtain the covariance matrix of the parameters `x`, `cov_x` must be
  317. multiplied by the variance of the residuals -- see curve_fit. Only
  318. returned if `full_output` is ``True``.
  319. infodict : dict
  320. a dictionary of optional outputs with the keys:
  321. ``nfev``
  322. The number of function calls
  323. ``fvec``
  324. The function evaluated at the output
  325. ``fjac``
  326. A permutation of the R matrix of a QR
  327. factorization of the final approximate
  328. Jacobian matrix, stored column wise.
  329. Together with ipvt, the covariance of the
  330. estimate can be approximated.
  331. ``ipvt``
  332. An integer array of length N which defines
  333. a permutation matrix, p, such that
  334. fjac*p = q*r, where r is upper triangular
  335. with diagonal elements of nonincreasing
  336. magnitude. Column j of p is column ipvt(j)
  337. of the identity matrix.
  338. ``qtf``
  339. The vector (transpose(q) * fvec).
  340. Only returned if `full_output` is ``True``.
  341. mesg : str
  342. A string message giving information about the cause of failure.
  343. Only returned if `full_output` is ``True``.
  344. ier : int
  345. An integer flag. If it is equal to 1, 2, 3 or 4, the solution was
  346. found. Otherwise, the solution was not found. In either case, the
  347. optional output variable 'mesg' gives more information.
  348. See Also
  349. --------
  350. least_squares : Newer interface to solve nonlinear least-squares problems
  351. with bounds on the variables. See ``method='lm'`` in particular.
  352. Notes
  353. -----
  354. "leastsq" is a wrapper around MINPACK's lmdif and lmder algorithms.
  355. cov_x is a Jacobian approximation to the Hessian of the least squares
  356. objective function.
  357. This approximation assumes that the objective function is based on the
  358. difference between some observed target data (ydata) and a (non-linear)
  359. function of the parameters `f(xdata, params)` ::
  360. func(params) = ydata - f(xdata, params)
  361. so that the objective function is ::
  362. min sum((ydata - f(xdata, params))**2, axis=0)
  363. params
  364. The solution, `x`, is always a 1-D array, regardless of the shape of `x0`,
  365. or whether `x0` is a scalar.
  366. Examples
  367. --------
  368. >>> from scipy.optimize import leastsq
  369. >>> def func(x):
  370. ... return 2*(x-3)**2+1
  371. >>> leastsq(func, 0)
  372. (array([2.99999999]), 1)
  373. """
  374. x0 = asarray(x0).flatten()
  375. n = len(x0)
  376. if not isinstance(args, tuple):
  377. args = (args,)
  378. shape, dtype = _check_func('leastsq', 'func', func, x0, args, n)
  379. m = shape[0]
  380. if n > m:
  381. raise TypeError(f"Improper input: func input vector length N={n} must"
  382. f" not exceed func output vector length M={m}")
  383. if epsfcn is None:
  384. epsfcn = finfo(dtype).eps
  385. if Dfun is None:
  386. if maxfev == 0:
  387. maxfev = 200*(n + 1)
  388. retval = _minpack._lmdif(func, x0, args, full_output, ftol, xtol,
  389. gtol, maxfev, epsfcn, factor, diag)
  390. else:
  391. if col_deriv:
  392. _check_func('leastsq', 'Dfun', Dfun, x0, args, n, (n, m))
  393. else:
  394. _check_func('leastsq', 'Dfun', Dfun, x0, args, n, (m, n))
  395. if maxfev == 0:
  396. maxfev = 100 * (n + 1)
  397. retval = _minpack._lmder(func, Dfun, x0, args, full_output,
  398. col_deriv, ftol, xtol, gtol, maxfev,
  399. factor, diag)
  400. errors = {0: ["Improper input parameters.", TypeError],
  401. 1: ["Both actual and predicted relative reductions "
  402. f"in the sum of squares\n are at most {ftol:f}", None],
  403. 2: ["The relative error between two consecutive "
  404. f"iterates is at most {xtol:f}", None],
  405. 3: ["Both actual and predicted relative reductions in "
  406. f"the sum of squares\n are at most {ftol:f} and the "
  407. "relative error between two consecutive "
  408. f"iterates is at \n most {xtol:f}", None],
  409. 4: ["The cosine of the angle between func(x) and any "
  410. f"column of the\n Jacobian is at most {gtol:f} in "
  411. "absolute value", None],
  412. 5: ["Number of calls to function has reached "
  413. f"maxfev = {maxfev}.", ValueError],
  414. 6: [f"ftol={ftol:f} is too small, no further reduction "
  415. "in the sum of squares\n is possible.",
  416. ValueError],
  417. 7: [f"xtol={xtol:f} is too small, no further improvement in "
  418. "the approximate\n solution is possible.",
  419. ValueError],
  420. 8: [f"gtol={gtol:f} is too small, func(x) is orthogonal to the "
  421. "columns of\n the Jacobian to machine precision.", ValueError]}
  422. # The FORTRAN return value (possible return values are >= 0 and <= 8)
  423. info = retval[-1]
  424. if full_output:
  425. cov_x = None
  426. if info in LEASTSQ_SUCCESS:
  427. # This was
  428. # perm = take(eye(n), retval[1]['ipvt'] - 1, 0)
  429. # r = triu(transpose(retval[1]['fjac'])[:n, :])
  430. # R = dot(r, perm)
  431. # cov_x = inv(dot(transpose(R), R))
  432. # but the explicit dot product was not necessary and sometimes
  433. # the result was not symmetric positive definite. See gh-4555.
  434. perm = retval[1]['ipvt']
  435. n = len(perm)
  436. r = triu(transpose(retval[1]['fjac'])[:n, :])
  437. inv_triu = linalg.get_lapack_funcs('trtri', (r,))
  438. try:
  439. # inverse of permuted matrix is a permutation of matrix inverse
  440. invR, trtri_info = inv_triu(r) # default: upper, non-unit diag
  441. if trtri_info != 0: # explicit comparison for readability
  442. raise LinAlgError(f'trtri returned info {trtri_info}')
  443. invR[perm] = invR.copy()
  444. cov_x = invR @ invR.T
  445. except (LinAlgError, ValueError):
  446. pass
  447. return (retval[0], cov_x) + retval[1:-1] + (errors[info][0], info)
  448. else:
  449. if info in LEASTSQ_FAILURE:
  450. warnings.warn(errors[info][0], RuntimeWarning, stacklevel=2)
  451. elif info == 0:
  452. raise errors[info][1](errors[info][0])
  453. return retval[0], info
  454. def _lightweight_memoizer(f):
  455. # very shallow memoization to address gh-13670: only remember the first set
  456. # of parameters and corresponding function value, and only attempt to use
  457. # them twice (the number of times the function is evaluated at x0).
  458. def _memoized_func(params):
  459. if _memoized_func.skip_lookup:
  460. return f(params)
  461. if np.all(_memoized_func.last_params == params):
  462. return _memoized_func.last_val
  463. elif _memoized_func.last_params is not None:
  464. _memoized_func.skip_lookup = True
  465. val = f(params)
  466. if _memoized_func.last_params is None:
  467. _memoized_func.last_params = np.copy(params)
  468. _memoized_func.last_val = val
  469. return val
  470. _memoized_func.last_params = None
  471. _memoized_func.last_val = None
  472. _memoized_func.skip_lookup = False
  473. return _memoized_func
  474. def _wrap_func(func, xdata, ydata, transform):
  475. if transform is None:
  476. def func_wrapped(params):
  477. return func(xdata, *params) - ydata
  478. elif transform.size == 1 or transform.ndim == 1:
  479. def func_wrapped(params):
  480. return transform * (func(xdata, *params) - ydata)
  481. else:
  482. # Chisq = (y - yd)^T C^{-1} (y-yd)
  483. # transform = L such that C = L L^T
  484. # C^{-1} = L^{-T} L^{-1}
  485. # Chisq = (y - yd)^T L^{-T} L^{-1} (y-yd)
  486. # Define (y-yd)' = L^{-1} (y-yd)
  487. # by solving
  488. # L (y-yd)' = (y-yd)
  489. # and minimize (y-yd)'^T (y-yd)'
  490. def func_wrapped(params):
  491. return solve_triangular(transform, func(xdata, *params) - ydata, lower=True)
  492. return func_wrapped
  493. def _wrap_jac(jac, xdata, transform):
  494. if transform is None:
  495. def jac_wrapped(params):
  496. return jac(xdata, *params)
  497. elif transform.ndim == 1:
  498. def jac_wrapped(params):
  499. return transform[:, np.newaxis] * np.asarray(jac(xdata, *params))
  500. else:
  501. def jac_wrapped(params):
  502. return solve_triangular(transform,
  503. np.asarray(jac(xdata, *params)),
  504. lower=True)
  505. return jac_wrapped
  506. def _initialize_feasible(lb, ub):
  507. p0 = np.ones_like(lb)
  508. lb_finite = np.isfinite(lb)
  509. ub_finite = np.isfinite(ub)
  510. mask = lb_finite & ub_finite
  511. p0[mask] = 0.5 * (lb[mask] + ub[mask])
  512. mask = lb_finite & ~ub_finite
  513. p0[mask] = lb[mask] + 1
  514. mask = ~lb_finite & ub_finite
  515. p0[mask] = ub[mask] - 1
  516. return p0
  517. def curve_fit(f, xdata, ydata, p0=None, sigma=None, absolute_sigma=False,
  518. check_finite=None, bounds=(-np.inf, np.inf), method=None,
  519. jac=None, *, full_output=False, nan_policy=None,
  520. **kwargs):
  521. """
  522. Use non-linear least squares to fit a function, f, to data.
  523. Assumes ``ydata = f(xdata, *params) + eps``.
  524. Parameters
  525. ----------
  526. f : callable
  527. The model function, f(x, ...). It must take the independent
  528. variable as the first argument and the parameters to fit as
  529. separate remaining arguments.
  530. xdata : array_like
  531. The independent variable where the data is measured.
  532. Should usually be an M-length sequence or an (k,M)-shaped array for
  533. functions with k predictors, and each element should be float
  534. convertible if it is an array like object.
  535. ydata : array_like
  536. The dependent data, a length M array - nominally ``f(xdata, ...)``.
  537. p0 : array_like, optional
  538. Initial guess for the parameters (length N). If None, then the
  539. initial values will all be 1 (if the number of parameters for the
  540. function can be determined using introspection, otherwise a
  541. ValueError is raised).
  542. sigma : None or scalar or M-length sequence or MxM array, optional
  543. Determines the uncertainty in `ydata`. If we define residuals as
  544. ``r = ydata - f(xdata, *popt)``, then the interpretation of `sigma`
  545. depends on its number of dimensions:
  546. - A scalar or 1-D `sigma` should contain values of standard deviations of
  547. errors in `ydata`. In this case, the optimized function is
  548. ``chisq = sum((r / sigma) ** 2)``.
  549. - A 2-D `sigma` should contain the covariance matrix of
  550. errors in `ydata`. In this case, the optimized function is
  551. ``chisq = r.T @ inv(sigma) @ r``.
  552. .. versionadded:: 0.19
  553. None (default) is equivalent of 1-D `sigma` filled with ones.
  554. absolute_sigma : bool, optional
  555. If True, `sigma` is used in an absolute sense and the estimated parameter
  556. covariance `pcov` reflects these absolute values.
  557. If False (default), only the relative magnitudes of the `sigma` values matter.
  558. The returned parameter covariance matrix `pcov` is based on scaling
  559. `sigma` by a constant factor. This constant is set by demanding that the
  560. reduced `chisq` for the optimal parameters `popt` when using the
  561. *scaled* `sigma` equals unity. In other words, `sigma` is scaled to
  562. match the sample variance of the residuals after the fit. Default is False.
  563. Mathematically,
  564. ``pcov(absolute_sigma=False) = pcov(absolute_sigma=True) * chisq(popt)/(M-N)``
  565. check_finite : bool, optional
  566. If True, check that the input arrays do not contain nans of infs,
  567. and raise a ValueError if they do. Setting this parameter to
  568. False may silently produce nonsensical results if the input arrays
  569. do contain nans. Default is True if `nan_policy` is not specified
  570. explicitly and False otherwise.
  571. bounds : 2-tuple of array_like or `Bounds`, optional
  572. Lower and upper bounds on parameters. Defaults to no bounds.
  573. There are two ways to specify the bounds:
  574. - Instance of `Bounds` class.
  575. - 2-tuple of array_like: Each element of the tuple must be either
  576. an array with the length equal to the number of parameters, or a
  577. scalar (in which case the bound is taken to be the same for all
  578. parameters). Use ``np.inf`` with an appropriate sign to disable
  579. bounds on all or some parameters.
  580. method : {'lm', 'trf', 'dogbox'}, optional
  581. Method to use for optimization. See `least_squares` for more details.
  582. Default is 'lm' for unconstrained problems and 'trf' if `bounds` are
  583. provided. The method 'lm' won't work when the number of observations
  584. is less than the number of variables, use 'trf' or 'dogbox' in this
  585. case.
  586. .. versionadded:: 0.17
  587. jac : callable, string or None, optional
  588. Function with signature ``jac(x, ...)`` which computes the Jacobian
  589. matrix of the model function with respect to parameters as a dense
  590. array_like structure. It will be scaled according to provided `sigma`.
  591. If None (default), the Jacobian will be estimated numerically.
  592. String keywords for 'trf' and 'dogbox' methods can be used to select
  593. a finite difference scheme, see `least_squares`.
  594. .. versionadded:: 0.18
  595. full_output : boolean, optional
  596. If True, this function returns additional information: `infodict`,
  597. `mesg`, and `ier`.
  598. .. versionadded:: 1.9
  599. nan_policy : {'raise', 'omit', None}, optional
  600. Defines how to handle when input contains nan.
  601. The following options are available (default is None):
  602. * 'raise': throws an error
  603. * 'omit': performs the calculations ignoring nan values
  604. * None: no special handling of NaNs is performed
  605. (except what is done by check_finite); the behavior when NaNs
  606. are present is implementation-dependent and may change.
  607. Note that if this value is specified explicitly (not None),
  608. `check_finite` will be set as False.
  609. .. versionadded:: 1.11
  610. **kwargs
  611. Keyword arguments passed to `leastsq` for ``method='lm'`` or
  612. `least_squares` otherwise.
  613. Returns
  614. -------
  615. popt : array
  616. Optimal values for the parameters so that the sum of the squared
  617. residuals of ``f(xdata, *popt) - ydata`` is minimized.
  618. pcov : 2-D array
  619. The estimated approximate covariance of popt. The diagonals provide
  620. the variance of the parameter estimate. To compute one standard
  621. deviation errors on the parameters, use
  622. ``perr = np.sqrt(np.diag(pcov))``. Note that the relationship between
  623. `cov` and parameter error estimates is derived based on a linear
  624. approximation to the model function around the optimum [1]_.
  625. When this approximation becomes inaccurate, `cov` may not provide an
  626. accurate measure of uncertainty.
  627. How the `sigma` parameter affects the estimated covariance
  628. depends on `absolute_sigma` argument, as described above.
  629. If the Jacobian matrix at the solution doesn't have a full rank, then
  630. 'lm' method returns a matrix filled with ``np.inf``, on the other hand
  631. 'trf' and 'dogbox' methods use Moore-Penrose pseudoinverse to compute
  632. the covariance matrix. Covariance matrices with large condition numbers
  633. (e.g. computed with `numpy.linalg.cond`) may indicate that results are
  634. unreliable.
  635. infodict : dict (returned only if `full_output` is True)
  636. a dictionary of optional outputs with the keys:
  637. ``nfev``
  638. The number of function calls. Methods 'trf' and 'dogbox' do not
  639. count function calls for numerical Jacobian approximation,
  640. as opposed to 'lm' method.
  641. ``fvec``
  642. The residual values evaluated at the solution, for a 1-D `sigma`
  643. this is ``(f(x, *popt) - ydata)/sigma``.
  644. ``fjac``
  645. A permutation of the R matrix of a QR
  646. factorization of the final approximate
  647. Jacobian matrix, stored column wise.
  648. Together with ipvt, the covariance of the
  649. estimate can be approximated.
  650. Method 'lm' only provides this information.
  651. ``ipvt``
  652. An integer array of length N which defines
  653. a permutation matrix, p, such that
  654. fjac*p = q*r, where r is upper triangular
  655. with diagonal elements of nonincreasing
  656. magnitude. Column j of p is column ipvt(j)
  657. of the identity matrix.
  658. Method 'lm' only provides this information.
  659. ``qtf``
  660. The vector (transpose(q) * fvec).
  661. Method 'lm' only provides this information.
  662. .. versionadded:: 1.9
  663. mesg : str (returned only if `full_output` is True)
  664. A string message giving information about the solution.
  665. .. versionadded:: 1.9
  666. ier : int (returned only if `full_output` is True)
  667. An integer flag. If it is equal to 1, 2, 3 or 4, the solution was
  668. found. Otherwise, the solution was not found. In either case, the
  669. optional output variable `mesg` gives more information.
  670. .. versionadded:: 1.9
  671. Raises
  672. ------
  673. ValueError
  674. if either `ydata` or `xdata` contain NaNs, or if incompatible options
  675. are used.
  676. RuntimeError
  677. if the least-squares minimization fails.
  678. OptimizeWarning
  679. if covariance of the parameters can not be estimated.
  680. See Also
  681. --------
  682. least_squares : Minimize the sum of squares of nonlinear functions.
  683. scipy.stats.linregress : Calculate a linear least squares regression for
  684. two sets of measurements.
  685. Notes
  686. -----
  687. Users should ensure that inputs `xdata`, `ydata`, and the output of `f`
  688. are ``float64``, or else the optimization may return incorrect results.
  689. With ``method='lm'``, the algorithm uses the Levenberg-Marquardt algorithm
  690. through `leastsq`. Note that this algorithm can only deal with
  691. unconstrained problems.
  692. Box constraints can be handled by methods 'trf' and 'dogbox'. Refer to
  693. the docstring of `least_squares` for more information.
  694. Parameters to be fitted must have similar scale. Differences of multiple
  695. orders of magnitude can lead to incorrect results. For the 'trf' and
  696. 'dogbox' methods, the `x_scale` keyword argument can be used to scale
  697. the parameters.
  698. `curve_fit` is for local optimization of parameters to minimize the sum of squares
  699. of residuals. For global optimization, other choices of objective function, and
  700. other advanced features, consider using SciPy's :ref:`tutorial_optimize_global`
  701. tools or the `LMFIT <https://lmfit.github.io/lmfit-py/index.html>`_ package.
  702. References
  703. ----------
  704. .. [1] K. Vugrin et al. Confidence region estimation techniques for nonlinear
  705. regression in groundwater flow: Three case studies. Water Resources
  706. Research, Vol. 43, W03423, :doi:`10.1029/2005WR004804`
  707. Examples
  708. --------
  709. >>> import numpy as np
  710. >>> import matplotlib.pyplot as plt
  711. >>> from scipy.optimize import curve_fit
  712. >>> def func(x, a, b, c):
  713. ... return a * np.exp(-b * x) + c
  714. Define the data to be fit with some noise:
  715. >>> xdata = np.linspace(0, 4, 50)
  716. >>> y = func(xdata, 2.5, 1.3, 0.5)
  717. >>> rng = np.random.default_rng()
  718. >>> y_noise = 0.2 * rng.normal(size=xdata.size)
  719. >>> ydata = y + y_noise
  720. >>> plt.plot(xdata, ydata, 'b-', label='data')
  721. Fit for the parameters a, b, c of the function `func`:
  722. >>> popt, pcov = curve_fit(func, xdata, ydata)
  723. >>> popt
  724. array([2.56274217, 1.37268521, 0.47427475])
  725. >>> plt.plot(xdata, func(xdata, *popt), 'r-',
  726. ... label='fit: a=%5.3f, b=%5.3f, c=%5.3f' % tuple(popt))
  727. Constrain the optimization to the region of ``0 <= a <= 3``,
  728. ``0 <= b <= 1`` and ``0 <= c <= 0.5``:
  729. >>> popt, pcov = curve_fit(func, xdata, ydata, bounds=(0, [3., 1., 0.5]))
  730. >>> popt
  731. array([2.43736712, 1. , 0.34463856])
  732. >>> plt.plot(xdata, func(xdata, *popt), 'g--',
  733. ... label='fit: a=%5.3f, b=%5.3f, c=%5.3f' % tuple(popt))
  734. >>> plt.xlabel('x')
  735. >>> plt.ylabel('y')
  736. >>> plt.legend()
  737. >>> plt.show()
  738. For reliable results, the model `func` should not be overparametrized;
  739. redundant parameters can cause unreliable covariance matrices and, in some
  740. cases, poorer quality fits. As a quick check of whether the model may be
  741. overparameterized, calculate the condition number of the covariance matrix:
  742. >>> np.linalg.cond(pcov)
  743. 34.571092161547405 # may vary
  744. The value is small, so it does not raise much concern. If, however, we were
  745. to add a fourth parameter ``d`` to `func` with the same effect as ``a``:
  746. >>> def func2(x, a, b, c, d):
  747. ... return a * d * np.exp(-b * x) + c # a and d are redundant
  748. >>> popt, pcov = curve_fit(func2, xdata, ydata)
  749. >>> np.linalg.cond(pcov)
  750. 1.13250718925596e+32 # may vary
  751. Such a large value is cause for concern. The diagonal elements of the
  752. covariance matrix, which is related to uncertainty of the fit, gives more
  753. information:
  754. >>> np.diag(pcov)
  755. array([1.48814742e+29, 3.78596560e-02, 5.39253738e-03, 2.76417220e+28]) # may vary
  756. Note that the first and last terms are much larger than the other elements,
  757. suggesting that the optimal values of these parameters are ambiguous and
  758. that only one of these parameters is needed in the model.
  759. If the optimal parameters of `f` differ by multiple orders of magnitude, the
  760. resulting fit can be inaccurate. Sometimes, `curve_fit` can fail to find any
  761. results:
  762. >>> ydata = func(xdata, 500000, 0.01, 15)
  763. >>> try:
  764. ... popt, pcov = curve_fit(func, xdata, ydata, method = 'trf')
  765. ... except RuntimeError as e:
  766. ... print(e)
  767. Optimal parameters not found: The maximum number of function evaluations is
  768. exceeded.
  769. If parameter scale is roughly known beforehand, it can be defined in
  770. `x_scale` argument:
  771. >>> popt, pcov = curve_fit(func, xdata, ydata, method = 'trf',
  772. ... x_scale = [1000, 1, 1])
  773. >>> popt
  774. array([5.00000000e+05, 1.00000000e-02, 1.49999999e+01])
  775. """
  776. if p0 is None:
  777. # determine number of parameters by inspecting the function
  778. sig = _getfullargspec(f)
  779. args = sig.args
  780. if len(args) < 2:
  781. raise ValueError("Unable to determine number of fit parameters.")
  782. n = len(args) - 1
  783. else:
  784. p0 = np.atleast_1d(p0)
  785. n = p0.size
  786. if isinstance(bounds, Bounds):
  787. lb, ub = bounds.lb, bounds.ub
  788. else:
  789. lb, ub = prepare_bounds(bounds, n)
  790. if p0 is None:
  791. p0 = _initialize_feasible(lb, ub)
  792. bounded_problem = np.any((lb > -np.inf) | (ub < np.inf))
  793. if method is None:
  794. if bounded_problem:
  795. method = 'trf'
  796. else:
  797. method = 'lm'
  798. if method == 'lm' and bounded_problem:
  799. raise ValueError("Method 'lm' only works for unconstrained problems. "
  800. "Use 'trf' or 'dogbox' instead.")
  801. if check_finite is None:
  802. check_finite = True if nan_policy is None else False
  803. # optimization may produce garbage for float32 inputs, cast them to float64
  804. if check_finite:
  805. ydata = np.asarray_chkfinite(ydata, float)
  806. else:
  807. ydata = np.asarray(ydata, float)
  808. if isinstance(xdata, list | tuple | np.ndarray):
  809. # `xdata` is passed straight to the user-defined `f`, so allow
  810. # non-array_like `xdata`.
  811. if check_finite:
  812. xdata = np.asarray_chkfinite(xdata, float)
  813. else:
  814. xdata = np.asarray(xdata, float)
  815. if ydata.size == 0:
  816. raise ValueError("`ydata` must not be empty!")
  817. # nan handling is needed only if check_finite is False because if True,
  818. # the x-y data are already checked, and they don't contain nans.
  819. if not check_finite and nan_policy is not None:
  820. if nan_policy == "propagate":
  821. msg = "`nan_policy='propagate'` is not supported by this function."
  822. raise ValueError(msg)
  823. if nan_policy not in ("raise", "omit"):
  824. # Override error message raised by _contains_nan
  825. msg = "nan_policy must be one of {None, 'raise', 'omit'}"
  826. raise ValueError(msg)
  827. x_contains_nan = _contains_nan(xdata, nan_policy)
  828. y_contains_nan = _contains_nan(ydata, nan_policy)
  829. if (x_contains_nan or y_contains_nan) and nan_policy == 'omit':
  830. # ignore NaNs for N dimensional arrays
  831. has_nan = np.isnan(xdata)
  832. has_nan = has_nan.any(axis=tuple(range(has_nan.ndim-1)))
  833. has_nan |= np.isnan(ydata)
  834. xdata = xdata[..., ~has_nan]
  835. ydata = ydata[~has_nan]
  836. # Also omit the corresponding entries from sigma
  837. if sigma is not None:
  838. sigma = np.asarray(sigma)
  839. if sigma.ndim == 1:
  840. sigma = sigma[~has_nan]
  841. elif sigma.ndim == 2:
  842. sigma = sigma[~has_nan, :]
  843. sigma = sigma[:, ~has_nan]
  844. # Determine type of sigma
  845. if sigma is not None:
  846. sigma = np.asarray(sigma)
  847. # if 1-D or a scalar, sigma are errors, define transform = 1/sigma
  848. if sigma.size == 1 or sigma.shape == (ydata.size,):
  849. transform = 1.0 / sigma
  850. # if 2-D, sigma is the covariance matrix,
  851. # define transform = L such that L L^T = C
  852. elif sigma.shape == (ydata.size, ydata.size):
  853. try:
  854. # scipy.linalg.cholesky requires lower=True to return L L^T = A
  855. transform = cholesky(sigma, lower=True)
  856. except LinAlgError as e:
  857. raise ValueError("`sigma` must be positive definite.") from e
  858. else:
  859. raise ValueError("`sigma` has incorrect shape.")
  860. else:
  861. transform = None
  862. func = _lightweight_memoizer(_wrap_func(f, xdata, ydata, transform))
  863. if callable(jac):
  864. jac = _lightweight_memoizer(_wrap_jac(jac, xdata, transform))
  865. elif jac is None and method != 'lm':
  866. jac = '2-point'
  867. if 'args' in kwargs:
  868. # The specification for the model function `f` does not support
  869. # additional arguments. Refer to the `curve_fit` docstring for
  870. # acceptable call signatures of `f`.
  871. raise ValueError("'args' is not a supported keyword argument.")
  872. if method == 'lm':
  873. # if ydata.size == 1, this might be used for broadcast.
  874. if ydata.size != 1 and n > ydata.size:
  875. raise TypeError(f"The number of func parameters={n} must not"
  876. f" exceed the number of data points={ydata.size}")
  877. res = leastsq(func, p0, Dfun=jac, full_output=1, **kwargs)
  878. popt, pcov, infodict, errmsg, ier = res
  879. ysize = len(infodict['fvec'])
  880. cost = np.sum(infodict['fvec'] ** 2)
  881. if ier not in [1, 2, 3, 4]:
  882. raise RuntimeError("Optimal parameters not found: " + errmsg)
  883. else:
  884. # Rename maxfev (leastsq) to max_nfev (least_squares), if specified.
  885. if 'max_nfev' not in kwargs:
  886. kwargs['max_nfev'] = kwargs.pop('maxfev', None)
  887. res = least_squares(func, p0, jac=jac, bounds=bounds, method=method,
  888. **kwargs)
  889. if not res.success:
  890. raise RuntimeError("Optimal parameters not found: " + res.message)
  891. infodict = dict(nfev=res.nfev, fvec=res.fun)
  892. ier = res.status
  893. errmsg = res.message
  894. ysize = len(res.fun)
  895. cost = 2 * res.cost # res.cost is half sum of squares!
  896. popt = res.x
  897. # Do Moore-Penrose inverse discarding zero singular values.
  898. _, s, VT = svd(res.jac, full_matrices=False)
  899. threshold = np.finfo(float).eps * max(res.jac.shape) * s[0]
  900. s = s[s > threshold]
  901. VT = VT[:s.size]
  902. pcov = np.dot(VT.T / s**2, VT)
  903. warn_cov = False
  904. if pcov is None or np.isnan(pcov).any():
  905. # indeterminate covariance
  906. pcov = zeros((len(popt), len(popt)), dtype=float)
  907. pcov.fill(inf)
  908. warn_cov = True
  909. elif not absolute_sigma:
  910. if ysize > p0.size:
  911. s_sq = cost / (ysize - p0.size)
  912. pcov = pcov * s_sq
  913. else:
  914. pcov.fill(inf)
  915. warn_cov = True
  916. if warn_cov:
  917. warnings.warn('Covariance of the parameters could not be estimated',
  918. category=OptimizeWarning, stacklevel=2)
  919. if full_output:
  920. return popt, pcov, infodict, errmsg, ier
  921. else:
  922. return popt, pcov
  923. def check_gradient(fcn, Dfcn, x0, args=(), col_deriv=0):
  924. """Perform a simple check on the gradient for correctness.
  925. """
  926. x = atleast_1d(x0)
  927. n = len(x)
  928. x = x.reshape((n,))
  929. fvec = atleast_1d(fcn(x, *args))
  930. m = len(fvec)
  931. fvec = fvec.reshape((m,))
  932. ldfjac = m
  933. fjac = atleast_1d(Dfcn(x, *args))
  934. fjac = fjac.reshape((m, n))
  935. if col_deriv == 0:
  936. fjac = transpose(fjac)
  937. xp = zeros((n,), float)
  938. err = zeros((m,), float)
  939. fvecp = None
  940. _minpack._chkder(m, n, x, fvec, fjac, ldfjac, xp, fvecp, 1, err)
  941. fvecp = atleast_1d(fcn(xp, *args))
  942. fvecp = fvecp.reshape((m,))
  943. _minpack._chkder(m, n, x, fvec, fjac, ldfjac, xp, fvecp, 2, err)
  944. good = (prod(greater(err, 0.5), axis=0))
  945. return (good, err)
  946. def _del2(p0, p1, d):
  947. return p0 - np.square(p1 - p0) / d
  948. def _relerr(actual, desired):
  949. return (actual - desired) / desired
  950. def _fixed_point_helper(func, x0, args, xtol, maxiter, use_accel):
  951. p0 = x0
  952. for i in range(maxiter):
  953. p1 = func(p0, *args)
  954. if use_accel:
  955. p2 = func(p1, *args)
  956. d = p2 - 2.0 * p1 + p0
  957. p = xpx.apply_where(d != 0, (p0, p1, d), _del2, fill_value=p2)
  958. else:
  959. p = p1
  960. relerr = xpx.apply_where(p0 != 0, (p, p0), _relerr, fill_value=p)
  961. if np.all(np.abs(relerr) < xtol):
  962. return p
  963. p0 = p
  964. msg = f"Failed to converge after {maxiter} iterations, value is {p}"
  965. raise RuntimeError(msg)
  966. def fixed_point(func, x0, args=(), xtol=1e-8, maxiter=500, method='del2'):
  967. """
  968. Find a fixed point of the function.
  969. Given a function of one or more variables and a starting point, find a
  970. fixed point of the function: i.e., where ``func(x0) == x0``.
  971. Parameters
  972. ----------
  973. func : function
  974. Function to evaluate.
  975. x0 : array_like
  976. Fixed point of function.
  977. args : tuple, optional
  978. Extra arguments to `func`.
  979. xtol : float, optional
  980. Convergence tolerance, defaults to 1e-08.
  981. maxiter : int, optional
  982. Maximum number of iterations, defaults to 500.
  983. method : {"del2", "iteration"}, optional
  984. Method of finding the fixed-point, defaults to "del2",
  985. which uses Steffensen's Method with Aitken's ``Del^2``
  986. convergence acceleration [1]_. The "iteration" method simply iterates
  987. the function until convergence is detected, without attempting to
  988. accelerate the convergence.
  989. References
  990. ----------
  991. .. [1] Burden, Faires, "Numerical Analysis", 5th edition, pg. 80
  992. Examples
  993. --------
  994. >>> import numpy as np
  995. >>> from scipy import optimize
  996. >>> def func(x, c1, c2):
  997. ... return np.sqrt(c1/(x+c2))
  998. >>> c1 = np.array([10,12.])
  999. >>> c2 = np.array([3, 5.])
  1000. >>> optimize.fixed_point(func, [1.2, 1.3], args=(c1,c2))
  1001. array([ 1.4920333 , 1.37228132])
  1002. """
  1003. use_accel = {'del2': True, 'iteration': False}[method]
  1004. x0 = _asarray_validated(x0, as_inexact=True)
  1005. return _fixed_point_helper(func, x0, args, xtol, maxiter, use_accel)