_root.py 28 KB

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