KrylovJacobian#
- class scipy.optimize.KrylovJacobian(rdiff=None, method='lgmres', inner_maxiter=20, inner_M=None, outer_k=10, **kw)[source]#
Find a root of a function, using Krylov approximation for inverse Jacobian.
This method is suitable for solving large-scale problems.
- Parameters:
- %(params_basic)s
- rdifffloat, optional
Relative step size to use in numerical differentiation.
- methodstr or callable, optional
Krylov 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_maxiterint, optional
Parameter 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_MLinearOperator or InverseJacobian
Preconditioner 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, withx
giving the current point, andf
the current function value.- outer_kint, optional
Size of the subspace kept across LGMRES nonlinear iterations. See
scipy.sparse.linalg.lgmres
for details.- inner_kwargskwargs
Keyword parameters for the “inner” Krylov solver (defined with method). Parameter names must start with the inner_ prefix which will be stripped before passing on the inner method. See, e.g.,
scipy.sparse.linalg.gmres
for details.- %(params_extra)s
Methods
aspreconditioner
matvec
setup
solve
update
See also
root
Interface to root finding algorithms for multivariate functions. See
method='krylov'
in particular.scipy.sparse.linalg.gmres
scipy.sparse.linalg.lgmres
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:
\[J v \approx (f(x + \omega*v/|v|) - f(x)) / \omega\]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].
References
[1]C. T. Kelley, Solving Nonlinear Equations with Newton’s Method, SIAM, pp.57-83, 2003. DOI:10.1137/1.9780898718898.ch3
[2]D.A. Knoll and D.E. Keyes, J. Comp. Phys. 193, 357 (2004). DOI:10.1016/j.jcp.2003.08.010
[3]A.H. Baker and E.R. Jessup and T. Manteuffel, SIAM J. Matrix Anal. Appl. 26, 962 (2005). DOI:10.1137/S0895479803422014
Examples
The following functions define a system of nonlinear equations
>>> def fun(x): ... return [x[0] + 0.5 * x[1] - 1.0, ... 0.5 * (x[1] - x[0]) ** 2]
A solution can be obtained as follows.
>>> from scipy import optimize >>> sol = optimize.newton_krylov(fun, [0, 0]) >>> sol array([0.66731771, 0.66536458])