_nonlin.py 51 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648
  1. # Copyright (C) 2009, Pauli Virtanen <pav@iki.fi>
  2. # Distributed under the same license as SciPy.
  3. import inspect
  4. import sys
  5. import warnings
  6. import numpy as np
  7. from numpy import asarray, dot, vdot
  8. from scipy.linalg import norm, solve, inv, qr, svd, LinAlgError
  9. import scipy.sparse.linalg
  10. import scipy.sparse
  11. from scipy.linalg import get_blas_funcs
  12. from scipy._lib._util import copy_if_needed, _dedent_for_py313
  13. from scipy._lib._util import getfullargspec_no_self as _getfullargspec
  14. from ._linesearch import scalar_search_wolfe1, scalar_search_armijo
  15. from inspect import signature
  16. from difflib import get_close_matches
  17. from types import GenericAlias
  18. __all__ = [
  19. 'broyden1', 'broyden2', 'anderson', 'linearmixing',
  20. 'diagbroyden', 'excitingmixing', 'newton_krylov',
  21. 'BroydenFirst', 'KrylovJacobian', 'InverseJacobian', 'NoConvergence']
  22. #------------------------------------------------------------------------------
  23. # Utility functions
  24. #------------------------------------------------------------------------------
  25. class NoConvergence(Exception):
  26. """Exception raised when nonlinear solver fails to converge within the specified
  27. `maxiter`."""
  28. pass
  29. def maxnorm(x):
  30. return np.absolute(x).max()
  31. def _as_inexact(x):
  32. """Return `x` as an array, of either floats or complex floats"""
  33. x = asarray(x)
  34. if not np.issubdtype(x.dtype, np.inexact):
  35. return asarray(x, dtype=np.float64)
  36. return x
  37. def _array_like(x, x0):
  38. """Return ndarray `x` as same array subclass and shape as `x0`"""
  39. x = np.reshape(x, np.shape(x0))
  40. wrap = getattr(x0, '__array_wrap__', x.__array_wrap__)
  41. return wrap(x)
  42. def _safe_norm(v):
  43. if not np.isfinite(v).all():
  44. return np.array(np.inf)
  45. return norm(v)
  46. #------------------------------------------------------------------------------
  47. # Generic nonlinear solver machinery
  48. #------------------------------------------------------------------------------
  49. _doc_parts = dict(
  50. params_basic=_dedent_for_py313("""
  51. F : function(x) -> f
  52. Function whose root to find; should take and return an array-like
  53. object.
  54. xin : array_like
  55. Initial guess for the solution
  56. """).strip(),
  57. params_extra=_dedent_for_py313("""
  58. iter : int, optional
  59. Number of iterations to make. If omitted (default), make as many
  60. as required to meet tolerances.
  61. verbose : bool, optional
  62. Print status to stdout on every iteration.
  63. maxiter : int, optional
  64. Maximum number of iterations to make. If more are needed to
  65. meet convergence, `NoConvergence` is raised.
  66. f_tol : float, optional
  67. Absolute tolerance (in max-norm) for the residual.
  68. If omitted, default is 6e-6.
  69. f_rtol : float, optional
  70. Relative tolerance for the residual. If omitted, not used.
  71. x_tol : float, optional
  72. Absolute minimum step size, as determined from the Jacobian
  73. approximation. If the step size is smaller than this, optimization
  74. is terminated as successful. If omitted, not used.
  75. x_rtol : float, optional
  76. Relative minimum step size. If omitted, not used.
  77. tol_norm : function(vector) -> scalar, optional
  78. Norm to use in convergence check. Default is the maximum norm.
  79. line_search : {None, 'armijo' (default), 'wolfe'}, optional
  80. Which type of a line search to use to determine the step size in the
  81. direction given by the Jacobian approximation. Defaults to 'armijo'.
  82. callback : function, optional
  83. Optional callback function. It is called on every iteration as
  84. ``callback(x, f)`` where `x` is the current solution and `f`
  85. the corresponding residual.
  86. Returns
  87. -------
  88. sol : ndarray
  89. An array (of similar array type as `x0`) containing the final solution.
  90. Raises
  91. ------
  92. NoConvergence
  93. When a solution was not found.
  94. """).strip()
  95. )
  96. def _set_doc(obj):
  97. if obj.__doc__:
  98. obj.__doc__ = obj.__doc__ % _doc_parts
  99. def nonlin_solve(F, x0, jacobian='krylov', iter=None, verbose=False,
  100. maxiter=None, f_tol=None, f_rtol=None, x_tol=None, x_rtol=None,
  101. tol_norm=None, line_search='armijo', callback=None,
  102. full_output=False, raise_exception=True):
  103. """
  104. Find a root of a function, in a way suitable for large-scale problems.
  105. Parameters
  106. ----------
  107. %(params_basic)s
  108. jacobian : Jacobian
  109. A Jacobian approximation: `Jacobian` object or something that
  110. `asjacobian` can transform to one. Alternatively, a string specifying
  111. which of the builtin Jacobian approximations to use:
  112. krylov, broyden1, broyden2, anderson
  113. diagbroyden, linearmixing, excitingmixing
  114. %(params_extra)s
  115. full_output : bool
  116. If true, returns a dictionary `info` containing convergence
  117. information.
  118. raise_exception : bool
  119. If True, a `NoConvergence` exception is raise if no solution is found.
  120. See Also
  121. --------
  122. asjacobian, Jacobian
  123. Notes
  124. -----
  125. This algorithm implements the inexact Newton method, with
  126. backtracking or full line searches. Several Jacobian
  127. approximations are available, including Krylov and Quasi-Newton
  128. methods.
  129. References
  130. ----------
  131. .. [KIM] C. T. Kelley, \"Iterative Methods for Linear and Nonlinear
  132. Equations\". Society for Industrial and Applied Mathematics. (1995)
  133. https://archive.siam.org/books/kelley/fr16/
  134. """
  135. # Can't use default parameters because it's being explicitly passed as None
  136. # from the calling function, so we need to set it here.
  137. tol_norm = maxnorm if tol_norm is None else tol_norm
  138. condition = TerminationCondition(f_tol=f_tol, f_rtol=f_rtol,
  139. x_tol=x_tol, x_rtol=x_rtol,
  140. iter=iter, norm=tol_norm)
  141. x0 = _as_inexact(x0)
  142. def func(z):
  143. return _as_inexact(F(_array_like(z, x0))).flatten()
  144. x = x0.flatten()
  145. dx = np.full_like(x, np.inf)
  146. Fx = func(x)
  147. Fx_norm = norm(Fx)
  148. jacobian = asjacobian(jacobian)
  149. jacobian.setup(x.copy(), Fx, func)
  150. if maxiter is None:
  151. if iter is not None:
  152. maxiter = iter + 1
  153. else:
  154. maxiter = 100*(x.size+1)
  155. if line_search is True:
  156. line_search = 'armijo'
  157. elif line_search is False:
  158. line_search = None
  159. if line_search not in (None, 'armijo', 'wolfe'):
  160. raise ValueError("Invalid line search")
  161. # Solver tolerance selection
  162. gamma = 0.9
  163. eta_max = 0.9999
  164. eta_treshold = 0.1
  165. eta = 1e-3
  166. for n in range(maxiter):
  167. status = condition.check(Fx, x, dx)
  168. if status:
  169. break
  170. # The tolerance, as computed for scipy.sparse.linalg.* routines
  171. tol = min(eta, eta*Fx_norm)
  172. dx = -jacobian.solve(Fx, tol=tol)
  173. if norm(dx) == 0:
  174. raise ValueError("Jacobian inversion yielded zero vector. "
  175. "This indicates a bug in the Jacobian "
  176. "approximation.")
  177. # Line search, or Newton step
  178. if line_search:
  179. s, x, Fx, Fx_norm_new = _nonlin_line_search(func, x, Fx, dx,
  180. line_search)
  181. else:
  182. s = 1.0
  183. x = x + dx
  184. Fx = func(x)
  185. Fx_norm_new = norm(Fx)
  186. jacobian.update(x.copy(), Fx)
  187. if callback:
  188. callback(x, Fx)
  189. # Adjust forcing parameters for inexact methods
  190. eta_A = gamma * Fx_norm_new**2 / Fx_norm**2
  191. if gamma * eta**2 < eta_treshold:
  192. eta = min(eta_max, eta_A)
  193. else:
  194. eta = min(eta_max, max(eta_A, gamma*eta**2))
  195. Fx_norm = Fx_norm_new
  196. # Print status
  197. if verbose:
  198. sys.stdout.write(f"{n}: |F(x)| = {tol_norm(Fx):g}; step {s:g}\n")
  199. sys.stdout.flush()
  200. else:
  201. if raise_exception:
  202. raise NoConvergence(_array_like(x, x0))
  203. else:
  204. status = 2
  205. if full_output:
  206. info = {'nit': condition.iteration,
  207. 'fun': Fx,
  208. 'status': status,
  209. 'success': status == 1,
  210. 'message': {1: 'A solution was found at the specified '
  211. 'tolerance.',
  212. 2: 'The maximum number of iterations allowed '
  213. 'has been reached.'
  214. }[status]
  215. }
  216. return _array_like(x, x0), info
  217. else:
  218. return _array_like(x, x0)
  219. _set_doc(nonlin_solve)
  220. def _nonlin_line_search(func, x, Fx, dx, search_type='armijo', rdiff=1e-8,
  221. smin=1e-2):
  222. tmp_s = [0]
  223. tmp_Fx = [Fx]
  224. tmp_phi = [norm(Fx)**2]
  225. s_norm = norm(x) / norm(dx)
  226. def phi(s, store=True):
  227. if s == tmp_s[0]:
  228. return tmp_phi[0]
  229. xt = x + s*dx
  230. v = func(xt)
  231. p = _safe_norm(v)**2
  232. if store:
  233. tmp_s[0] = s
  234. tmp_phi[0] = p
  235. tmp_Fx[0] = v
  236. return p
  237. def derphi(s):
  238. ds = (abs(s) + s_norm + 1) * rdiff
  239. return (phi(s+ds, store=False) - phi(s)) / ds
  240. if search_type == 'wolfe':
  241. s, phi1, phi0 = scalar_search_wolfe1(phi, derphi, tmp_phi[0],
  242. xtol=1e-2, amin=smin)
  243. elif search_type == 'armijo':
  244. s, phi1 = scalar_search_armijo(phi, tmp_phi[0], -tmp_phi[0],
  245. amin=smin)
  246. if s is None:
  247. # XXX: No suitable step length found. Take the full Newton step,
  248. # and hope for the best.
  249. s = 1.0
  250. x = x + s*dx
  251. if s == tmp_s[0]:
  252. Fx = tmp_Fx[0]
  253. else:
  254. Fx = func(x)
  255. Fx_norm = norm(Fx)
  256. return s, x, Fx, Fx_norm
  257. class TerminationCondition:
  258. """
  259. Termination condition for an iteration. It is terminated if
  260. - |F| < f_rtol*|F_0|, AND
  261. - |F| < f_tol
  262. AND
  263. - |dx| < x_rtol*|x|, AND
  264. - |dx| < x_tol
  265. """
  266. def __init__(self, f_tol=None, f_rtol=None, x_tol=None, x_rtol=None,
  267. iter=None, norm=maxnorm):
  268. if f_tol is None:
  269. f_tol = np.finfo(np.float64).eps ** (1./3)
  270. if f_rtol is None:
  271. f_rtol = np.inf
  272. if x_tol is None:
  273. x_tol = np.inf
  274. if x_rtol is None:
  275. x_rtol = np.inf
  276. self.x_tol = x_tol
  277. self.x_rtol = x_rtol
  278. self.f_tol = f_tol
  279. self.f_rtol = f_rtol
  280. self.norm = norm
  281. self.iter = iter
  282. self.f0_norm = None
  283. self.iteration = 0
  284. def check(self, f, x, dx):
  285. self.iteration += 1
  286. f_norm = self.norm(f)
  287. x_norm = self.norm(x)
  288. dx_norm = self.norm(dx)
  289. if self.f0_norm is None:
  290. self.f0_norm = f_norm
  291. if f_norm == 0:
  292. return 1
  293. if self.iter is not None:
  294. # backwards compatibility with SciPy 0.6.0
  295. return 2 * (self.iteration > self.iter)
  296. # NB: condition must succeed for rtol=inf even if norm == 0
  297. return int((f_norm <= self.f_tol
  298. and f_norm/self.f_rtol <= self.f0_norm)
  299. and (dx_norm <= self.x_tol
  300. and dx_norm/self.x_rtol <= x_norm))
  301. #------------------------------------------------------------------------------
  302. # Generic Jacobian approximation
  303. #------------------------------------------------------------------------------
  304. class Jacobian:
  305. """
  306. Common interface for Jacobians or Jacobian approximations.
  307. The optional methods come useful when implementing trust region
  308. etc., algorithms that often require evaluating transposes of the
  309. Jacobian.
  310. Methods
  311. -------
  312. solve
  313. Returns J^-1 * v
  314. update
  315. Updates Jacobian to point `x` (where the function has residual `Fx`)
  316. matvec : optional
  317. Returns J * v
  318. rmatvec : optional
  319. Returns A^H * v
  320. rsolve : optional
  321. Returns A^-H * v
  322. matmat : optional
  323. Returns A * V, where V is a dense matrix with dimensions (N,K).
  324. todense : optional
  325. Form the dense Jacobian matrix. Necessary for dense trust region
  326. algorithms, and useful for testing.
  327. Attributes
  328. ----------
  329. shape
  330. Matrix dimensions (M, N)
  331. dtype
  332. Data type of the matrix.
  333. func : callable, optional
  334. Function the Jacobian corresponds to
  335. """
  336. # generic type compatibility with scipy-stubs
  337. __class_getitem__ = classmethod(GenericAlias)
  338. def __init__(self, **kw):
  339. names = ["solve", "update", "matvec", "rmatvec", "rsolve",
  340. "matmat", "todense", "shape", "dtype"]
  341. for name, value in kw.items():
  342. if name not in names:
  343. raise ValueError(f"Unknown keyword argument {name}")
  344. if value is not None:
  345. setattr(self, name, kw[name])
  346. if hasattr(self, "todense"):
  347. def __array__(self, dtype=None, copy=None):
  348. if dtype is not None:
  349. raise ValueError(f"`dtype` must be None, was {dtype}")
  350. return self.todense()
  351. def aspreconditioner(self):
  352. return InverseJacobian(self)
  353. def solve(self, v, tol=0):
  354. raise NotImplementedError
  355. def update(self, x, F):
  356. pass
  357. def setup(self, x, F, func):
  358. self.func = func
  359. self.shape = (F.size, x.size)
  360. self.dtype = F.dtype
  361. if self.__class__.setup is Jacobian.setup:
  362. # Call on the first point unless overridden
  363. self.update(x, F)
  364. class InverseJacobian:
  365. """
  366. A simple wrapper that inverts the Jacobian using the `solve` method.
  367. .. legacy:: class
  368. See the newer, more consistent interfaces in :mod:`scipy.optimize`.
  369. Parameters
  370. ----------
  371. jacobian : Jacobian
  372. The Jacobian to invert.
  373. Attributes
  374. ----------
  375. shape
  376. Matrix dimensions (M, N)
  377. dtype
  378. Data type of the matrix.
  379. """
  380. # generic type compatibility with scipy-stubs
  381. __class_getitem__ = classmethod(GenericAlias)
  382. def __init__(self, jacobian):
  383. self.jacobian = jacobian
  384. self.matvec = jacobian.solve
  385. self.update = jacobian.update
  386. if hasattr(jacobian, 'setup'):
  387. self.setup = jacobian.setup
  388. if hasattr(jacobian, 'rsolve'):
  389. self.rmatvec = jacobian.rsolve
  390. @property
  391. def shape(self):
  392. return self.jacobian.shape
  393. @property
  394. def dtype(self):
  395. return self.jacobian.dtype
  396. def asjacobian(J):
  397. """
  398. Convert given object to one suitable for use as a Jacobian.
  399. """
  400. spsolve = scipy.sparse.linalg.spsolve
  401. if isinstance(J, Jacobian):
  402. return J
  403. elif inspect.isclass(J) and issubclass(J, Jacobian):
  404. return J()
  405. elif isinstance(J, np.ndarray):
  406. if J.ndim > 2:
  407. raise ValueError('array must have rank <= 2')
  408. J = np.atleast_2d(np.asarray(J))
  409. if J.shape[0] != J.shape[1]:
  410. raise ValueError('array must be square')
  411. return Jacobian(matvec=lambda v: dot(J, v),
  412. rmatvec=lambda v: dot(J.conj().T, v),
  413. solve=lambda v, tol=0: solve(J, v),
  414. rsolve=lambda v, tol=0: solve(J.conj().T, v),
  415. dtype=J.dtype, shape=J.shape)
  416. elif scipy.sparse.issparse(J):
  417. if J.shape[0] != J.shape[1]:
  418. raise ValueError('matrix must be square')
  419. return Jacobian(matvec=lambda v: J @ v,
  420. rmatvec=lambda v: J.conj().T @ v,
  421. solve=lambda v, tol=0: spsolve(J, v),
  422. rsolve=lambda v, tol=0: spsolve(J.conj().T, v),
  423. dtype=J.dtype, shape=J.shape)
  424. elif hasattr(J, 'shape') and hasattr(J, 'dtype') and hasattr(J, 'solve'):
  425. return Jacobian(matvec=getattr(J, 'matvec'),
  426. rmatvec=getattr(J, 'rmatvec'),
  427. solve=J.solve,
  428. rsolve=getattr(J, 'rsolve'),
  429. update=getattr(J, 'update'),
  430. setup=getattr(J, 'setup'),
  431. dtype=J.dtype,
  432. shape=J.shape)
  433. elif callable(J):
  434. # Assume it's a function J(x) that returns the Jacobian
  435. class Jac(Jacobian):
  436. def update(self, x, F):
  437. self.x = x
  438. def solve(self, v, tol=0):
  439. m = J(self.x)
  440. if isinstance(m, np.ndarray):
  441. return solve(m, v)
  442. elif scipy.sparse.issparse(m):
  443. return spsolve(m, v)
  444. else:
  445. raise ValueError("Unknown matrix type")
  446. def matvec(self, v):
  447. m = J(self.x)
  448. if isinstance(m, np.ndarray):
  449. return dot(m, v)
  450. elif scipy.sparse.issparse(m):
  451. return m @ v
  452. else:
  453. raise ValueError("Unknown matrix type")
  454. def rsolve(self, v, tol=0):
  455. m = J(self.x)
  456. if isinstance(m, np.ndarray):
  457. return solve(m.conj().T, v)
  458. elif scipy.sparse.issparse(m):
  459. return spsolve(m.conj().T, v)
  460. else:
  461. raise ValueError("Unknown matrix type")
  462. def rmatvec(self, v):
  463. m = J(self.x)
  464. if isinstance(m, np.ndarray):
  465. return dot(m.conj().T, v)
  466. elif scipy.sparse.issparse(m):
  467. return m.conj().T @ v
  468. else:
  469. raise ValueError("Unknown matrix type")
  470. return Jac()
  471. elif isinstance(J, str):
  472. return dict(broyden1=BroydenFirst,
  473. broyden2=BroydenSecond,
  474. anderson=Anderson,
  475. diagbroyden=DiagBroyden,
  476. linearmixing=LinearMixing,
  477. excitingmixing=ExcitingMixing,
  478. krylov=KrylovJacobian)[J]()
  479. else:
  480. raise TypeError('Cannot convert object to a Jacobian')
  481. #------------------------------------------------------------------------------
  482. # Broyden
  483. #------------------------------------------------------------------------------
  484. class GenericBroyden(Jacobian):
  485. # generic type compatibility with scipy-stubs
  486. __class_getitem__ = classmethod(GenericAlias)
  487. def setup(self, x0, f0, func):
  488. Jacobian.setup(self, x0, f0, func)
  489. self.last_f = f0
  490. self.last_x = x0
  491. if hasattr(self, 'alpha') and self.alpha is None:
  492. # Autoscale the initial Jacobian parameter
  493. # unless we have already guessed the solution.
  494. normf0 = norm(f0)
  495. if normf0:
  496. self.alpha = 0.5*max(norm(x0), 1) / normf0
  497. else:
  498. self.alpha = 1.0
  499. def _update(self, x, f, dx, df, dx_norm, df_norm):
  500. raise NotImplementedError
  501. def update(self, x, f):
  502. df = f - self.last_f
  503. dx = x - self.last_x
  504. self._update(x, f, dx, df, norm(dx), norm(df))
  505. self.last_f = f
  506. self.last_x = x
  507. class LowRankMatrix:
  508. r"""
  509. A matrix represented as
  510. .. math:: \alpha I + \sum_{n=0}^{n=M} c_n d_n^\dagger
  511. However, if the rank of the matrix reaches the dimension of the vectors,
  512. full matrix representation will be used thereon.
  513. """
  514. # generic type compatibility with scipy-stubs
  515. __class_getitem__ = classmethod(GenericAlias)
  516. def __init__(self, alpha, n, dtype):
  517. self.alpha = alpha
  518. self.cs = []
  519. self.ds = []
  520. self.n = n
  521. self.dtype = dtype
  522. self.collapsed = None
  523. @staticmethod
  524. def _matvec(v, alpha, cs, ds):
  525. axpy, scal, dotc = get_blas_funcs(['axpy', 'scal', 'dotc'],
  526. cs[:1] + [v])
  527. w = alpha * v
  528. for c, d in zip(cs, ds):
  529. a = dotc(d, v)
  530. w = axpy(c, w, w.size, a)
  531. return w
  532. @staticmethod
  533. def _solve(v, alpha, cs, ds):
  534. """Evaluate w = M^-1 v"""
  535. if len(cs) == 0:
  536. return v/alpha
  537. # (B + C D^H)^-1 = B^-1 - B^-1 C (I + D^H B^-1 C)^-1 D^H B^-1
  538. axpy, dotc = get_blas_funcs(['axpy', 'dotc'], cs[:1] + [v])
  539. c0 = cs[0]
  540. A = alpha * np.identity(len(cs), dtype=c0.dtype)
  541. for i, d in enumerate(ds):
  542. for j, c in enumerate(cs):
  543. A[i,j] += dotc(d, c)
  544. q = np.zeros(len(cs), dtype=c0.dtype)
  545. for j, d in enumerate(ds):
  546. q[j] = dotc(d, v)
  547. q /= alpha
  548. q = solve(A, q)
  549. w = v/alpha
  550. for c, qc in zip(cs, q):
  551. w = axpy(c, w, w.size, -qc)
  552. return w
  553. def matvec(self, v):
  554. """Evaluate w = M v"""
  555. if self.collapsed is not None:
  556. return np.dot(self.collapsed, v)
  557. return LowRankMatrix._matvec(v, self.alpha, self.cs, self.ds)
  558. def rmatvec(self, v):
  559. """Evaluate w = M^H v"""
  560. if self.collapsed is not None:
  561. return np.dot(self.collapsed.T.conj(), v)
  562. return LowRankMatrix._matvec(v, np.conj(self.alpha), self.ds, self.cs)
  563. def solve(self, v, tol=0):
  564. """Evaluate w = M^-1 v"""
  565. if self.collapsed is not None:
  566. return solve(self.collapsed, v)
  567. return LowRankMatrix._solve(v, self.alpha, self.cs, self.ds)
  568. def rsolve(self, v, tol=0):
  569. """Evaluate w = M^-H v"""
  570. if self.collapsed is not None:
  571. return solve(self.collapsed.T.conj(), v)
  572. return LowRankMatrix._solve(v, np.conj(self.alpha), self.ds, self.cs)
  573. def append(self, c, d):
  574. if self.collapsed is not None:
  575. self.collapsed += c[:,None] * d[None,:].conj()
  576. return
  577. self.cs.append(c)
  578. self.ds.append(d)
  579. if len(self.cs) > c.size:
  580. self.collapse()
  581. def __array__(self, dtype=None, copy=None):
  582. if dtype is not None:
  583. warnings.warn("LowRankMatrix is scipy-internal code, `dtype` "
  584. f"should only be None but was {dtype} (not handled)",
  585. stacklevel=3)
  586. if copy is not None:
  587. warnings.warn("LowRankMatrix is scipy-internal code, `copy` "
  588. f"should only be None but was {copy} (not handled)",
  589. stacklevel=3)
  590. if self.collapsed is not None:
  591. return self.collapsed
  592. Gm = self.alpha*np.identity(self.n, dtype=self.dtype)
  593. for c, d in zip(self.cs, self.ds):
  594. Gm += c[:,None]*d[None,:].conj()
  595. return Gm
  596. def collapse(self):
  597. """Collapse the low-rank matrix to a full-rank one."""
  598. self.collapsed = np.array(self, copy=copy_if_needed)
  599. self.cs = None
  600. self.ds = None
  601. self.alpha = None
  602. def restart_reduce(self, rank):
  603. """
  604. Reduce the rank of the matrix by dropping all vectors.
  605. """
  606. if self.collapsed is not None:
  607. return
  608. assert rank > 0
  609. if len(self.cs) > rank:
  610. del self.cs[:]
  611. del self.ds[:]
  612. def simple_reduce(self, rank):
  613. """
  614. Reduce the rank of the matrix by dropping oldest vectors.
  615. """
  616. if self.collapsed is not None:
  617. return
  618. assert rank > 0
  619. while len(self.cs) > rank:
  620. del self.cs[0]
  621. del self.ds[0]
  622. def svd_reduce(self, max_rank, to_retain=None):
  623. """
  624. Reduce the rank of the matrix by retaining some SVD components.
  625. This corresponds to the \"Broyden Rank Reduction Inverse\"
  626. algorithm described in [1]_.
  627. Note that the SVD decomposition can be done by solving only a
  628. problem whose size is the effective rank of this matrix, which
  629. is viable even for large problems.
  630. Parameters
  631. ----------
  632. max_rank : int
  633. Maximum rank of this matrix after reduction.
  634. to_retain : int, optional
  635. Number of SVD components to retain when reduction is done
  636. (ie. rank > max_rank). Default is ``max_rank - 2``.
  637. References
  638. ----------
  639. .. [1] B.A. van der Rotten, PhD thesis,
  640. \"A limited memory Broyden method to solve high-dimensional
  641. systems of nonlinear equations\". Mathematisch Instituut,
  642. Universiteit Leiden, The Netherlands (2003).
  643. https://web.archive.org/web/20161022015821/http://www.math.leidenuniv.nl/scripties/Rotten.pdf
  644. """
  645. if self.collapsed is not None:
  646. return
  647. p = max_rank
  648. if to_retain is not None:
  649. q = to_retain
  650. else:
  651. q = p - 2
  652. if self.cs:
  653. p = min(p, len(self.cs[0]))
  654. q = max(0, min(q, p-1))
  655. m = len(self.cs)
  656. if m < p:
  657. # nothing to do
  658. return
  659. C = np.array(self.cs).T
  660. D = np.array(self.ds).T
  661. D, R = qr(D, mode='economic')
  662. C = dot(C, R.T.conj())
  663. U, S, WH = svd(C, full_matrices=False)
  664. C = dot(C, inv(WH))
  665. D = dot(D, WH.T.conj())
  666. for k in range(q):
  667. self.cs[k] = C[:,k].copy()
  668. self.ds[k] = D[:,k].copy()
  669. del self.cs[q:]
  670. del self.ds[q:]
  671. _doc_parts['broyden_params'] = _dedent_for_py313("""
  672. alpha : float, optional
  673. Initial guess for the Jacobian is ``(-1/alpha)``.
  674. reduction_method : str or tuple, optional
  675. Method used in ensuring that the rank of the Broyden matrix
  676. stays low. Can either be a string giving the name of the method,
  677. or a tuple of the form ``(method, param1, param2, ...)``
  678. that gives the name of the method and values for additional parameters.
  679. Methods available:
  680. - ``restart``: drop all matrix columns. Has no extra parameters.
  681. - ``simple``: drop oldest matrix column. Has no extra parameters.
  682. - ``svd``: keep only the most significant SVD components.
  683. Takes an extra parameter, ``to_retain``, which determines the
  684. number of SVD components to retain when rank reduction is done.
  685. Default is ``max_rank - 2``.
  686. max_rank : int, optional
  687. Maximum rank for the Broyden matrix.
  688. Default is infinity (i.e., no rank reduction).
  689. """).strip()
  690. class BroydenFirst(GenericBroyden):
  691. """
  692. Find a root of a function, using Broyden's first Jacobian approximation.
  693. This method is also known as "Broyden's good method".
  694. Parameters
  695. ----------
  696. %(params_basic)s
  697. %(broyden_params)s
  698. %(params_extra)s
  699. See Also
  700. --------
  701. root : Interface to root finding algorithms for multivariate
  702. functions. See ``method='broyden1'`` in particular.
  703. Notes
  704. -----
  705. This algorithm implements the inverse Jacobian Quasi-Newton update
  706. .. math:: H_+ = H + (dx - H df) dx^\\dagger H / ( dx^\\dagger H df)
  707. which corresponds to Broyden's first Jacobian update
  708. .. math:: J_+ = J + (df - J dx) dx^\\dagger / dx^\\dagger dx
  709. References
  710. ----------
  711. .. [1] B.A. van der Rotten, PhD thesis,
  712. "A limited memory Broyden method to solve high-dimensional
  713. systems of nonlinear equations". Mathematisch Instituut,
  714. Universiteit Leiden, The Netherlands (2003).
  715. https://math.leidenuniv.nl/scripties/Rotten.pdf
  716. Examples
  717. --------
  718. The following functions define a system of nonlinear equations
  719. >>> def fun(x):
  720. ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
  721. ... 0.5 * (x[1] - x[0])**3 + x[1]]
  722. A solution can be obtained as follows.
  723. >>> from scipy import optimize
  724. >>> sol = optimize.broyden1(fun, [0, 0])
  725. >>> sol
  726. array([0.84116396, 0.15883641])
  727. """
  728. def __init__(self, alpha=None, reduction_method='restart', max_rank=None):
  729. GenericBroyden.__init__(self)
  730. self.alpha = alpha
  731. self.Gm = None
  732. if max_rank is None:
  733. max_rank = np.inf
  734. self.max_rank = max_rank
  735. if isinstance(reduction_method, str):
  736. reduce_params = ()
  737. else:
  738. reduce_params = reduction_method[1:]
  739. reduction_method = reduction_method[0]
  740. reduce_params = (max_rank - 1,) + reduce_params
  741. if reduction_method == 'svd':
  742. self._reduce = lambda: self.Gm.svd_reduce(*reduce_params)
  743. elif reduction_method == 'simple':
  744. self._reduce = lambda: self.Gm.simple_reduce(*reduce_params)
  745. elif reduction_method == 'restart':
  746. self._reduce = lambda: self.Gm.restart_reduce(*reduce_params)
  747. else:
  748. raise ValueError(f"Unknown rank reduction method '{reduction_method}'")
  749. def setup(self, x, F, func):
  750. GenericBroyden.setup(self, x, F, func)
  751. self.Gm = LowRankMatrix(-self.alpha, self.shape[0], self.dtype)
  752. def todense(self):
  753. return inv(self.Gm)
  754. def solve(self, f, tol=0):
  755. r = self.Gm.matvec(f)
  756. if not np.isfinite(r).all():
  757. # singular; reset the Jacobian approximation
  758. self.setup(self.last_x, self.last_f, self.func)
  759. return self.Gm.matvec(f)
  760. return r
  761. def matvec(self, f):
  762. return self.Gm.solve(f)
  763. def rsolve(self, f, tol=0):
  764. return self.Gm.rmatvec(f)
  765. def rmatvec(self, f):
  766. return self.Gm.rsolve(f)
  767. def _update(self, x, f, dx, df, dx_norm, df_norm):
  768. self._reduce() # reduce first to preserve secant condition
  769. v = self.Gm.rmatvec(dx)
  770. c = dx - self.Gm.matvec(df)
  771. d = v / vdot(df, v)
  772. self.Gm.append(c, d)
  773. class BroydenSecond(BroydenFirst):
  774. """
  775. Find a root of a function, using Broyden\'s second Jacobian approximation.
  776. This method is also known as \"Broyden's bad method\".
  777. Parameters
  778. ----------
  779. %(params_basic)s
  780. %(broyden_params)s
  781. %(params_extra)s
  782. See Also
  783. --------
  784. root : Interface to root finding algorithms for multivariate
  785. functions. See ``method='broyden2'`` in particular.
  786. Notes
  787. -----
  788. This algorithm implements the inverse Jacobian Quasi-Newton update
  789. .. math:: H_+ = H + (dx - H df) df^\\dagger / ( df^\\dagger df)
  790. corresponding to Broyden's second method.
  791. References
  792. ----------
  793. .. [1] B.A. van der Rotten, PhD thesis,
  794. \"A limited memory Broyden method to solve high-dimensional
  795. systems of nonlinear equations\". Mathematisch Instituut,
  796. Universiteit Leiden, The Netherlands (2003).
  797. https://web.archive.org/web/20161022015821/http://www.math.leidenuniv.nl/scripties/Rotten.pdf
  798. Examples
  799. --------
  800. The following functions define a system of nonlinear equations
  801. >>> def fun(x):
  802. ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
  803. ... 0.5 * (x[1] - x[0])**3 + x[1]]
  804. A solution can be obtained as follows.
  805. >>> from scipy import optimize
  806. >>> sol = optimize.broyden2(fun, [0, 0])
  807. >>> sol
  808. array([0.84116365, 0.15883529])
  809. """
  810. def _update(self, x, f, dx, df, dx_norm, df_norm):
  811. self._reduce() # reduce first to preserve secant condition
  812. v = df
  813. c = dx - self.Gm.matvec(df)
  814. d = v / df_norm**2
  815. self.Gm.append(c, d)
  816. #------------------------------------------------------------------------------
  817. # Broyden-like (restricted memory)
  818. #------------------------------------------------------------------------------
  819. class Anderson(GenericBroyden):
  820. """
  821. Find a root of a function, using (extended) Anderson mixing.
  822. The Jacobian is formed by for a 'best' solution in the space
  823. spanned by last `M` vectors. As a result, only a MxM matrix
  824. inversions and MxN multiplications are required. [Ey]_
  825. Parameters
  826. ----------
  827. %(params_basic)s
  828. alpha : float, optional
  829. Initial guess for the Jacobian is (-1/alpha).
  830. M : float, optional
  831. Number of previous vectors to retain. Defaults to 5.
  832. w0 : float, optional
  833. Regularization parameter for numerical stability.
  834. Compared to unity, good values of the order of 0.01.
  835. %(params_extra)s
  836. See Also
  837. --------
  838. root : Interface to root finding algorithms for multivariate
  839. functions. See ``method='anderson'`` in particular.
  840. References
  841. ----------
  842. .. [Ey] V. Eyert, J. Comp. Phys., 124, 271 (1996).
  843. Examples
  844. --------
  845. The following functions define a system of nonlinear equations
  846. >>> def fun(x):
  847. ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
  848. ... 0.5 * (x[1] - x[0])**3 + x[1]]
  849. A solution can be obtained as follows.
  850. >>> from scipy import optimize
  851. >>> sol = optimize.anderson(fun, [0, 0])
  852. >>> sol
  853. array([0.84116588, 0.15883789])
  854. """
  855. # Note:
  856. #
  857. # Anderson method maintains a rank M approximation of the inverse Jacobian,
  858. #
  859. # J^-1 v ~ -v*alpha + (dX + alpha dF) A^-1 dF^H v
  860. # A = W + dF^H dF
  861. # W = w0^2 diag(dF^H dF)
  862. #
  863. # so that for w0 = 0 the secant condition applies for last M iterates, i.e.,
  864. #
  865. # J^-1 df_j = dx_j
  866. #
  867. # for all j = 0 ... M-1.
  868. #
  869. # Moreover, (from Sherman-Morrison-Woodbury formula)
  870. #
  871. # J v ~ [ b I - b^2 C (I + b dF^H A^-1 C)^-1 dF^H ] v
  872. # C = (dX + alpha dF) A^-1
  873. # b = -1/alpha
  874. #
  875. # and after simplification
  876. #
  877. # J v ~ -v/alpha + (dX/alpha + dF) (dF^H dX - alpha W)^-1 dF^H v
  878. #
  879. def __init__(self, alpha=None, w0=0.01, M=5):
  880. GenericBroyden.__init__(self)
  881. self.alpha = alpha
  882. self.M = M
  883. self.dx = []
  884. self.df = []
  885. self.gamma = None
  886. self.w0 = w0
  887. def solve(self, f, tol=0):
  888. dx = -self.alpha*f
  889. n = len(self.dx)
  890. if n == 0:
  891. return dx
  892. df_f = np.empty(n, dtype=f.dtype)
  893. for k in range(n):
  894. df_f[k] = vdot(self.df[k], f)
  895. try:
  896. gamma = solve(self.a, df_f)
  897. except LinAlgError:
  898. # singular; reset the Jacobian approximation
  899. del self.dx[:]
  900. del self.df[:]
  901. return dx
  902. for m in range(n):
  903. dx += gamma[m]*(self.dx[m] + self.alpha*self.df[m])
  904. return dx
  905. def matvec(self, f):
  906. dx = -f/self.alpha
  907. n = len(self.dx)
  908. if n == 0:
  909. return dx
  910. df_f = np.empty(n, dtype=f.dtype)
  911. for k in range(n):
  912. df_f[k] = vdot(self.df[k], f)
  913. b = np.empty((n, n), dtype=f.dtype)
  914. for i in range(n):
  915. for j in range(n):
  916. b[i,j] = vdot(self.df[i], self.dx[j])
  917. if i == j and self.w0 != 0:
  918. b[i,j] -= vdot(self.df[i], self.df[i])*self.w0**2*self.alpha
  919. gamma = solve(b, df_f)
  920. for m in range(n):
  921. dx += gamma[m]*(self.df[m] + self.dx[m]/self.alpha)
  922. return dx
  923. def _update(self, x, f, dx, df, dx_norm, df_norm):
  924. if self.M == 0:
  925. return
  926. self.dx.append(dx)
  927. self.df.append(df)
  928. while len(self.dx) > self.M:
  929. self.dx.pop(0)
  930. self.df.pop(0)
  931. n = len(self.dx)
  932. a = np.zeros((n, n), dtype=f.dtype)
  933. for i in range(n):
  934. for j in range(i, n):
  935. if i == j:
  936. wd = self.w0**2
  937. else:
  938. wd = 0
  939. a[i,j] = (1+wd)*vdot(self.df[i], self.df[j])
  940. a += np.triu(a, 1).T.conj()
  941. self.a = a
  942. #------------------------------------------------------------------------------
  943. # Simple iterations
  944. #------------------------------------------------------------------------------
  945. class DiagBroyden(GenericBroyden):
  946. """
  947. Find a root of a function, using diagonal Broyden Jacobian approximation.
  948. The Jacobian approximation is derived from previous iterations, by
  949. retaining only the diagonal of Broyden matrices.
  950. .. warning::
  951. This algorithm may be useful for specific problems, but whether
  952. it will work may depend strongly on the problem.
  953. Parameters
  954. ----------
  955. %(params_basic)s
  956. alpha : float, optional
  957. Initial guess for the Jacobian is (-1/alpha).
  958. %(params_extra)s
  959. See Also
  960. --------
  961. root : Interface to root finding algorithms for multivariate
  962. functions. See ``method='diagbroyden'`` in particular.
  963. Examples
  964. --------
  965. The following functions define a system of nonlinear equations
  966. >>> def fun(x):
  967. ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
  968. ... 0.5 * (x[1] - x[0])**3 + x[1]]
  969. A solution can be obtained as follows.
  970. >>> from scipy import optimize
  971. >>> sol = optimize.diagbroyden(fun, [0, 0])
  972. >>> sol
  973. array([0.84116403, 0.15883384])
  974. """
  975. def __init__(self, alpha=None):
  976. GenericBroyden.__init__(self)
  977. self.alpha = alpha
  978. def setup(self, x, F, func):
  979. GenericBroyden.setup(self, x, F, func)
  980. self.d = np.full((self.shape[0],), 1 / self.alpha, dtype=self.dtype)
  981. def solve(self, f, tol=0):
  982. return -f / self.d
  983. def matvec(self, f):
  984. return -f * self.d
  985. def rsolve(self, f, tol=0):
  986. return -f / self.d.conj()
  987. def rmatvec(self, f):
  988. return -f * self.d.conj()
  989. def todense(self):
  990. return np.diag(-self.d)
  991. def _update(self, x, f, dx, df, dx_norm, df_norm):
  992. self.d -= (df + self.d*dx)*dx/dx_norm**2
  993. class LinearMixing(GenericBroyden):
  994. """
  995. Find a root of a function, using a scalar Jacobian approximation.
  996. .. warning::
  997. This algorithm may be useful for specific problems, but whether
  998. it will work may depend strongly on the problem.
  999. Parameters
  1000. ----------
  1001. %(params_basic)s
  1002. alpha : float, optional
  1003. The Jacobian approximation is (-1/alpha).
  1004. %(params_extra)s
  1005. See Also
  1006. --------
  1007. root : Interface to root finding algorithms for multivariate
  1008. functions. See ``method='linearmixing'`` in particular.
  1009. """
  1010. def __init__(self, alpha=None):
  1011. GenericBroyden.__init__(self)
  1012. self.alpha = alpha
  1013. def solve(self, f, tol=0):
  1014. return -f*self.alpha
  1015. def matvec(self, f):
  1016. return -f/self.alpha
  1017. def rsolve(self, f, tol=0):
  1018. return -f*np.conj(self.alpha)
  1019. def rmatvec(self, f):
  1020. return -f/np.conj(self.alpha)
  1021. def todense(self):
  1022. return np.diag(np.full(self.shape[0], -1/self.alpha))
  1023. def _update(self, x, f, dx, df, dx_norm, df_norm):
  1024. pass
  1025. class ExcitingMixing(GenericBroyden):
  1026. """
  1027. Find a root of a function, using a tuned diagonal Jacobian approximation.
  1028. The Jacobian matrix is diagonal and is tuned on each iteration.
  1029. .. warning::
  1030. This algorithm may be useful for specific problems, but whether
  1031. it will work may depend strongly on the problem.
  1032. See Also
  1033. --------
  1034. root : Interface to root finding algorithms for multivariate
  1035. functions. See ``method='excitingmixing'`` in particular.
  1036. Parameters
  1037. ----------
  1038. %(params_basic)s
  1039. alpha : float, optional
  1040. Initial Jacobian approximation is (-1/alpha).
  1041. alphamax : float, optional
  1042. The entries of the diagonal Jacobian are kept in the range
  1043. ``[alpha, alphamax]``.
  1044. %(params_extra)s
  1045. """
  1046. def __init__(self, alpha=None, alphamax=1.0):
  1047. GenericBroyden.__init__(self)
  1048. self.alpha = alpha
  1049. self.alphamax = alphamax
  1050. self.beta = None
  1051. def setup(self, x, F, func):
  1052. GenericBroyden.setup(self, x, F, func)
  1053. self.beta = np.full((self.shape[0],), self.alpha, dtype=self.dtype)
  1054. def solve(self, f, tol=0):
  1055. return -f*self.beta
  1056. def matvec(self, f):
  1057. return -f/self.beta
  1058. def rsolve(self, f, tol=0):
  1059. return -f*self.beta.conj()
  1060. def rmatvec(self, f):
  1061. return -f/self.beta.conj()
  1062. def todense(self):
  1063. return np.diag(-1/self.beta)
  1064. def _update(self, x, f, dx, df, dx_norm, df_norm):
  1065. incr = f*self.last_f > 0
  1066. self.beta[incr] += self.alpha
  1067. self.beta[~incr] = self.alpha
  1068. np.clip(self.beta, 0, self.alphamax, out=self.beta)
  1069. #------------------------------------------------------------------------------
  1070. # Iterative/Krylov approximated Jacobians
  1071. #------------------------------------------------------------------------------
  1072. class KrylovJacobian(Jacobian):
  1073. """
  1074. Find a root of a function, using Krylov approximation for inverse Jacobian.
  1075. This method is suitable for solving large-scale problems.
  1076. Parameters
  1077. ----------
  1078. %(params_basic)s
  1079. rdiff : float, optional
  1080. Relative step size to use in numerical differentiation.
  1081. method : str or callable, optional
  1082. Krylov method to use to approximate the Jacobian. Can be a string,
  1083. or a function implementing the same interface as the iterative
  1084. solvers in `scipy.sparse.linalg`. If a string, needs to be one of:
  1085. ``'lgmres'``, ``'gmres'``, ``'bicgstab'``, ``'cgs'``, ``'minres'``,
  1086. ``'tfqmr'``.
  1087. The default is `scipy.sparse.linalg.lgmres`.
  1088. inner_maxiter : int, optional
  1089. Parameter to pass to the "inner" Krylov solver: maximum number of
  1090. iterations. Iteration will stop after maxiter steps even if the
  1091. specified tolerance has not been achieved.
  1092. inner_M : LinearOperator or InverseJacobian
  1093. Preconditioner for the inner Krylov iteration.
  1094. Note that you can use also inverse Jacobians as (adaptive)
  1095. preconditioners. For example,
  1096. >>> from scipy.optimize import BroydenFirst, KrylovJacobian
  1097. >>> from scipy.optimize import InverseJacobian
  1098. >>> jac = BroydenFirst()
  1099. >>> kjac = KrylovJacobian(inner_M=InverseJacobian(jac))
  1100. If the preconditioner has a method named 'update', it will be called
  1101. as ``update(x, f)`` after each nonlinear step, with ``x`` giving
  1102. the current point, and ``f`` the current function value.
  1103. outer_k : int, optional
  1104. Size of the subspace kept across LGMRES nonlinear iterations.
  1105. See `scipy.sparse.linalg.lgmres` for details.
  1106. inner_kwargs : kwargs
  1107. Keyword parameters for the "inner" Krylov solver
  1108. (defined with `method`). Parameter names must start with
  1109. the `inner_` prefix which will be stripped before passing on
  1110. the inner method. See, e.g., `scipy.sparse.linalg.gmres` for details.
  1111. %(params_extra)s
  1112. See Also
  1113. --------
  1114. root : Interface to root finding algorithms for multivariate
  1115. functions. See ``method='krylov'`` in particular.
  1116. scipy.sparse.linalg.gmres
  1117. scipy.sparse.linalg.lgmres
  1118. Notes
  1119. -----
  1120. This function implements a Newton-Krylov solver. The basic idea is
  1121. to compute the inverse of the Jacobian with an iterative Krylov
  1122. method. These methods require only evaluating the Jacobian-vector
  1123. products, which are conveniently approximated by a finite difference:
  1124. .. math:: J v \\approx (f(x + \\omega*v/|v|) - f(x)) / \\omega
  1125. Due to the use of iterative matrix inverses, these methods can
  1126. deal with large nonlinear problems.
  1127. SciPy's `scipy.sparse.linalg` module offers a selection of Krylov
  1128. solvers to choose from. The default here is `lgmres`, which is a
  1129. variant of restarted GMRES iteration that reuses some of the
  1130. information obtained in the previous Newton steps to invert
  1131. Jacobians in subsequent steps.
  1132. For a review on Newton-Krylov methods, see for example [1]_,
  1133. and for the LGMRES sparse inverse method, see [2]_.
  1134. References
  1135. ----------
  1136. .. [1] C. T. Kelley, Solving Nonlinear Equations with Newton's Method,
  1137. SIAM, pp.57-83, 2003.
  1138. :doi:`10.1137/1.9780898718898.ch3`
  1139. .. [2] D.A. Knoll and D.E. Keyes, J. Comp. Phys. 193, 357 (2004).
  1140. :doi:`10.1016/j.jcp.2003.08.010`
  1141. .. [3] A.H. Baker and E.R. Jessup and T. Manteuffel,
  1142. SIAM J. Matrix Anal. Appl. 26, 962 (2005).
  1143. :doi:`10.1137/S0895479803422014`
  1144. Examples
  1145. --------
  1146. The following functions define a system of nonlinear equations
  1147. >>> def fun(x):
  1148. ... return [x[0] + 0.5 * x[1] - 1.0,
  1149. ... 0.5 * (x[1] - x[0]) ** 2]
  1150. A solution can be obtained as follows.
  1151. >>> from scipy import optimize
  1152. >>> sol = optimize.newton_krylov(fun, [0, 0])
  1153. >>> sol
  1154. array([0.66731771, 0.66536458])
  1155. """
  1156. def __init__(self, rdiff=None, method='lgmres', inner_maxiter=20,
  1157. inner_M=None, outer_k=10, **kw):
  1158. self.preconditioner = inner_M
  1159. self.rdiff = rdiff
  1160. # Note that this retrieves one of the named functions, or otherwise
  1161. # uses `method` as is (i.e., for a user-provided callable).
  1162. self.method = dict(
  1163. bicgstab=scipy.sparse.linalg.bicgstab,
  1164. gmres=scipy.sparse.linalg.gmres,
  1165. lgmres=scipy.sparse.linalg.lgmres,
  1166. cgs=scipy.sparse.linalg.cgs,
  1167. minres=scipy.sparse.linalg.minres,
  1168. tfqmr=scipy.sparse.linalg.tfqmr,
  1169. ).get(method, method)
  1170. self.method_kw = dict(maxiter=inner_maxiter, M=self.preconditioner)
  1171. if self.method is scipy.sparse.linalg.gmres:
  1172. # Replace GMRES's outer iteration with Newton steps
  1173. self.method_kw['restart'] = inner_maxiter
  1174. self.method_kw['maxiter'] = 1
  1175. self.method_kw.setdefault('atol', 0)
  1176. elif self.method in (scipy.sparse.linalg.gcrotmk,
  1177. scipy.sparse.linalg.bicgstab,
  1178. scipy.sparse.linalg.cgs):
  1179. self.method_kw.setdefault('atol', 0)
  1180. elif self.method is scipy.sparse.linalg.lgmres:
  1181. self.method_kw['outer_k'] = outer_k
  1182. # Replace LGMRES's outer iteration with Newton steps
  1183. self.method_kw['maxiter'] = 1
  1184. # Carry LGMRES's `outer_v` vectors across nonlinear iterations
  1185. self.method_kw.setdefault('outer_v', [])
  1186. self.method_kw.setdefault('prepend_outer_v', True)
  1187. # But don't carry the corresponding Jacobian*v products, in case
  1188. # the Jacobian changes a lot in the nonlinear step
  1189. #
  1190. # XXX: some trust-region inspired ideas might be more efficient...
  1191. # See e.g., Brown & Saad. But needs to be implemented separately
  1192. # since it's not an inexact Newton method.
  1193. self.method_kw.setdefault('store_outer_Av', False)
  1194. self.method_kw.setdefault('atol', 0)
  1195. # Retrieve the signature of the method to find the valid parameters
  1196. valid_inner_params = [
  1197. k for k in signature(self.method).parameters
  1198. if k not in ('self', 'args', 'kwargs')
  1199. ]
  1200. for key, value in kw.items():
  1201. if not key.startswith("inner_"):
  1202. raise ValueError(f"Unknown parameter {key}")
  1203. if key[6:] not in valid_inner_params:
  1204. # Use difflib to find close matches to the invalid key
  1205. inner_param_suggestions = get_close_matches(key[6:],
  1206. valid_inner_params,
  1207. n=1)
  1208. if inner_param_suggestions:
  1209. suggestion_msg = (f" Did you mean '"
  1210. f"{inner_param_suggestions[0]}'?")
  1211. else:
  1212. suggestion_msg = ""
  1213. # warn user that the parameter is not valid for the inner method
  1214. warnings.warn(
  1215. f"Option '{key}' is invalid for the inner method: {method}."
  1216. " It will be ignored."
  1217. "Please check inner method documentation for valid options."
  1218. + suggestion_msg,
  1219. stacklevel=3,
  1220. category=UserWarning,
  1221. # using `skip_file_prefixes` would be a good idea
  1222. # and should be added once we drop support for Python 3.11
  1223. )
  1224. # ignore this parameter and continue
  1225. continue
  1226. self.method_kw[key[6:]] = value
  1227. def _update_diff_step(self):
  1228. mx = abs(self.x0).max()
  1229. mf = abs(self.f0).max()
  1230. self.omega = self.rdiff * max(1, mx) / max(1, mf)
  1231. def matvec(self, v):
  1232. nv = norm(v)
  1233. if nv == 0:
  1234. return 0*v
  1235. sc = self.omega / nv
  1236. r = (self.func(self.x0 + sc*v) - self.f0) / sc
  1237. if not np.all(np.isfinite(r)) and np.all(np.isfinite(v)):
  1238. raise ValueError('Function returned non-finite results')
  1239. return r
  1240. def solve(self, rhs, tol=0):
  1241. if 'rtol' in self.method_kw:
  1242. sol, info = self.method(self.op, rhs, **self.method_kw)
  1243. else:
  1244. sol, info = self.method(self.op, rhs, rtol=tol, **self.method_kw)
  1245. return sol
  1246. def update(self, x, f):
  1247. self.x0 = x
  1248. self.f0 = f
  1249. self._update_diff_step()
  1250. # Update also the preconditioner, if possible
  1251. if self.preconditioner is not None:
  1252. if hasattr(self.preconditioner, 'update'):
  1253. self.preconditioner.update(x, f)
  1254. def setup(self, x, f, func):
  1255. Jacobian.setup(self, x, f, func)
  1256. self.x0 = x
  1257. self.f0 = f
  1258. self.op = scipy.sparse.linalg.aslinearoperator(self)
  1259. if self.rdiff is None:
  1260. self.rdiff = np.finfo(x.dtype).eps ** (1./2)
  1261. self._update_diff_step()
  1262. # Setup also the preconditioner, if possible
  1263. if self.preconditioner is not None:
  1264. if hasattr(self.preconditioner, 'setup'):
  1265. self.preconditioner.setup(x, f, func)
  1266. #------------------------------------------------------------------------------
  1267. # Wrapper functions
  1268. #------------------------------------------------------------------------------
  1269. def _nonlin_wrapper(name, jac):
  1270. """
  1271. Construct a solver wrapper with given name and Jacobian approx.
  1272. It inspects the keyword arguments of ``jac.__init__``, and allows to
  1273. use the same arguments in the wrapper function, in addition to the
  1274. keyword arguments of `nonlin_solve`
  1275. """
  1276. signature = _getfullargspec(jac.__init__)
  1277. args, varargs, varkw, defaults, kwonlyargs, kwdefaults, _ = signature
  1278. kwargs = list(zip(args[-len(defaults):], defaults))
  1279. kw_str = ", ".join([f"{k}={v!r}" for k, v in kwargs])
  1280. if kw_str:
  1281. kw_str = ", " + kw_str
  1282. kwkw_str = ", ".join([f"{k}={k}" for k, v in kwargs])
  1283. if kwkw_str:
  1284. kwkw_str = kwkw_str + ", "
  1285. if kwonlyargs:
  1286. raise ValueError(f'Unexpected signature {signature}')
  1287. # Construct the wrapper function so that its keyword arguments
  1288. # are visible in pydoc.help etc.
  1289. wrapper = """
  1290. def %(name)s(F, xin, iter=None %(kw)s, verbose=False, maxiter=None,
  1291. f_tol=None, f_rtol=None, x_tol=None, x_rtol=None,
  1292. tol_norm=None, line_search='armijo', callback=None, **kw):
  1293. jac = %(jac)s(%(kwkw)s **kw)
  1294. return nonlin_solve(F, xin, jac, iter, verbose, maxiter,
  1295. f_tol, f_rtol, x_tol, x_rtol, tol_norm, line_search,
  1296. callback)
  1297. """
  1298. wrapper = wrapper % dict(name=name, kw=kw_str, jac=jac.__name__,
  1299. kwkw=kwkw_str)
  1300. ns = {}
  1301. ns.update(globals())
  1302. exec(wrapper, ns)
  1303. func = ns[name]
  1304. func.__doc__ = jac.__doc__
  1305. _set_doc(func)
  1306. return func
  1307. broyden1 = _nonlin_wrapper('broyden1', BroydenFirst)
  1308. broyden2 = _nonlin_wrapper('broyden2', BroydenSecond)
  1309. anderson = _nonlin_wrapper('anderson', Anderson)
  1310. linearmixing = _nonlin_wrapper('linearmixing', LinearMixing)
  1311. diagbroyden = _nonlin_wrapper('diagbroyden', DiagBroyden)
  1312. excitingmixing = _nonlin_wrapper('excitingmixing', ExcitingMixing)
  1313. newton_krylov = _nonlin_wrapper('newton_krylov', KrylovJacobian)