_nonlin.py 48 KB

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