_root.py 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733
  1. """
  2. Unified interfaces to root finding algorithms.
  3. Functions
  4. ---------
  5. - root : find a root of a vector function.
  6. """
  7. __all__ = ['root']
  8. import numpy as np
  9. from warnings import warn
  10. from ._optimize import MemoizeJac, OptimizeResult, _check_unknown_options
  11. from ._minpack_py import _root_hybr, leastsq
  12. from ._spectral import _root_df_sane
  13. from . import _nonlin as nonlin
  14. ROOT_METHODS = ['hybr', 'lm', 'broyden1', 'broyden2', 'anderson',
  15. 'linearmixing', 'diagbroyden', 'excitingmixing', 'krylov',
  16. 'df-sane']
  17. def root(fun, x0, args=(), method='hybr', jac=None, tol=None, callback=None,
  18. options=None):
  19. r"""
  20. Find a root of a vector function.
  21. Parameters
  22. ----------
  23. fun : callable
  24. A vector function to find a root of.
  25. Suppose the callable has signature ``f0(x, *my_args, **my_kwargs)``, where
  26. ``my_args`` and ``my_kwargs`` are required positional and keyword arguments.
  27. Rather than passing ``f0`` as the callable, wrap it to accept
  28. only ``x``; e.g., pass ``fun=lambda x: f0(x, *my_args, **my_kwargs)`` as the
  29. callable, where ``my_args`` (tuple) and ``my_kwargs`` (dict) have been
  30. gathered before invoking this function.
  31. x0 : ndarray
  32. Initial guess.
  33. args : tuple, optional
  34. Extra arguments passed to the objective function and its Jacobian.
  35. method : str, optional
  36. Type of solver. Should be one of
  37. - 'hybr' :ref:`(see here) <optimize.root-hybr>`
  38. - 'lm' :ref:`(see here) <optimize.root-lm>`
  39. - 'broyden1' :ref:`(see here) <optimize.root-broyden1>`
  40. - 'broyden2' :ref:`(see here) <optimize.root-broyden2>`
  41. - 'anderson' :ref:`(see here) <optimize.root-anderson>`
  42. - 'linearmixing' :ref:`(see here) <optimize.root-linearmixing>`
  43. - 'diagbroyden' :ref:`(see here) <optimize.root-diagbroyden>`
  44. - 'excitingmixing' :ref:`(see here) <optimize.root-excitingmixing>`
  45. - 'krylov' :ref:`(see here) <optimize.root-krylov>`
  46. - 'df-sane' :ref:`(see here) <optimize.root-dfsane>`
  47. jac : bool or callable, optional
  48. If `jac` is a Boolean and is True, `fun` is assumed to return the
  49. value of Jacobian along with the objective function. If False, the
  50. Jacobian will be estimated numerically.
  51. `jac` can also be a callable returning the Jacobian of `fun`. In
  52. this case, it must accept the same arguments as `fun`.
  53. tol : float, optional
  54. Tolerance for termination. For detailed control, use solver-specific
  55. options.
  56. callback : function, optional
  57. Optional callback function. It is called on every iteration as
  58. ``callback(x, f)`` where `x` is the current solution and `f`
  59. the corresponding residual. For all methods but 'hybr' and 'lm'.
  60. options : dict, optional
  61. A dictionary of solver options. E.g., `xtol` or `maxiter`, see
  62. :obj:`show_options()` for details.
  63. Returns
  64. -------
  65. sol : OptimizeResult
  66. The solution represented as a ``OptimizeResult`` object.
  67. Important attributes are: ``x`` the solution array, ``success`` a
  68. Boolean flag indicating if the algorithm exited successfully and
  69. ``message`` which describes the cause of the termination. See
  70. `OptimizeResult` for a description of other attributes.
  71. See also
  72. --------
  73. show_options : Additional options accepted by the solvers
  74. Notes
  75. -----
  76. This section describes the available solvers that can be selected by the
  77. 'method' parameter. The default method is *hybr*.
  78. Method *hybr* uses a modification of the Powell hybrid method as
  79. implemented in MINPACK [1]_.
  80. Method *lm* solves the system of nonlinear equations in a least squares
  81. sense using a modification of the Levenberg-Marquardt algorithm as
  82. implemented in MINPACK [1]_.
  83. Method *df-sane* is a derivative-free spectral method. [3]_
  84. Methods *broyden1*, *broyden2*, *anderson*, *linearmixing*,
  85. *diagbroyden*, *excitingmixing*, *krylov* are inexact Newton methods,
  86. with backtracking or full line searches [2]_. Each method corresponds
  87. to a particular Jacobian approximations.
  88. - Method *broyden1* uses Broyden's first Jacobian approximation, it is
  89. known as Broyden's good method.
  90. - Method *broyden2* uses Broyden's second Jacobian approximation, it
  91. is known as Broyden's bad method.
  92. - Method *anderson* uses (extended) Anderson mixing.
  93. - Method *Krylov* uses Krylov approximation for inverse Jacobian. It
  94. is suitable for large-scale problem.
  95. - Method *diagbroyden* uses diagonal Broyden Jacobian approximation.
  96. - Method *linearmixing* uses a scalar Jacobian approximation.
  97. - Method *excitingmixing* uses a tuned diagonal Jacobian
  98. approximation.
  99. .. warning::
  100. The algorithms implemented for methods *diagbroyden*,
  101. *linearmixing* and *excitingmixing* may be useful for specific
  102. problems, but whether they will work may depend strongly on the
  103. problem.
  104. .. versionadded:: 0.11.0
  105. References
  106. ----------
  107. .. [1] More, Jorge J., Burton S. Garbow, and Kenneth E. Hillstrom.
  108. 1980. User Guide for MINPACK-1.
  109. .. [2] C. T. Kelley. 1995. Iterative Methods for Linear and Nonlinear
  110. Equations. Society for Industrial and Applied Mathematics.
  111. :doi:`10.1137/1.9781611970944`
  112. .. [3] W. La Cruz, J.M. Martinez, M. Raydan.
  113. Math. Comp. 75, 1429 (2006).
  114. Examples
  115. --------
  116. The following functions define a system of nonlinear equations and its
  117. jacobian.
  118. >>> import numpy as np
  119. >>> def fun(x):
  120. ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
  121. ... 0.5 * (x[1] - x[0])**3 + x[1]]
  122. >>> def jac(x):
  123. ... return np.array([[1 + 1.5 * (x[0] - x[1])**2,
  124. ... -1.5 * (x[0] - x[1])**2],
  125. ... [-1.5 * (x[1] - x[0])**2,
  126. ... 1 + 1.5 * (x[1] - x[0])**2]])
  127. A solution can be obtained as follows.
  128. >>> from scipy import optimize
  129. >>> sol = optimize.root(fun, [0, 0], jac=jac, method='hybr')
  130. >>> sol.x
  131. array([ 0.8411639, 0.1588361])
  132. **Large problem**
  133. Suppose that we needed to solve the following integrodifferential
  134. equation on the square :math:`[0,1]\times[0,1]`:
  135. .. math::
  136. \nabla^2 P = 10 \left(\int_0^1\int_0^1\cosh(P)\,dx\,dy\right)^2
  137. with :math:`P(x,1) = 1` and :math:`P=0` elsewhere on the boundary of
  138. the square.
  139. The solution can be found using the ``method='krylov'`` solver:
  140. >>> from scipy import optimize
  141. >>> # parameters
  142. >>> nx, ny = 75, 75
  143. >>> hx, hy = 1./(nx-1), 1./(ny-1)
  144. >>> P_left, P_right = 0, 0
  145. >>> P_top, P_bottom = 1, 0
  146. >>> def residual(P):
  147. ... d2x = np.zeros_like(P)
  148. ... d2y = np.zeros_like(P)
  149. ...
  150. ... d2x[1:-1] = (P[2:] - 2*P[1:-1] + P[:-2]) / hx/hx
  151. ... d2x[0] = (P[1] - 2*P[0] + P_left)/hx/hx
  152. ... d2x[-1] = (P_right - 2*P[-1] + P[-2])/hx/hx
  153. ...
  154. ... d2y[:,1:-1] = (P[:,2:] - 2*P[:,1:-1] + P[:,:-2])/hy/hy
  155. ... d2y[:,0] = (P[:,1] - 2*P[:,0] + P_bottom)/hy/hy
  156. ... d2y[:,-1] = (P_top - 2*P[:,-1] + P[:,-2])/hy/hy
  157. ...
  158. ... return d2x + d2y - 10*np.cosh(P).mean()**2
  159. >>> guess = np.zeros((nx, ny), float)
  160. >>> sol = optimize.root(residual, guess, method='krylov')
  161. >>> print('Residual: %g' % abs(residual(sol.x)).max())
  162. Residual: 5.7972e-06 # may vary
  163. >>> import matplotlib.pyplot as plt
  164. >>> x, y = np.mgrid[0:1:(nx*1j), 0:1:(ny*1j)]
  165. >>> plt.pcolormesh(x, y, sol.x, shading='gouraud')
  166. >>> plt.colorbar()
  167. >>> plt.show()
  168. """
  169. def _wrapped_fun(*fargs):
  170. """
  171. Wrapped `func` to track the number of times
  172. the function has been called.
  173. """
  174. _wrapped_fun.nfev += 1
  175. return fun(*fargs)
  176. _wrapped_fun.nfev = 0
  177. if not isinstance(args, tuple):
  178. args = (args,)
  179. meth = method.lower()
  180. if options is None:
  181. options = {}
  182. if callback is not None and meth in ('hybr', 'lm'):
  183. warn(f'Method {method} does not accept callback.',
  184. RuntimeWarning, stacklevel=2)
  185. # fun also returns the Jacobian
  186. if not callable(jac) and meth in ('hybr', 'lm'):
  187. if bool(jac):
  188. fun = MemoizeJac(fun)
  189. jac = fun.derivative
  190. else:
  191. jac = None
  192. # set default tolerances
  193. if tol is not None:
  194. options = dict(options)
  195. if meth in ('hybr', 'lm'):
  196. options.setdefault('xtol', tol)
  197. elif meth in ('df-sane',):
  198. options.setdefault('ftol', tol)
  199. elif meth in ('broyden1', 'broyden2', 'anderson', 'linearmixing',
  200. 'diagbroyden', 'excitingmixing', 'krylov'):
  201. options.setdefault('xtol', tol)
  202. options.setdefault('xatol', np.inf)
  203. options.setdefault('ftol', np.inf)
  204. options.setdefault('fatol', np.inf)
  205. if meth == 'hybr':
  206. sol = _root_hybr(_wrapped_fun, x0, args=args, jac=jac, **options)
  207. elif meth == 'lm':
  208. sol = _root_leastsq(_wrapped_fun, x0, args=args, jac=jac, **options)
  209. elif meth == 'df-sane':
  210. _warn_jac_unused(jac, method)
  211. sol = _root_df_sane(_wrapped_fun, x0, args=args, callback=callback,
  212. **options)
  213. elif meth in ('broyden1', 'broyden2', 'anderson', 'linearmixing',
  214. 'diagbroyden', 'excitingmixing', 'krylov'):
  215. _warn_jac_unused(jac, method)
  216. sol = _root_nonlin_solve(_wrapped_fun, x0, args=args, jac=jac,
  217. _method=meth, _callback=callback,
  218. **options)
  219. else:
  220. raise ValueError(f'Unknown solver {method}')
  221. sol.nfev = _wrapped_fun.nfev
  222. return sol
  223. def _warn_jac_unused(jac, method):
  224. if jac is not None:
  225. warn(f'Method {method} does not use the jacobian (jac).',
  226. RuntimeWarning, stacklevel=2)
  227. def _root_leastsq(fun, x0, args=(), jac=None,
  228. col_deriv=0, xtol=1.49012e-08, ftol=1.49012e-08,
  229. gtol=0.0, maxiter=0, eps=0.0, factor=100, diag=None,
  230. **unknown_options):
  231. """
  232. Solve for least squares with Levenberg-Marquardt
  233. Options
  234. -------
  235. col_deriv : bool
  236. non-zero to specify that the Jacobian function computes derivatives
  237. down the columns (faster, because there is no transpose operation).
  238. ftol : float
  239. Relative error desired in the sum of squares.
  240. xtol : float
  241. Relative error desired in the approximate solution.
  242. gtol : float
  243. Orthogonality desired between the function vector and the columns
  244. of the Jacobian.
  245. maxiter : int
  246. The maximum number of calls to the function. If zero, then
  247. 100*(N+1) is the maximum where N is the number of elements in x0.
  248. eps : float
  249. A suitable step length for the forward-difference approximation of
  250. the Jacobian (for Dfun=None). If `eps` is less than the machine
  251. precision, it is assumed that the relative errors in the functions
  252. are of the order of the machine precision.
  253. factor : float
  254. A parameter determining the initial step bound
  255. (``factor * || diag * x||``). Should be in interval ``(0.1, 100)``.
  256. diag : sequence
  257. N positive entries that serve as a scale factors for the variables.
  258. """
  259. nfev = 0
  260. def _wrapped_fun(*fargs):
  261. """
  262. Wrapped `func` to track the number of times
  263. the function has been called.
  264. """
  265. nonlocal nfev
  266. nfev += 1
  267. return fun(*fargs)
  268. _check_unknown_options(unknown_options)
  269. x, cov_x, info, msg, ier = leastsq(_wrapped_fun, x0, args=args,
  270. Dfun=jac, full_output=True,
  271. col_deriv=col_deriv, xtol=xtol,
  272. ftol=ftol, gtol=gtol,
  273. maxfev=maxiter, epsfcn=eps,
  274. factor=factor, diag=diag)
  275. sol = OptimizeResult(x=x, message=msg, status=ier,
  276. success=ier in (1, 2, 3, 4), cov_x=cov_x,
  277. fun=info.pop('fvec'), method="lm")
  278. sol.update(info)
  279. sol.nfev = nfev
  280. return sol
  281. def _root_nonlin_solve(fun, x0, args=(), jac=None,
  282. _callback=None, _method=None,
  283. nit=None, disp=False, maxiter=None,
  284. ftol=None, fatol=None, xtol=None, xatol=None,
  285. tol_norm=None, line_search='armijo', jac_options=None,
  286. **unknown_options):
  287. _check_unknown_options(unknown_options)
  288. f_tol = fatol
  289. f_rtol = ftol
  290. x_tol = xatol
  291. x_rtol = xtol
  292. verbose = disp
  293. if jac_options is None:
  294. jac_options = dict()
  295. jacobian = {'broyden1': nonlin.BroydenFirst,
  296. 'broyden2': nonlin.BroydenSecond,
  297. 'anderson': nonlin.Anderson,
  298. 'linearmixing': nonlin.LinearMixing,
  299. 'diagbroyden': nonlin.DiagBroyden,
  300. 'excitingmixing': nonlin.ExcitingMixing,
  301. 'krylov': nonlin.KrylovJacobian
  302. }[_method]
  303. if args:
  304. if jac is True:
  305. def f(x):
  306. return fun(x, *args)[0]
  307. else:
  308. def f(x):
  309. return fun(x, *args)
  310. else:
  311. f = fun
  312. x, info = nonlin.nonlin_solve(f, x0, jacobian=jacobian(**jac_options),
  313. iter=nit, verbose=verbose,
  314. maxiter=maxiter, f_tol=f_tol,
  315. f_rtol=f_rtol, x_tol=x_tol,
  316. x_rtol=x_rtol, tol_norm=tol_norm,
  317. line_search=line_search,
  318. callback=_callback, full_output=True,
  319. raise_exception=False)
  320. sol = OptimizeResult(x=x, method=_method)
  321. sol.update(info)
  322. return sol
  323. def _root_broyden1_doc():
  324. """
  325. Options
  326. -------
  327. nit : int, optional
  328. Number of iterations to make. If omitted (default), make as many
  329. as required to meet tolerances.
  330. disp : bool, optional
  331. Print status to stdout on every iteration.
  332. maxiter : int, optional
  333. Maximum number of iterations to make.
  334. ftol : float, optional
  335. Relative tolerance for the residual. If omitted, not used.
  336. fatol : float, optional
  337. Absolute tolerance (in max-norm) for the residual.
  338. If omitted, default is 6e-6.
  339. xtol : float, optional
  340. Relative minimum step size. If omitted, not used.
  341. xatol : float, optional
  342. Absolute minimum step size, as determined from the Jacobian
  343. approximation. If the step size is smaller than this, optimization
  344. is terminated as successful. If omitted, not used.
  345. tol_norm : function(vector) -> scalar, optional
  346. Norm to use in convergence check. Default is the maximum norm.
  347. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  348. Which type of a line search to use to determine the step size in
  349. the direction given by the Jacobian approximation. Defaults to
  350. 'armijo'.
  351. jac_options : dict, optional
  352. Options for the respective Jacobian approximation.
  353. alpha : float, optional
  354. Initial guess for the Jacobian is (-1/alpha).
  355. reduction_method : str or tuple, optional
  356. Method used in ensuring that the rank of the Broyden
  357. matrix stays low. Can either be a string giving the
  358. name of the method, or a tuple of the form ``(method,
  359. param1, param2, ...)`` that gives the name of the
  360. method and values for additional parameters.
  361. Methods available:
  362. - ``restart``: drop all matrix columns. Has no extra parameters.
  363. - ``simple``: drop oldest matrix column. Has no extra parameters.
  364. - ``svd``: keep only the most significant SVD components.
  365. Takes an extra parameter, ``to_retain``, which determines the
  366. number of SVD components to retain when rank reduction is done.
  367. Default is ``max_rank - 2``.
  368. max_rank : int, optional
  369. Maximum rank for the Broyden matrix.
  370. Default is infinity (i.e., no rank reduction).
  371. Examples
  372. --------
  373. >>> def func(x):
  374. ... return np.cos(x) + x[::-1] - [1, 2, 3, 4]
  375. ...
  376. >>> from scipy import optimize
  377. >>> res = optimize.root(func, [1, 1, 1, 1], method='broyden1', tol=1e-14)
  378. >>> x = res.x
  379. >>> x
  380. array([4.04674914, 3.91158389, 2.71791677, 1.61756251])
  381. >>> np.cos(x) + x[::-1]
  382. array([1., 2., 3., 4.])
  383. """
  384. pass
  385. def _root_broyden2_doc():
  386. """
  387. Options
  388. -------
  389. nit : int, optional
  390. Number of iterations to make. If omitted (default), make as many
  391. as required to meet tolerances.
  392. disp : bool, optional
  393. Print status to stdout on every iteration.
  394. maxiter : int, optional
  395. Maximum number of iterations to make.
  396. ftol : float, optional
  397. Relative tolerance for the residual. If omitted, not used.
  398. fatol : float, optional
  399. Absolute tolerance (in max-norm) for the residual.
  400. If omitted, default is 6e-6.
  401. xtol : float, optional
  402. Relative minimum step size. If omitted, not used.
  403. xatol : float, optional
  404. Absolute minimum step size, as determined from the Jacobian
  405. approximation. If the step size is smaller than this, optimization
  406. is terminated as successful. If omitted, not used.
  407. tol_norm : function(vector) -> scalar, optional
  408. Norm to use in convergence check. Default is the maximum norm.
  409. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  410. Which type of a line search to use to determine the step size in
  411. the direction given by the Jacobian approximation. Defaults to
  412. 'armijo'.
  413. jac_options : dict, optional
  414. Options for the respective Jacobian approximation.
  415. alpha : float, optional
  416. Initial guess for the Jacobian is (-1/alpha).
  417. reduction_method : str or tuple, optional
  418. Method used in ensuring that the rank of the Broyden
  419. matrix stays low. Can either be a string giving the
  420. name of the method, or a tuple of the form ``(method,
  421. param1, param2, ...)`` that gives the name of the
  422. method and values for additional parameters.
  423. Methods available:
  424. - ``restart``: drop all matrix columns. Has no extra parameters.
  425. - ``simple``: drop oldest matrix column. Has no extra parameters.
  426. - ``svd``: keep only the most significant SVD components.
  427. Takes an extra parameter, ``to_retain``, which determines the
  428. number of SVD components to retain when rank reduction is done.
  429. Default is ``max_rank - 2``.
  430. max_rank : int, optional
  431. Maximum rank for the Broyden matrix.
  432. Default is infinity (i.e., no rank reduction).
  433. """
  434. pass
  435. def _root_anderson_doc():
  436. """
  437. Options
  438. -------
  439. nit : int, optional
  440. Number of iterations to make. If omitted (default), make as many
  441. as required to meet tolerances.
  442. disp : bool, optional
  443. Print status to stdout on every iteration.
  444. maxiter : int, optional
  445. Maximum number of iterations to make.
  446. ftol : float, optional
  447. Relative tolerance for the residual. If omitted, not used.
  448. fatol : float, optional
  449. Absolute tolerance (in max-norm) for the residual.
  450. If omitted, default is 6e-6.
  451. xtol : float, optional
  452. Relative minimum step size. If omitted, not used.
  453. xatol : float, optional
  454. Absolute minimum step size, as determined from the Jacobian
  455. approximation. If the step size is smaller than this, optimization
  456. is terminated as successful. If omitted, not used.
  457. tol_norm : function(vector) -> scalar, optional
  458. Norm to use in convergence check. Default is the maximum norm.
  459. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  460. Which type of a line search to use to determine the step size in
  461. the direction given by the Jacobian approximation. Defaults to
  462. 'armijo'.
  463. jac_options : dict, optional
  464. Options for the respective Jacobian approximation.
  465. alpha : float, optional
  466. Initial guess for the Jacobian is (-1/alpha).
  467. M : float, optional
  468. Number of previous vectors to retain. Defaults to 5.
  469. w0 : float, optional
  470. Regularization parameter for numerical stability.
  471. Compared to unity, good values of the order of 0.01.
  472. """
  473. pass
  474. def _root_linearmixing_doc():
  475. """
  476. Options
  477. -------
  478. nit : int, optional
  479. Number of iterations to make. If omitted (default), make as many
  480. as required to meet tolerances.
  481. disp : bool, optional
  482. Print status to stdout on every iteration.
  483. maxiter : int, optional
  484. Maximum number of iterations to make.
  485. ftol : float, optional
  486. Relative tolerance for the residual. If omitted, not used.
  487. fatol : float, optional
  488. Absolute tolerance (in max-norm) for the residual.
  489. If omitted, default is 6e-6.
  490. xtol : float, optional
  491. Relative minimum step size. If omitted, not used.
  492. xatol : float, optional
  493. Absolute minimum step size, as determined from the Jacobian
  494. approximation. If the step size is smaller than this, optimization
  495. is terminated as successful. If omitted, not used.
  496. tol_norm : function(vector) -> scalar, optional
  497. Norm to use in convergence check. Default is the maximum norm.
  498. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  499. Which type of a line search to use to determine the step size in
  500. the direction given by the Jacobian approximation. Defaults to
  501. 'armijo'.
  502. jac_options : dict, optional
  503. Options for the respective Jacobian approximation.
  504. alpha : float, optional
  505. initial guess for the jacobian is (-1/alpha).
  506. """
  507. pass
  508. def _root_diagbroyden_doc():
  509. """
  510. Options
  511. -------
  512. nit : int, optional
  513. Number of iterations to make. If omitted (default), make as many
  514. as required to meet tolerances.
  515. disp : bool, optional
  516. Print status to stdout on every iteration.
  517. maxiter : int, optional
  518. Maximum number of iterations to make.
  519. ftol : float, optional
  520. Relative tolerance for the residual. If omitted, not used.
  521. fatol : float, optional
  522. Absolute tolerance (in max-norm) for the residual.
  523. If omitted, default is 6e-6.
  524. xtol : float, optional
  525. Relative minimum step size. If omitted, not used.
  526. xatol : float, optional
  527. Absolute minimum step size, as determined from the Jacobian
  528. approximation. If the step size is smaller than this, optimization
  529. is terminated as successful. If omitted, not used.
  530. tol_norm : function(vector) -> scalar, optional
  531. Norm to use in convergence check. Default is the maximum norm.
  532. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  533. Which type of a line search to use to determine the step size in
  534. the direction given by the Jacobian approximation. Defaults to
  535. 'armijo'.
  536. jac_options : dict, optional
  537. Options for the respective Jacobian approximation.
  538. alpha : float, optional
  539. initial guess for the jacobian is (-1/alpha).
  540. """
  541. pass
  542. def _root_excitingmixing_doc():
  543. """
  544. Options
  545. -------
  546. nit : int, optional
  547. Number of iterations to make. If omitted (default), make as many
  548. as required to meet tolerances.
  549. disp : bool, optional
  550. Print status to stdout on every iteration.
  551. maxiter : int, optional
  552. Maximum number of iterations to make.
  553. ftol : float, optional
  554. Relative tolerance for the residual. If omitted, not used.
  555. fatol : float, optional
  556. Absolute tolerance (in max-norm) for the residual.
  557. If omitted, default is 6e-6.
  558. xtol : float, optional
  559. Relative minimum step size. If omitted, not used.
  560. xatol : float, optional
  561. Absolute minimum step size, as determined from the Jacobian
  562. approximation. If the step size is smaller than this, optimization
  563. is terminated as successful. If omitted, not used.
  564. tol_norm : function(vector) -> scalar, optional
  565. Norm to use in convergence check. Default is the maximum norm.
  566. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  567. Which type of a line search to use to determine the step size in
  568. the direction given by the Jacobian approximation. Defaults to
  569. 'armijo'.
  570. jac_options : dict, optional
  571. Options for the respective Jacobian approximation.
  572. alpha : float, optional
  573. Initial Jacobian approximation is (-1/alpha).
  574. alphamax : float, optional
  575. The entries of the diagonal Jacobian are kept in the range
  576. ``[alpha, alphamax]``.
  577. """
  578. pass
  579. def _root_krylov_doc():
  580. """
  581. Options
  582. -------
  583. nit : int, optional
  584. Number of iterations to make. If omitted (default), make as many
  585. as required to meet tolerances.
  586. disp : bool, optional
  587. Print status to stdout on every iteration.
  588. maxiter : int, optional
  589. Maximum number of iterations to make.
  590. ftol : float, optional
  591. Relative tolerance for the residual. If omitted, not used.
  592. fatol : float, optional
  593. Absolute tolerance (in max-norm) for the residual.
  594. If omitted, default is 6e-6.
  595. xtol : float, optional
  596. Relative minimum step size. If omitted, not used.
  597. xatol : float, optional
  598. Absolute minimum step size, as determined from the Jacobian
  599. approximation. If the step size is smaller than this, optimization
  600. is terminated as successful. If omitted, not used.
  601. tol_norm : function(vector) -> scalar, optional
  602. Norm to use in convergence check. Default is the maximum norm.
  603. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  604. Which type of a line search to use to determine the step size in
  605. the direction given by the Jacobian approximation. Defaults to
  606. 'armijo'.
  607. jac_options : dict, optional
  608. Options for the respective Jacobian approximation.
  609. rdiff : float, optional
  610. Relative step size to use in numerical differentiation.
  611. method : str or callable, optional
  612. Krylov method to use to approximate the Jacobian. Can be a string,
  613. or a function implementing the same interface as the iterative
  614. solvers in `scipy.sparse.linalg`. If a string, needs to be one of:
  615. ``'lgmres'``, ``'gmres'``, ``'bicgstab'``, ``'cgs'``, ``'minres'``,
  616. ``'tfqmr'``.
  617. The default is `scipy.sparse.linalg.lgmres`.
  618. inner_M : LinearOperator or InverseJacobian
  619. Preconditioner for the inner Krylov iteration.
  620. Note that you can use also inverse Jacobians as (adaptive)
  621. preconditioners. For example,
  622. >>> jac = BroydenFirst()
  623. >>> kjac = KrylovJacobian(inner_M=jac.inverse).
  624. If the preconditioner has a method named 'update', it will
  625. be called as ``update(x, f)`` after each nonlinear step,
  626. with ``x`` giving the current point, and ``f`` the current
  627. function value.
  628. inner_rtol, inner_atol, inner_callback, ...
  629. Parameters to pass on to the "inner" Krylov solver.
  630. For a full list of options, see the documentation for the
  631. solver you are using. By default this is `scipy.sparse.linalg.lgmres`.
  632. If the solver has been overridden through `method`, see the documentation
  633. for that solver instead.
  634. To use an option for that solver, prepend ``inner_`` to it.
  635. For example, to control the ``rtol`` argument to the solver,
  636. set the `inner_rtol` option here.
  637. outer_k : int, optional
  638. Size of the subspace kept across LGMRES nonlinear
  639. iterations.
  640. See `scipy.sparse.linalg.lgmres` for details.
  641. """
  642. pass