bundles / scipy 1.17.1 / scipy / optimize / _nonlin / newton_krylov
function
scipy.optimize._nonlin:newton_krylov
Signature
def newton_krylov ( F , xin , iter = None , rdiff = None , method = lgmres , inner_maxiter = 20 , inner_M = None , outer_k = 10 , verbose = False , maxiter = None , f_tol = None , f_rtol = None , x_tol = None , x_rtol = None , tol_norm = None , line_search = armijo , callback = None , ** kw ) Summary
Find a root of a function, using Krylov approximation for inverse Jacobian.
Extended Summary
This method is suitable for solving large-scale problems.
Parameters
F: function(x) -> fFunction whose root to find; should take and return an array-like object.
xin: array_likeInitial guess for the solution
rdiff: float, optionalRelative step size to use in numerical differentiation.
method: str or callable, optionalKrylov method to use to approximate the Jacobian. Can be a string, or a function implementing the same interface as the iterative solvers in scipy.sparse.linalg. If a string, needs to be one of:
'lgmres','gmres','bicgstab','cgs','minres','tfqmr'.The default is scipy.sparse.linalg.lgmres.
inner_maxiter: int, optionalParameter to pass to the "inner" Krylov solver: maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved.
inner_M: LinearOperator or InverseJacobianPreconditioner for the inner Krylov iteration. Note that you can use also inverse Jacobians as (adaptive) preconditioners. For example,
>>> from scipy.optimize import BroydenFirst, KrylovJacobian >>> from scipy.optimize import InverseJacobian >>> jac = BroydenFirst() >>> kjac = KrylovJacobian(inner_M=InverseJacobian(jac))
If the preconditioner has a method named 'update', it will be called as
update(x, f)after each nonlinear step, withxgiving the current point, andfthe current function value.outer_k: int, optionalSize of the subspace kept across LGMRES nonlinear iterations. See scipy.sparse.linalg.lgmres for details.
inner_kwargs: kwargsKeyword parameters for the "inner" Krylov solver (defined with
method). Parameter names must start with theinner_prefix which will be stripped before passing on the inner method. See, e.g., scipy.sparse.linalg.gmres for details.iter: int, optionalNumber of iterations to make. If omitted (default), make as many as required to meet tolerances.
verbose: bool, optionalPrint status to stdout on every iteration.
maxiter: int, optionalMaximum number of iterations to make. If more are needed to meet convergence, NoConvergence is raised.
f_tol: float, optionalAbsolute tolerance (in max-norm) for the residual. If omitted, default is 6e-6.
f_rtol: float, optionalRelative tolerance for the residual. If omitted, not used.
x_tol: float, optionalAbsolute minimum step size, as determined from the Jacobian approximation. If the step size is smaller than this, optimization is terminated as successful. If omitted, not used.
x_rtol: float, optionalRelative minimum step size. If omitted, not used.
tol_norm: function(vector) -> scalar, optionalNorm to use in convergence check. Default is the maximum norm.
line_search: {None, 'armijo' (default), 'wolfe'}, optionalWhich type of a line search to use to determine the step size in the direction given by the Jacobian approximation. Defaults to 'armijo'.
callback: function, optionalOptional callback function. It is called on every iteration as
callback(x, f)wherexis the current solution andfthe corresponding residual.
Returns
sol: ndarrayAn array (of similar array type as
x0) containing the final solution.
Raises
: NoConvergenceWhen a solution was not found.
Notes
This function implements a Newton-Krylov solver. The basic idea is to compute the inverse of the Jacobian with an iterative Krylov method. These methods require only evaluating the Jacobian-vector products, which are conveniently approximated by a finite difference:
Due to the use of iterative matrix inverses, these methods can deal with large nonlinear problems.
SciPy's scipy.sparse.linalg module offers a selection of Krylov solvers to choose from. The default here is lgmres, which is a variant of restarted GMRES iteration that reuses some of the information obtained in the previous Newton steps to invert Jacobians in subsequent steps.
For a review on Newton-Krylov methods, see for example [1], and for the LGMRES sparse inverse method, see [2].
Examples
The following functions define a system of nonlinear equationsdef fun(x): return [x[0] + 0.5 * x[1] - 1.0, 0.5 * (x[1] - x[0]) ** 2]✓
from scipy import optimize sol = optimize.newton_krylov(fun, [0, 0])✓
sol
✗See also
- root
Interface to root finding algorithms for multivariate functions. See
method='krylov'in particular.- scipy.sparse.linalg.gmres
- scipy.sparse.linalg.lgmres
Aliases
-
scipy.optimize.newton_krylov