iterative.py 33 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051
  1. import warnings
  2. import numpy as np
  3. from scipy.sparse.linalg._interface import LinearOperator
  4. from .utils import make_system
  5. from scipy.linalg import get_lapack_funcs
  6. __all__ = ['bicg', 'bicgstab', 'cg', 'cgs', 'gmres', 'qmr']
  7. def _get_atol_rtol(name, b_norm, atol=0., rtol=1e-5):
  8. """
  9. A helper function to handle tolerance normalization
  10. """
  11. if atol == 'legacy' or atol is None or atol < 0:
  12. msg = (f"'scipy.sparse.linalg.{name}' called with invalid `atol`={atol}; "
  13. "if set, `atol` must be a real, non-negative number.")
  14. raise ValueError(msg)
  15. atol = max(float(atol), float(rtol) * float(b_norm))
  16. return atol, rtol
  17. def bicg(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M=None, callback=None):
  18. """
  19. Solve ``Ax = b`` with the BIConjugate Gradient method.
  20. Parameters
  21. ----------
  22. A : {sparse array, ndarray, LinearOperator}
  23. The real or complex N-by-N matrix of the linear system.
  24. Alternatively, `A` can be a linear operator which can
  25. produce ``Ax`` and ``A^T x`` using, e.g.,
  26. ``scipy.sparse.linalg.LinearOperator``.
  27. b : ndarray
  28. Right hand side of the linear system. Has shape (N,) or (N,1).
  29. x0 : ndarray
  30. Starting guess for the solution.
  31. rtol, atol : float, optional
  32. Parameters for the convergence test. For convergence,
  33. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  34. The default is ``atol=0.`` and ``rtol=1e-5``.
  35. maxiter : integer
  36. Maximum number of iterations. Iteration will stop after maxiter
  37. steps even if the specified tolerance has not been achieved.
  38. M : {sparse array, ndarray, LinearOperator}
  39. Preconditioner for `A`. It should approximate the
  40. inverse of `A` (see Notes). Effective preconditioning dramatically improves the
  41. rate of convergence, which implies that fewer iterations are needed
  42. to reach a given error tolerance.
  43. callback : function
  44. User-supplied function to call after each iteration. It is called
  45. as ``callback(xk)``, where ``xk`` is the current solution vector.
  46. Returns
  47. -------
  48. x : ndarray
  49. The converged solution.
  50. info : integer
  51. Provides convergence information:
  52. 0 : successful exit
  53. >0 : convergence to tolerance not achieved, number of iterations
  54. <0 : parameter breakdown
  55. Notes
  56. -----
  57. The preconditioner `M` should be a matrix such that ``M @ A`` has a smaller
  58. condition number than `A`, see [1]_ .
  59. References
  60. ----------
  61. .. [1] "Preconditioner", Wikipedia,
  62. https://en.wikipedia.org/wiki/Preconditioner
  63. .. [2] "Biconjugate gradient method", Wikipedia,
  64. https://en.wikipedia.org/wiki/Biconjugate_gradient_method
  65. Examples
  66. --------
  67. >>> import numpy as np
  68. >>> from scipy.sparse import csc_array
  69. >>> from scipy.sparse.linalg import bicg
  70. >>> A = csc_array([[3, 2, 0], [1, -1, 0], [0, 5, 1.]])
  71. >>> b = np.array([2., 4., -1.])
  72. >>> x, exitCode = bicg(A, b, atol=1e-5)
  73. >>> print(exitCode) # 0 indicates successful convergence
  74. 0
  75. >>> np.allclose(A.dot(x), b)
  76. True
  77. """
  78. A, M, x, b = make_system(A, M, x0, b)
  79. bnrm2 = np.linalg.norm(b)
  80. atol, _ = _get_atol_rtol('bicg', bnrm2, atol, rtol)
  81. if bnrm2 == 0:
  82. return b, 0
  83. n = len(b)
  84. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  85. if maxiter is None:
  86. maxiter = n*10
  87. matvec, rmatvec = A.matvec, A.rmatvec
  88. psolve, rpsolve = M.matvec, M.rmatvec
  89. rhotol = np.finfo(x.dtype.char).eps**2
  90. # Dummy values to initialize vars, silence linter warnings
  91. rho_prev, p, ptilde = None, None, None
  92. r = b - matvec(x) if x.any() else b.copy()
  93. rtilde = r.copy()
  94. for iteration in range(maxiter):
  95. if np.linalg.norm(r) < atol: # Are we done?
  96. return x, 0
  97. z = psolve(r)
  98. ztilde = rpsolve(rtilde)
  99. # order matters in this dot product
  100. rho_cur = dotprod(rtilde, z)
  101. if np.abs(rho_cur) < rhotol: # Breakdown case
  102. return x, -10
  103. if iteration > 0:
  104. beta = rho_cur / rho_prev
  105. p *= beta
  106. p += z
  107. ptilde *= beta.conj()
  108. ptilde += ztilde
  109. else: # First spin
  110. p = z.copy()
  111. ptilde = ztilde.copy()
  112. q = matvec(p)
  113. qtilde = rmatvec(ptilde)
  114. rv = dotprod(ptilde, q)
  115. if rv == 0:
  116. return x, -11
  117. alpha = rho_cur / rv
  118. x += alpha*p
  119. r -= alpha*q
  120. rtilde -= alpha.conj()*qtilde
  121. rho_prev = rho_cur
  122. if callback:
  123. callback(x)
  124. else: # for loop exhausted
  125. # Return incomplete progress
  126. return x, maxiter
  127. def bicgstab(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M=None,
  128. callback=None):
  129. """
  130. Solve ``Ax = b`` with the BIConjugate Gradient STABilized method.
  131. Parameters
  132. ----------
  133. A : {sparse array, ndarray, LinearOperator}
  134. The real or complex N-by-N matrix of the linear system.
  135. Alternatively, `A` can be a linear operator which can
  136. produce ``Ax`` and ``A^T x`` using, e.g.,
  137. ``scipy.sparse.linalg.LinearOperator``.
  138. b : ndarray
  139. Right hand side of the linear system. Has shape (N,) or (N,1).
  140. x0 : ndarray
  141. Starting guess for the solution.
  142. rtol, atol : float, optional
  143. Parameters for the convergence test. For convergence,
  144. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  145. The default is ``atol=0.`` and ``rtol=1e-5``.
  146. maxiter : integer
  147. Maximum number of iterations. Iteration will stop after maxiter
  148. steps even if the specified tolerance has not been achieved.
  149. M : {sparse array, ndarray, LinearOperator}
  150. Preconditioner for `A`. It should approximate the
  151. inverse of `A` (see Notes). Effective preconditioning dramatically improves the
  152. rate of convergence, which implies that fewer iterations are needed
  153. to reach a given error tolerance.
  154. callback : function
  155. User-supplied function to call after each iteration. It is called
  156. as ``callback(xk)``, where ``xk`` is the current solution vector.
  157. Returns
  158. -------
  159. x : ndarray
  160. The converged solution.
  161. info : integer
  162. Provides convergence information:
  163. 0 : successful exit
  164. >0 : convergence to tolerance not achieved, number of iterations
  165. <0 : parameter breakdown
  166. Notes
  167. -----
  168. The preconditioner `M` should be a matrix such that ``M @ A`` has a smaller
  169. condition number than `A`, see [1]_ .
  170. References
  171. ----------
  172. .. [1] "Preconditioner", Wikipedia,
  173. https://en.wikipedia.org/wiki/Preconditioner
  174. .. [2] "Biconjugate gradient stabilized method",
  175. Wikipedia, https://en.wikipedia.org/wiki/Biconjugate_gradient_stabilized_method
  176. Examples
  177. --------
  178. >>> import numpy as np
  179. >>> from scipy.sparse import csc_array
  180. >>> from scipy.sparse.linalg import bicgstab
  181. >>> R = np.array([[4, 2, 0, 1],
  182. ... [3, 0, 0, 2],
  183. ... [0, 1, 1, 1],
  184. ... [0, 2, 1, 0]])
  185. >>> A = csc_array(R)
  186. >>> b = np.array([-1, -0.5, -1, 2])
  187. >>> x, exit_code = bicgstab(A, b, atol=1e-5)
  188. >>> print(exit_code) # 0 indicates successful convergence
  189. 0
  190. >>> np.allclose(A.dot(x), b)
  191. True
  192. """
  193. A, M, x, b = make_system(A, M, x0, b)
  194. bnrm2 = np.linalg.norm(b)
  195. atol, _ = _get_atol_rtol('bicgstab', bnrm2, atol, rtol)
  196. if bnrm2 == 0:
  197. return b, 0
  198. n = len(b)
  199. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  200. if maxiter is None:
  201. maxiter = n*10
  202. matvec = A.matvec
  203. psolve = M.matvec
  204. # These values make no sense but coming from original Fortran code
  205. # sqrt might have been meant instead.
  206. rhotol = np.finfo(x.dtype.char).eps**2
  207. omegatol = rhotol
  208. # Dummy values to initialize vars, silence linter warnings
  209. rho_prev, omega, alpha, p, v = None, None, None, None, None
  210. r = b - matvec(x) if x.any() else b.copy()
  211. rtilde = r.copy()
  212. for iteration in range(maxiter):
  213. if np.linalg.norm(r) < atol: # Are we done?
  214. return x, 0
  215. rho = dotprod(rtilde, r)
  216. if np.abs(rho) < rhotol: # rho breakdown
  217. return x, -10
  218. if iteration > 0:
  219. if np.abs(omega) < omegatol: # omega breakdown
  220. return x, -11
  221. beta = (rho / rho_prev) * (alpha / omega)
  222. p -= omega*v
  223. p *= beta
  224. p += r
  225. else: # First spin
  226. s = np.empty_like(r)
  227. p = r.copy()
  228. phat = psolve(p)
  229. v = matvec(phat)
  230. rv = dotprod(rtilde, v)
  231. if rv == 0:
  232. return x, -11
  233. alpha = rho / rv
  234. r -= alpha*v
  235. s[:] = r[:]
  236. if np.linalg.norm(s) < atol:
  237. x += alpha*phat
  238. return x, 0
  239. shat = psolve(s)
  240. t = matvec(shat)
  241. omega = dotprod(t, s) / dotprod(t, t)
  242. x += alpha*phat
  243. x += omega*shat
  244. r -= omega*t
  245. rho_prev = rho
  246. if callback:
  247. callback(x)
  248. else: # for loop exhausted
  249. # Return incomplete progress
  250. return x, maxiter
  251. def cg(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M=None, callback=None):
  252. """
  253. Solve ``Ax = b`` with the Conjugate Gradient method, for a symmetric,
  254. positive-definite `A`.
  255. Parameters
  256. ----------
  257. A : {sparse array, ndarray, LinearOperator}
  258. The real or complex N-by-N matrix of the linear system.
  259. `A` must represent a hermitian, positive definite matrix.
  260. Alternatively, `A` can be a linear operator which can
  261. produce ``Ax`` using, e.g.,
  262. ``scipy.sparse.linalg.LinearOperator``.
  263. b : ndarray
  264. Right hand side of the linear system. Has shape (N,) or (N,1).
  265. x0 : ndarray
  266. Starting guess for the solution.
  267. rtol, atol : float, optional
  268. Parameters for the convergence test. For convergence,
  269. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  270. The default is ``atol=0.`` and ``rtol=1e-5``.
  271. maxiter : integer
  272. Maximum number of iterations. Iteration will stop after maxiter
  273. steps even if the specified tolerance has not been achieved.
  274. M : {sparse array, ndarray, LinearOperator}
  275. Preconditioner for `A`. `M` must represent a hermitian, positive definite
  276. matrix. It should approximate the inverse of `A` (see Notes).
  277. Effective preconditioning dramatically improves the
  278. rate of convergence, which implies that fewer iterations are needed
  279. to reach a given error tolerance.
  280. callback : function
  281. User-supplied function to call after each iteration. It is called
  282. as ``callback(xk)``, where ``xk`` is the current solution vector.
  283. Returns
  284. -------
  285. x : ndarray
  286. The converged solution.
  287. info : integer
  288. Provides convergence information:
  289. 0 : successful exit
  290. >0 : convergence to tolerance not achieved, number of iterations
  291. Notes
  292. -----
  293. The preconditioner `M` should be a matrix such that ``M @ A`` has a smaller
  294. condition number than `A`, see [2]_.
  295. References
  296. ----------
  297. .. [1] "Conjugate Gradient Method, Wikipedia,
  298. https://en.wikipedia.org/wiki/Conjugate_gradient_method
  299. .. [2] "Preconditioner",
  300. Wikipedia, https://en.wikipedia.org/wiki/Preconditioner
  301. Examples
  302. --------
  303. >>> import numpy as np
  304. >>> from scipy.sparse import csc_array
  305. >>> from scipy.sparse.linalg import cg
  306. >>> P = np.array([[4, 0, 1, 0],
  307. ... [0, 5, 0, 0],
  308. ... [1, 0, 3, 2],
  309. ... [0, 0, 2, 4]])
  310. >>> A = csc_array(P)
  311. >>> b = np.array([-1, -0.5, -1, 2])
  312. >>> x, exit_code = cg(A, b, atol=1e-5)
  313. >>> print(exit_code) # 0 indicates successful convergence
  314. 0
  315. >>> np.allclose(A.dot(x), b)
  316. True
  317. """
  318. A, M, x, b = make_system(A, M, x0, b)
  319. bnrm2 = np.linalg.norm(b)
  320. atol, _ = _get_atol_rtol('cg', bnrm2, atol, rtol)
  321. if bnrm2 == 0:
  322. return b, 0
  323. n = len(b)
  324. if maxiter is None:
  325. maxiter = n*10
  326. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  327. matvec = A.matvec
  328. psolve = M.matvec
  329. r = b - matvec(x) if x.any() else b.copy()
  330. # Dummy value to initialize var, silences warnings
  331. rho_prev, p = None, None
  332. for iteration in range(maxiter):
  333. if np.linalg.norm(r) < atol: # Are we done?
  334. return x, 0
  335. z = psolve(r)
  336. rho_cur = dotprod(r, z)
  337. if iteration > 0:
  338. beta = rho_cur / rho_prev
  339. p *= beta
  340. p += z
  341. else: # First spin
  342. p = np.empty_like(r)
  343. p[:] = z[:]
  344. q = matvec(p)
  345. alpha = rho_cur / dotprod(p, q)
  346. x += alpha*p
  347. r -= alpha*q
  348. rho_prev = rho_cur
  349. if callback:
  350. callback(x)
  351. else: # for loop exhausted
  352. # Return incomplete progress
  353. return x, maxiter
  354. def cgs(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M=None, callback=None):
  355. """
  356. Solve ``Ax = b`` with the Conjugate Gradient Squared method.
  357. Parameters
  358. ----------
  359. A : {sparse array, ndarray, LinearOperator}
  360. The real-valued N-by-N matrix of the linear system.
  361. Alternatively, `A` can be a linear operator which can
  362. produce ``Ax`` using, e.g.,
  363. ``scipy.sparse.linalg.LinearOperator``.
  364. b : ndarray
  365. Right hand side of the linear system. Has shape (N,) or (N,1).
  366. x0 : ndarray
  367. Starting guess for the solution.
  368. rtol, atol : float, optional
  369. Parameters for the convergence test. For convergence,
  370. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  371. The default is ``atol=0.`` and ``rtol=1e-5``.
  372. maxiter : integer
  373. Maximum number of iterations. Iteration will stop after maxiter
  374. steps even if the specified tolerance has not been achieved.
  375. M : {sparse array, ndarray, LinearOperator}
  376. Preconditioner for ``A``. It should approximate the
  377. inverse of `A` (see Notes). Effective preconditioning dramatically improves the
  378. rate of convergence, which implies that fewer iterations are needed
  379. to reach a given error tolerance.
  380. callback : function
  381. User-supplied function to call after each iteration. It is called
  382. as ``callback(xk)``, where ``xk`` is the current solution vector.
  383. Returns
  384. -------
  385. x : ndarray
  386. The converged solution.
  387. info : integer
  388. Provides convergence information:
  389. 0 : successful exit
  390. >0 : convergence to tolerance not achieved, number of iterations
  391. <0 : parameter breakdown
  392. Notes
  393. -----
  394. The preconditioner `M` should be a matrix such that ``M @ A`` has a smaller
  395. condition number than `A`, see [1]_.
  396. References
  397. ----------
  398. .. [1] "Preconditioner", Wikipedia,
  399. https://en.wikipedia.org/wiki/Preconditioner
  400. .. [2] "Conjugate gradient squared", Wikipedia,
  401. https://en.wikipedia.org/wiki/Conjugate_gradient_squared_method
  402. Examples
  403. --------
  404. >>> import numpy as np
  405. >>> from scipy.sparse import csc_array
  406. >>> from scipy.sparse.linalg import cgs
  407. >>> R = np.array([[4, 2, 0, 1],
  408. ... [3, 0, 0, 2],
  409. ... [0, 1, 1, 1],
  410. ... [0, 2, 1, 0]])
  411. >>> A = csc_array(R)
  412. >>> b = np.array([-1, -0.5, -1, 2])
  413. >>> x, exit_code = cgs(A, b)
  414. >>> print(exit_code) # 0 indicates successful convergence
  415. 0
  416. >>> np.allclose(A.dot(x), b)
  417. True
  418. """
  419. A, M, x, b = make_system(A, M, x0, b)
  420. bnrm2 = np.linalg.norm(b)
  421. atol, _ = _get_atol_rtol('cgs', bnrm2, atol, rtol)
  422. if bnrm2 == 0:
  423. return b, 0
  424. n = len(b)
  425. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  426. if maxiter is None:
  427. maxiter = n*10
  428. matvec = A.matvec
  429. psolve = M.matvec
  430. rhotol = np.finfo(x.dtype.char).eps**2
  431. r = b - matvec(x) if x.any() else b.copy()
  432. rtilde = r.copy()
  433. bnorm = np.linalg.norm(b)
  434. if bnorm == 0:
  435. bnorm = 1
  436. # Dummy values to initialize vars, silence linter warnings
  437. rho_prev, p, u, q = None, None, None, None
  438. for iteration in range(maxiter):
  439. rnorm = np.linalg.norm(r)
  440. if rnorm < atol: # Are we done?
  441. return x, 0
  442. rho_cur = dotprod(rtilde, r)
  443. if np.abs(rho_cur) < rhotol: # Breakdown case
  444. return x, -10
  445. if iteration > 0:
  446. beta = rho_cur / rho_prev
  447. # u = r + beta * q
  448. # p = u + beta * (q + beta * p);
  449. u[:] = r[:]
  450. u += beta*q
  451. p *= beta
  452. p += q
  453. p *= beta
  454. p += u
  455. else: # First spin
  456. p = r.copy()
  457. u = r.copy()
  458. q = np.empty_like(r)
  459. phat = psolve(p)
  460. vhat = matvec(phat)
  461. rv = dotprod(rtilde, vhat)
  462. if rv == 0: # Dot product breakdown
  463. return x, -11
  464. alpha = rho_cur / rv
  465. q[:] = u[:]
  466. q -= alpha*vhat
  467. uhat = psolve(u + q)
  468. x += alpha*uhat
  469. # Due to numerical error build-up the actual residual is computed
  470. # instead of the following two lines that were in the original
  471. # FORTRAN templates, still using a single matvec.
  472. # qhat = matvec(uhat)
  473. # r -= alpha*qhat
  474. r = b - matvec(x)
  475. rho_prev = rho_cur
  476. if callback:
  477. callback(x)
  478. else: # for loop exhausted
  479. # Return incomplete progress
  480. return x, maxiter
  481. def gmres(A, b, x0=None, *, rtol=1e-5, atol=0., restart=None, maxiter=None, M=None,
  482. callback=None, callback_type=None):
  483. """
  484. Solve ``Ax = b`` with the Generalized Minimal RESidual method.
  485. Parameters
  486. ----------
  487. A : {sparse array, ndarray, LinearOperator}
  488. The real or complex N-by-N matrix of the linear system.
  489. Alternatively, `A` can be a linear operator which can
  490. produce ``Ax`` using, e.g.,
  491. ``scipy.sparse.linalg.LinearOperator``.
  492. b : ndarray
  493. Right hand side of the linear system. Has shape (N,) or (N,1).
  494. x0 : ndarray
  495. Starting guess for the solution (a vector of zeros by default).
  496. atol, rtol : float
  497. Parameters for the convergence test. For convergence,
  498. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  499. The default is ``atol=0.`` and ``rtol=1e-5``.
  500. restart : int, optional
  501. Number of iterations between restarts. Larger values increase
  502. iteration cost, but may be necessary for convergence.
  503. If omitted, ``min(20, n)`` is used.
  504. maxiter : int, optional
  505. Maximum number of iterations (restart cycles). Iteration will stop
  506. after maxiter steps even if the specified tolerance has not been
  507. achieved. See `callback_type`.
  508. M : {sparse array, ndarray, LinearOperator}
  509. Inverse of the preconditioner of `A`. `M` should approximate the
  510. inverse of `A` and be easy to solve for (see Notes). Effective
  511. preconditioning dramatically improves the rate of convergence,
  512. which implies that fewer iterations are needed to reach a given
  513. error tolerance. By default, no preconditioner is used.
  514. In this implementation, left preconditioning is used,
  515. and the preconditioned residual is minimized. However, the final
  516. convergence is tested with respect to the ``b - A @ x`` residual.
  517. callback : function
  518. User-supplied function to call after each iteration. It is called
  519. as ``callback(args)``, where ``args`` are selected by `callback_type`.
  520. callback_type : {'x', 'pr_norm', 'legacy'}, optional
  521. Callback function argument requested:
  522. - ``x``: current iterate (ndarray), called on every restart
  523. - ``pr_norm``: relative (preconditioned) residual norm (float),
  524. called on every inner iteration
  525. - ``legacy`` (default): same as ``pr_norm``, but also changes the
  526. meaning of `maxiter` to count inner iterations instead of restart
  527. cycles.
  528. This keyword has no effect if `callback` is not set.
  529. Returns
  530. -------
  531. x : ndarray
  532. The converged solution.
  533. info : int
  534. Provides convergence information:
  535. 0 : successful exit
  536. >0 : convergence to tolerance not achieved, number of iterations
  537. See Also
  538. --------
  539. LinearOperator
  540. Notes
  541. -----
  542. A preconditioner, P, is chosen such that P is close to A but easy to solve
  543. for. The preconditioner parameter required by this routine is
  544. ``M = P^-1``. The inverse should preferably not be calculated
  545. explicitly. Rather, use the following template to produce M::
  546. # Construct a linear operator that computes P^-1 @ x.
  547. import scipy.sparse.linalg as spla
  548. M_x = lambda x: spla.spsolve(P, x)
  549. M = spla.LinearOperator((n, n), M_x)
  550. Examples
  551. --------
  552. >>> import numpy as np
  553. >>> from scipy.sparse import csc_array
  554. >>> from scipy.sparse.linalg import gmres
  555. >>> A = csc_array([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
  556. >>> b = np.array([2, 4, -1], dtype=float)
  557. >>> x, exitCode = gmres(A, b, atol=1e-5)
  558. >>> print(exitCode) # 0 indicates successful convergence
  559. 0
  560. >>> np.allclose(A.dot(x), b)
  561. True
  562. """
  563. if callback is not None and callback_type is None:
  564. # Warn about 'callback_type' semantic changes.
  565. # Probably should be removed only in far future, Scipy 2.0 or so.
  566. msg = ("scipy.sparse.linalg.gmres called without specifying "
  567. "`callback_type`. The default value will be changed in"
  568. " a future release. For compatibility, specify a value "
  569. "for `callback_type` explicitly, e.g., "
  570. "``gmres(..., callback_type='pr_norm')``, or to retain the "
  571. "old behavior ``gmres(..., callback_type='legacy')``"
  572. )
  573. warnings.warn(msg, category=DeprecationWarning, stacklevel=3)
  574. if callback_type is None:
  575. callback_type = 'legacy'
  576. if callback_type not in ('x', 'pr_norm', 'legacy'):
  577. raise ValueError(f"Unknown callback_type: {callback_type!r}")
  578. if callback is None:
  579. callback_type = None
  580. A, M, x, b = make_system(A, M, x0, b)
  581. matvec = A.matvec
  582. psolve = M.matvec
  583. n = len(b)
  584. bnrm2 = np.linalg.norm(b)
  585. atol, _ = _get_atol_rtol('gmres', bnrm2, atol, rtol)
  586. if bnrm2 == 0:
  587. return b, 0
  588. eps = np.finfo(x.dtype.char).eps
  589. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  590. if maxiter is None:
  591. maxiter = n*10
  592. if restart is None:
  593. restart = 20
  594. restart = min(restart, n)
  595. Mb_nrm2 = np.linalg.norm(psolve(b))
  596. # ====================================================
  597. # =========== Tolerance control from gh-8400 =========
  598. # ====================================================
  599. # Tolerance passed to GMRESREVCOM applies to the inner
  600. # iteration and deals with the left-preconditioned
  601. # residual.
  602. ptol_max_factor = 1.
  603. ptol = Mb_nrm2 * min(ptol_max_factor, atol / bnrm2)
  604. presid = 0.
  605. # ====================================================
  606. lartg = get_lapack_funcs('lartg', dtype=x.dtype)
  607. # allocate internal variables
  608. v = np.empty([restart+1, n], dtype=x.dtype)
  609. h = np.zeros([restart, restart+1], dtype=x.dtype)
  610. givens = np.zeros([restart, 2], dtype=x.dtype)
  611. # legacy iteration count
  612. inner_iter = 0
  613. for iteration in range(maxiter):
  614. if iteration == 0:
  615. r = b - matvec(x) if x.any() else b.copy()
  616. if np.linalg.norm(r) < atol: # Are we done?
  617. return x, 0
  618. v[0, :] = psolve(r)
  619. tmp = np.linalg.norm(v[0, :])
  620. v[0, :] *= (1 / tmp)
  621. # RHS of the Hessenberg problem
  622. S = np.zeros(restart+1, dtype=x.dtype)
  623. S[0] = tmp
  624. breakdown = False
  625. for col in range(restart):
  626. av = matvec(v[col, :])
  627. w = psolve(av)
  628. # Modified Gram-Schmidt
  629. h0 = np.linalg.norm(w)
  630. for k in range(col+1):
  631. tmp = dotprod(v[k, :], w)
  632. h[col, k] = tmp
  633. w -= tmp*v[k, :]
  634. h1 = np.linalg.norm(w)
  635. h[col, col + 1] = h1
  636. v[col + 1, :] = w[:]
  637. # Exact solution indicator
  638. if h1 <= eps*h0:
  639. h[col, col + 1] = 0
  640. breakdown = True
  641. else:
  642. v[col + 1, :] *= (1 / h1)
  643. # apply past Givens rotations to current h column
  644. for k in range(col):
  645. c, s = givens[k, 0], givens[k, 1]
  646. n0, n1 = h[col, [k, k+1]]
  647. h[col, [k, k + 1]] = [c*n0 + s*n1, -s.conj()*n0 + c*n1]
  648. # get and apply current rotation to h and S
  649. c, s, mag = lartg(h[col, col], h[col, col+1])
  650. givens[col, :] = [c, s]
  651. h[col, [col, col+1]] = mag, 0
  652. # S[col+1] component is always 0
  653. tmp = -np.conjugate(s)*S[col]
  654. S[[col, col + 1]] = [c*S[col], tmp]
  655. presid = np.abs(tmp)
  656. inner_iter += 1
  657. if callback_type in ('legacy', 'pr_norm'):
  658. callback(presid / bnrm2)
  659. # Legacy behavior
  660. if callback_type == 'legacy' and inner_iter == maxiter:
  661. break
  662. if presid <= ptol or breakdown:
  663. break
  664. # Solve h(col, col) upper triangular system and allow pseudo-solve
  665. # singular cases as in (but without the f2py copies):
  666. # y = trsv(h[:col+1, :col+1].T, S[:col+1])
  667. if h[col, col] == 0:
  668. S[col] = 0
  669. y = np.zeros([col+1], dtype=x.dtype)
  670. y[:] = S[:col+1]
  671. for k in range(col, 0, -1):
  672. if y[k] != 0:
  673. y[k] /= h[k, k]
  674. tmp = y[k]
  675. y[:k] -= tmp*h[k, :k]
  676. if y[0] != 0:
  677. y[0] /= h[0, 0]
  678. x += y @ v[:col+1, :]
  679. r = b - matvec(x)
  680. rnorm = np.linalg.norm(r)
  681. # Legacy exit
  682. if callback_type == 'legacy' and inner_iter == maxiter:
  683. return x, 0 if rnorm <= atol else maxiter
  684. if callback_type == 'x':
  685. callback(x)
  686. if rnorm <= atol:
  687. break
  688. elif breakdown:
  689. # Reached breakdown (= exact solution), but the external
  690. # tolerance check failed. Bail out with failure.
  691. break
  692. elif presid <= ptol:
  693. # Inner loop passed but outer didn't
  694. ptol_max_factor = max(eps, 0.25 * ptol_max_factor)
  695. else:
  696. ptol_max_factor = min(1.0, 1.5 * ptol_max_factor)
  697. ptol = presid * min(ptol_max_factor, atol / rnorm)
  698. info = 0 if (rnorm <= atol) else maxiter
  699. return x, info
  700. def qmr(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M1=None, M2=None,
  701. callback=None):
  702. """
  703. Solve ``Ax = b`` with the Quasi-Minimal Residual method.
  704. Parameters
  705. ----------
  706. A : {sparse array, ndarray, LinearOperator}
  707. The real-valued N-by-N matrix of the linear system.
  708. Alternatively, ``A`` can be a linear operator which can
  709. produce ``Ax`` and ``A^T x`` using, e.g.,
  710. ``scipy.sparse.linalg.LinearOperator``.
  711. b : ndarray
  712. Right hand side of the linear system. Has shape (N,) or (N,1).
  713. x0 : ndarray
  714. Starting guess for the solution.
  715. atol, rtol : float, optional
  716. Parameters for the convergence test. For convergence,
  717. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  718. The default is ``atol=0.`` and ``rtol=1e-5``.
  719. maxiter : integer
  720. Maximum number of iterations. Iteration will stop after maxiter
  721. steps even if the specified tolerance has not been achieved.
  722. M1 : {sparse array, ndarray, LinearOperator}
  723. Left preconditioner for A.
  724. M2 : {sparse array, ndarray, LinearOperator}
  725. Right preconditioner for A. Used together with the left
  726. preconditioner M1. The matrix M1@A@M2 should have better
  727. conditioned than A alone.
  728. callback : function
  729. User-supplied function to call after each iteration. It is called
  730. as callback(xk), where xk is the current solution vector.
  731. Returns
  732. -------
  733. x : ndarray
  734. The converged solution.
  735. info : integer
  736. Provides convergence information:
  737. 0 : successful exit
  738. >0 : convergence to tolerance not achieved, number of iterations
  739. <0 : parameter breakdown
  740. See Also
  741. --------
  742. LinearOperator
  743. Examples
  744. --------
  745. >>> import numpy as np
  746. >>> from scipy.sparse import csc_array
  747. >>> from scipy.sparse.linalg import qmr
  748. >>> A = csc_array([[3., 2., 0.], [1., -1., 0.], [0., 5., 1.]])
  749. >>> b = np.array([2., 4., -1.])
  750. >>> x, exitCode = qmr(A, b, atol=1e-5)
  751. >>> print(exitCode) # 0 indicates successful convergence
  752. 0
  753. >>> np.allclose(A.dot(x), b)
  754. True
  755. """
  756. A_ = A
  757. A, M, x, b = make_system(A, None, x0, b)
  758. bnrm2 = np.linalg.norm(b)
  759. atol, _ = _get_atol_rtol('qmr', bnrm2, atol, rtol)
  760. if bnrm2 == 0:
  761. return b, 0
  762. if M1 is None and M2 is None:
  763. if hasattr(A_, 'psolve'):
  764. def left_psolve(b):
  765. return A_.psolve(b, 'left')
  766. def right_psolve(b):
  767. return A_.psolve(b, 'right')
  768. def left_rpsolve(b):
  769. return A_.rpsolve(b, 'left')
  770. def right_rpsolve(b):
  771. return A_.rpsolve(b, 'right')
  772. M1 = LinearOperator(A.shape,
  773. matvec=left_psolve,
  774. rmatvec=left_rpsolve)
  775. M2 = LinearOperator(A.shape,
  776. matvec=right_psolve,
  777. rmatvec=right_rpsolve)
  778. else:
  779. def id(b):
  780. return b
  781. M1 = LinearOperator(A.shape, matvec=id, rmatvec=id)
  782. M2 = LinearOperator(A.shape, matvec=id, rmatvec=id)
  783. n = len(b)
  784. if maxiter is None:
  785. maxiter = n*10
  786. dotprod = np.vdot if np.iscomplexobj(x) else np.dot
  787. rhotol = np.finfo(x.dtype.char).eps
  788. betatol = rhotol
  789. gammatol = rhotol
  790. deltatol = rhotol
  791. epsilontol = rhotol
  792. xitol = rhotol
  793. r = b - A.matvec(x) if x.any() else b.copy()
  794. vtilde = r.copy()
  795. y = M1.matvec(vtilde)
  796. rho = np.linalg.norm(y)
  797. wtilde = r.copy()
  798. z = M2.rmatvec(wtilde)
  799. xi = np.linalg.norm(z)
  800. gamma, eta, theta = 1, -1, 0
  801. v = np.empty_like(vtilde)
  802. w = np.empty_like(wtilde)
  803. # Dummy values to initialize vars, silence linter warnings
  804. epsilon, q, d, p, s = None, None, None, None, None
  805. for iteration in range(maxiter):
  806. if np.linalg.norm(r) < atol: # Are we done?
  807. return x, 0
  808. if np.abs(rho) < rhotol: # rho breakdown
  809. return x, -10
  810. if np.abs(xi) < xitol: # xi breakdown
  811. return x, -15
  812. v[:] = vtilde[:]
  813. v *= (1 / rho)
  814. y *= (1 / rho)
  815. w[:] = wtilde[:]
  816. w *= (1 / xi)
  817. z *= (1 / xi)
  818. delta = dotprod(z, y)
  819. if np.abs(delta) < deltatol: # delta breakdown
  820. return x, -13
  821. ytilde = M2.matvec(y)
  822. ztilde = M1.rmatvec(z)
  823. if iteration > 0:
  824. ytilde -= (xi * delta / epsilon) * p
  825. p[:] = ytilde[:]
  826. ztilde -= (rho * (delta / epsilon).conj()) * q
  827. q[:] = ztilde[:]
  828. else: # First spin
  829. p = ytilde.copy()
  830. q = ztilde.copy()
  831. ptilde = A.matvec(p)
  832. epsilon = dotprod(q, ptilde)
  833. if np.abs(epsilon) < epsilontol: # epsilon breakdown
  834. return x, -14
  835. beta = epsilon / delta
  836. if np.abs(beta) < betatol: # beta breakdown
  837. return x, -11
  838. vtilde[:] = ptilde[:]
  839. vtilde -= beta*v
  840. y = M1.matvec(vtilde)
  841. rho_prev = rho
  842. rho = np.linalg.norm(y)
  843. wtilde[:] = w[:]
  844. wtilde *= - beta.conj()
  845. wtilde += A.rmatvec(q)
  846. z = M2.rmatvec(wtilde)
  847. xi = np.linalg.norm(z)
  848. gamma_prev = gamma
  849. theta_prev = theta
  850. theta = rho / (gamma_prev * np.abs(beta))
  851. gamma = 1 / np.sqrt(1 + theta**2)
  852. if np.abs(gamma) < gammatol: # gamma breakdown
  853. return x, -12
  854. eta *= -(rho_prev / beta) * (gamma / gamma_prev)**2
  855. if iteration > 0:
  856. d *= (theta_prev * gamma) ** 2
  857. d += eta*p
  858. s *= (theta_prev * gamma) ** 2
  859. s += eta*ptilde
  860. else:
  861. d = p.copy()
  862. d *= eta
  863. s = ptilde.copy()
  864. s *= eta
  865. x += d
  866. r -= s
  867. if callback:
  868. callback(x)
  869. else: # for loop exhausted
  870. # Return incomplete progress
  871. return x, maxiter