iterative.py 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753
  1. """Iterative methods for solving linear systems"""
  2. from __future__ import division, print_function, absolute_import
  3. __all__ = ['bicg','bicgstab','cg','cgs','gmres','qmr']
  4. import warnings
  5. import numpy as np
  6. from . import _iterative
  7. from scipy.sparse.linalg.interface import LinearOperator
  8. from .utils import make_system
  9. from scipy._lib._util import _aligned_zeros
  10. from scipy._lib._threadsafety import non_reentrant
  11. _type_conv = {'f':'s', 'd':'d', 'F':'c', 'D':'z'}
  12. # Part of the docstring common to all iterative solvers
  13. common_doc1 = \
  14. """
  15. Parameters
  16. ----------
  17. A : {sparse matrix, dense matrix, LinearOperator}"""
  18. common_doc2 = \
  19. """b : {array, matrix}
  20. Right hand side of the linear system. Has shape (N,) or (N,1).
  21. Returns
  22. -------
  23. x : {array, matrix}
  24. The converged solution.
  25. info : integer
  26. Provides convergence information:
  27. 0 : successful exit
  28. >0 : convergence to tolerance not achieved, number of iterations
  29. <0 : illegal input or breakdown
  30. Other Parameters
  31. ----------------
  32. x0 : {array, matrix}
  33. Starting guess for the solution.
  34. tol, atol : float, optional
  35. Tolerances for convergence, ``norm(residual) <= max(tol*norm(b), atol)``.
  36. The default for ``atol`` is ``'legacy'``, which emulates
  37. a different legacy behavior.
  38. .. warning::
  39. The default value for `atol` will be changed in a future release.
  40. For future compatibility, specify `atol` explicitly.
  41. maxiter : integer
  42. Maximum number of iterations. Iteration will stop after maxiter
  43. steps even if the specified tolerance has not been achieved.
  44. M : {sparse matrix, dense matrix, LinearOperator}
  45. Preconditioner for A. The preconditioner should approximate the
  46. inverse of A. Effective preconditioning dramatically improves the
  47. rate of convergence, which implies that fewer iterations are needed
  48. to reach a given error tolerance.
  49. callback : function
  50. User-supplied function to call after each iteration. It is called
  51. as callback(xk), where xk is the current solution vector.
  52. """
  53. def _stoptest(residual, atol):
  54. """
  55. Successful termination condition for the solvers.
  56. """
  57. resid = np.linalg.norm(residual)
  58. if resid <= atol:
  59. return resid, 1
  60. else:
  61. return resid, 0
  62. def _get_atol(tol, atol, bnrm2, get_residual, routine_name):
  63. """
  64. Parse arguments for absolute tolerance in termination condition.
  65. Parameters
  66. ----------
  67. tol, atol : object
  68. The arguments passed into the solver routine by user.
  69. bnrm2 : float
  70. 2-norm of the rhs vector.
  71. get_residual : callable
  72. Callable ``get_residual()`` that returns the initial value of
  73. the residual.
  74. routine_name : str
  75. Name of the routine.
  76. """
  77. if atol is None:
  78. warnings.warn("scipy.sparse.linalg.{name} called without specifying `atol`. "
  79. "The default value will be changed in a future release. "
  80. "For compatibility, specify a value for `atol` explicitly, e.g., "
  81. "``{name}(..., atol=0)``, or to retain the old behavior "
  82. "``{name}(..., atol='legacy')``".format(name=routine_name),
  83. category=DeprecationWarning, stacklevel=4)
  84. atol = 'legacy'
  85. tol = float(tol)
  86. if atol == 'legacy':
  87. # emulate old legacy behavior
  88. resid = get_residual()
  89. if resid <= tol:
  90. return 'exit'
  91. if bnrm2 == 0:
  92. return tol
  93. else:
  94. return tol * float(bnrm2)
  95. else:
  96. return max(float(atol), tol * float(bnrm2))
  97. def set_docstring(header, Ainfo, footer='', atol_default='0'):
  98. def combine(fn):
  99. fn.__doc__ = '\n'.join((header, common_doc1,
  100. ' ' + Ainfo.replace('\n', '\n '),
  101. common_doc2, footer))
  102. return fn
  103. return combine
  104. @set_docstring('Use BIConjugate Gradient iteration to solve ``Ax = b``.',
  105. 'The real or complex N-by-N matrix of the linear system.\n'
  106. 'It is required that the linear operator can produce\n'
  107. '``Ax`` and ``A^T x``.')
  108. @non_reentrant()
  109. def bicg(A, b, x0=None, tol=1e-5, maxiter=None, M=None, callback=None, atol=None):
  110. A,M,x,b,postprocess = make_system(A, M, x0, b)
  111. n = len(b)
  112. if maxiter is None:
  113. maxiter = n*10
  114. matvec, rmatvec = A.matvec, A.rmatvec
  115. psolve, rpsolve = M.matvec, M.rmatvec
  116. ltr = _type_conv[x.dtype.char]
  117. revcom = getattr(_iterative, ltr + 'bicgrevcom')
  118. get_residual = lambda: np.linalg.norm(matvec(x) - b)
  119. atol = _get_atol(tol, atol, np.linalg.norm(b), get_residual, 'bicg')
  120. if atol == 'exit':
  121. return postprocess(x), 0
  122. resid = atol
  123. ndx1 = 1
  124. ndx2 = -1
  125. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  126. work = _aligned_zeros(6*n,dtype=x.dtype)
  127. ijob = 1
  128. info = 0
  129. ftflag = True
  130. iter_ = maxiter
  131. while True:
  132. olditer = iter_
  133. x, iter_, resid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  134. revcom(b, x, work, iter_, resid, info, ndx1, ndx2, ijob)
  135. if callback is not None and iter_ > olditer:
  136. callback(x)
  137. slice1 = slice(ndx1-1, ndx1-1+n)
  138. slice2 = slice(ndx2-1, ndx2-1+n)
  139. if (ijob == -1):
  140. if callback is not None:
  141. callback(x)
  142. break
  143. elif (ijob == 1):
  144. work[slice2] *= sclr2
  145. work[slice2] += sclr1*matvec(work[slice1])
  146. elif (ijob == 2):
  147. work[slice2] *= sclr2
  148. work[slice2] += sclr1*rmatvec(work[slice1])
  149. elif (ijob == 3):
  150. work[slice1] = psolve(work[slice2])
  151. elif (ijob == 4):
  152. work[slice1] = rpsolve(work[slice2])
  153. elif (ijob == 5):
  154. work[slice2] *= sclr2
  155. work[slice2] += sclr1*matvec(x)
  156. elif (ijob == 6):
  157. if ftflag:
  158. info = -1
  159. ftflag = False
  160. resid, info = _stoptest(work[slice1], atol)
  161. ijob = 2
  162. if info > 0 and iter_ == maxiter and not (resid <= atol):
  163. # info isn't set appropriately otherwise
  164. info = iter_
  165. return postprocess(x), info
  166. @set_docstring('Use BIConjugate Gradient STABilized iteration to solve '
  167. '``Ax = b``.',
  168. 'The real or complex N-by-N matrix of the linear system.')
  169. @non_reentrant()
  170. def bicgstab(A, b, x0=None, tol=1e-5, maxiter=None, M=None, callback=None, atol=None):
  171. A, M, x, b, postprocess = make_system(A, M, x0, b)
  172. n = len(b)
  173. if maxiter is None:
  174. maxiter = n*10
  175. matvec = A.matvec
  176. psolve = M.matvec
  177. ltr = _type_conv[x.dtype.char]
  178. revcom = getattr(_iterative, ltr + 'bicgstabrevcom')
  179. get_residual = lambda: np.linalg.norm(matvec(x) - b)
  180. atol = _get_atol(tol, atol, np.linalg.norm(b), get_residual, 'bicgstab')
  181. if atol == 'exit':
  182. return postprocess(x), 0
  183. resid = atol
  184. ndx1 = 1
  185. ndx2 = -1
  186. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  187. work = _aligned_zeros(7*n,dtype=x.dtype)
  188. ijob = 1
  189. info = 0
  190. ftflag = True
  191. iter_ = maxiter
  192. while True:
  193. olditer = iter_
  194. x, iter_, resid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  195. revcom(b, x, work, iter_, resid, info, ndx1, ndx2, ijob)
  196. if callback is not None and iter_ > olditer:
  197. callback(x)
  198. slice1 = slice(ndx1-1, ndx1-1+n)
  199. slice2 = slice(ndx2-1, ndx2-1+n)
  200. if (ijob == -1):
  201. if callback is not None:
  202. callback(x)
  203. break
  204. elif (ijob == 1):
  205. work[slice2] *= sclr2
  206. work[slice2] += sclr1*matvec(work[slice1])
  207. elif (ijob == 2):
  208. work[slice1] = psolve(work[slice2])
  209. elif (ijob == 3):
  210. work[slice2] *= sclr2
  211. work[slice2] += sclr1*matvec(x)
  212. elif (ijob == 4):
  213. if ftflag:
  214. info = -1
  215. ftflag = False
  216. resid, info = _stoptest(work[slice1], atol)
  217. ijob = 2
  218. if info > 0 and iter_ == maxiter and not (resid <= atol):
  219. # info isn't set appropriately otherwise
  220. info = iter_
  221. return postprocess(x), info
  222. @set_docstring('Use Conjugate Gradient iteration to solve ``Ax = b``.',
  223. 'The real or complex N-by-N matrix of the linear system.\n'
  224. '``A`` must represent a hermitian, positive definite matrix.')
  225. @non_reentrant()
  226. def cg(A, b, x0=None, tol=1e-5, maxiter=None, M=None, callback=None, atol=None):
  227. A, M, x, b, postprocess = make_system(A, M, x0, b)
  228. n = len(b)
  229. if maxiter is None:
  230. maxiter = n*10
  231. matvec = A.matvec
  232. psolve = M.matvec
  233. ltr = _type_conv[x.dtype.char]
  234. revcom = getattr(_iterative, ltr + 'cgrevcom')
  235. get_residual = lambda: np.linalg.norm(matvec(x) - b)
  236. atol = _get_atol(tol, atol, np.linalg.norm(b), get_residual, 'cg')
  237. if atol == 'exit':
  238. return postprocess(x), 0
  239. resid = atol
  240. ndx1 = 1
  241. ndx2 = -1
  242. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  243. work = _aligned_zeros(4*n,dtype=x.dtype)
  244. ijob = 1
  245. info = 0
  246. ftflag = True
  247. iter_ = maxiter
  248. while True:
  249. olditer = iter_
  250. x, iter_, resid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  251. revcom(b, x, work, iter_, resid, info, ndx1, ndx2, ijob)
  252. if callback is not None and iter_ > olditer:
  253. callback(x)
  254. slice1 = slice(ndx1-1, ndx1-1+n)
  255. slice2 = slice(ndx2-1, ndx2-1+n)
  256. if (ijob == -1):
  257. if callback is not None:
  258. callback(x)
  259. break
  260. elif (ijob == 1):
  261. work[slice2] *= sclr2
  262. work[slice2] += sclr1*matvec(work[slice1])
  263. elif (ijob == 2):
  264. work[slice1] = psolve(work[slice2])
  265. elif (ijob == 3):
  266. work[slice2] *= sclr2
  267. work[slice2] += sclr1*matvec(x)
  268. elif (ijob == 4):
  269. if ftflag:
  270. info = -1
  271. ftflag = False
  272. resid, info = _stoptest(work[slice1], atol)
  273. if info == 1 and iter_ > 1:
  274. # recompute residual and recheck, to avoid
  275. # accumulating rounding error
  276. work[slice1] = b - matvec(x)
  277. resid, info = _stoptest(work[slice1], atol)
  278. ijob = 2
  279. if info > 0 and iter_ == maxiter and not (resid <= atol):
  280. # info isn't set appropriately otherwise
  281. info = iter_
  282. return postprocess(x), info
  283. @set_docstring('Use Conjugate Gradient Squared iteration to solve ``Ax = b``.',
  284. 'The real-valued N-by-N matrix of the linear system.')
  285. @non_reentrant()
  286. def cgs(A, b, x0=None, tol=1e-5, maxiter=None, M=None, callback=None, atol=None):
  287. A, M, x, b, postprocess = make_system(A, M, x0, b)
  288. n = len(b)
  289. if maxiter is None:
  290. maxiter = n*10
  291. matvec = A.matvec
  292. psolve = M.matvec
  293. ltr = _type_conv[x.dtype.char]
  294. revcom = getattr(_iterative, ltr + 'cgsrevcom')
  295. get_residual = lambda: np.linalg.norm(matvec(x) - b)
  296. atol = _get_atol(tol, atol, np.linalg.norm(b), get_residual, 'cgs')
  297. if atol == 'exit':
  298. return postprocess(x), 0
  299. resid = atol
  300. ndx1 = 1
  301. ndx2 = -1
  302. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  303. work = _aligned_zeros(7*n,dtype=x.dtype)
  304. ijob = 1
  305. info = 0
  306. ftflag = True
  307. iter_ = maxiter
  308. while True:
  309. olditer = iter_
  310. x, iter_, resid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  311. revcom(b, x, work, iter_, resid, info, ndx1, ndx2, ijob)
  312. if callback is not None and iter_ > olditer:
  313. callback(x)
  314. slice1 = slice(ndx1-1, ndx1-1+n)
  315. slice2 = slice(ndx2-1, ndx2-1+n)
  316. if (ijob == -1):
  317. if callback is not None:
  318. callback(x)
  319. break
  320. elif (ijob == 1):
  321. work[slice2] *= sclr2
  322. work[slice2] += sclr1*matvec(work[slice1])
  323. elif (ijob == 2):
  324. work[slice1] = psolve(work[slice2])
  325. elif (ijob == 3):
  326. work[slice2] *= sclr2
  327. work[slice2] += sclr1*matvec(x)
  328. elif (ijob == 4):
  329. if ftflag:
  330. info = -1
  331. ftflag = False
  332. resid, info = _stoptest(work[slice1], atol)
  333. if info == 1 and iter_ > 1:
  334. # recompute residual and recheck, to avoid
  335. # accumulating rounding error
  336. work[slice1] = b - matvec(x)
  337. resid, info = _stoptest(work[slice1], atol)
  338. ijob = 2
  339. if info == -10:
  340. # termination due to breakdown: check for convergence
  341. resid, ok = _stoptest(b - matvec(x), atol)
  342. if ok:
  343. info = 0
  344. if info > 0 and iter_ == maxiter and not (resid <= atol):
  345. # info isn't set appropriately otherwise
  346. info = iter_
  347. return postprocess(x), info
  348. @non_reentrant()
  349. def gmres(A, b, x0=None, tol=1e-5, restart=None, maxiter=None, M=None, callback=None,
  350. restrt=None, atol=None):
  351. """
  352. Use Generalized Minimal RESidual iteration to solve ``Ax = b``.
  353. Parameters
  354. ----------
  355. A : {sparse matrix, dense matrix, LinearOperator}
  356. The real or complex N-by-N matrix of the linear system.
  357. b : {array, matrix}
  358. Right hand side of the linear system. Has shape (N,) or (N,1).
  359. Returns
  360. -------
  361. x : {array, matrix}
  362. The converged solution.
  363. info : int
  364. Provides convergence information:
  365. * 0 : successful exit
  366. * >0 : convergence to tolerance not achieved, number of iterations
  367. * <0 : illegal input or breakdown
  368. Other parameters
  369. ----------------
  370. x0 : {array, matrix}
  371. Starting guess for the solution (a vector of zeros by default).
  372. tol, atol : float, optional
  373. Tolerances for convergence, ``norm(residual) <= max(tol*norm(b), atol)``.
  374. The default for ``atol`` is ``'legacy'``, which emulates
  375. a different legacy behavior.
  376. .. warning::
  377. The default value for `atol` will be changed in a future release.
  378. For future compatibility, specify `atol` explicitly.
  379. restart : int, optional
  380. Number of iterations between restarts. Larger values increase
  381. iteration cost, but may be necessary for convergence.
  382. Default is 20.
  383. maxiter : int, optional
  384. Maximum number of iterations (restart cycles). Iteration will stop
  385. after maxiter steps even if the specified tolerance has not been
  386. achieved.
  387. M : {sparse matrix, dense matrix, LinearOperator}
  388. Inverse of the preconditioner of A. M should approximate the
  389. inverse of A and be easy to solve for (see Notes). Effective
  390. preconditioning dramatically improves the rate of convergence,
  391. which implies that fewer iterations are needed to reach a given
  392. error tolerance. By default, no preconditioner is used.
  393. callback : function
  394. User-supplied function to call after each iteration. It is called
  395. as callback(rk), where rk is the current residual vector.
  396. restrt : int, optional
  397. DEPRECATED - use `restart` instead.
  398. See Also
  399. --------
  400. LinearOperator
  401. Notes
  402. -----
  403. A preconditioner, P, is chosen such that P is close to A but easy to solve
  404. for. The preconditioner parameter required by this routine is
  405. ``M = P^-1``. The inverse should preferably not be calculated
  406. explicitly. Rather, use the following template to produce M::
  407. # Construct a linear operator that computes P^-1 * x.
  408. import scipy.sparse.linalg as spla
  409. M_x = lambda x: spla.spsolve(P, x)
  410. M = spla.LinearOperator((n, n), M_x)
  411. Examples
  412. --------
  413. >>> from scipy.sparse import csc_matrix
  414. >>> from scipy.sparse.linalg import gmres
  415. >>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
  416. >>> b = np.array([2, 4, -1], dtype=float)
  417. >>> x, exitCode = gmres(A, b)
  418. >>> print(exitCode) # 0 indicates successful convergence
  419. 0
  420. >>> np.allclose(A.dot(x), b)
  421. True
  422. """
  423. # Change 'restrt' keyword to 'restart'
  424. if restrt is None:
  425. restrt = restart
  426. elif restart is not None:
  427. raise ValueError("Cannot specify both restart and restrt keywords. "
  428. "Preferably use 'restart' only.")
  429. A, M, x, b,postprocess = make_system(A, M, x0, b)
  430. n = len(b)
  431. if maxiter is None:
  432. maxiter = n*10
  433. if restrt is None:
  434. restrt = 20
  435. restrt = min(restrt, n)
  436. matvec = A.matvec
  437. psolve = M.matvec
  438. ltr = _type_conv[x.dtype.char]
  439. revcom = getattr(_iterative, ltr + 'gmresrevcom')
  440. bnrm2 = np.linalg.norm(b)
  441. Mb_nrm2 = np.linalg.norm(psolve(b))
  442. get_residual = lambda: np.linalg.norm(matvec(x) - b)
  443. atol = _get_atol(tol, atol, bnrm2, get_residual, 'gmres')
  444. if atol == 'exit':
  445. return postprocess(x), 0
  446. if bnrm2 == 0:
  447. return postprocess(b), 0
  448. # Tolerance passed to GMRESREVCOM applies to the inner iteration
  449. # and deals with the left-preconditioned residual.
  450. ptol_max_factor = 1.0
  451. ptol = Mb_nrm2 * min(ptol_max_factor, atol / bnrm2)
  452. resid = np.nan
  453. presid = np.nan
  454. ndx1 = 1
  455. ndx2 = -1
  456. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  457. work = _aligned_zeros((6+restrt)*n,dtype=x.dtype)
  458. work2 = _aligned_zeros((restrt+1)*(2*restrt+2),dtype=x.dtype)
  459. ijob = 1
  460. info = 0
  461. ftflag = True
  462. iter_ = maxiter
  463. old_ijob = ijob
  464. first_pass = True
  465. resid_ready = False
  466. iter_num = 1
  467. while True:
  468. x, iter_, presid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  469. revcom(b, x, restrt, work, work2, iter_, presid, info, ndx1, ndx2, ijob, ptol)
  470. slice1 = slice(ndx1-1, ndx1-1+n)
  471. slice2 = slice(ndx2-1, ndx2-1+n)
  472. if (ijob == -1): # gmres success, update last residual
  473. if resid_ready and callback is not None:
  474. callback(presid / bnrm2)
  475. resid_ready = False
  476. break
  477. elif (ijob == 1):
  478. work[slice2] *= sclr2
  479. work[slice2] += sclr1*matvec(x)
  480. elif (ijob == 2):
  481. work[slice1] = psolve(work[slice2])
  482. if not first_pass and old_ijob == 3:
  483. resid_ready = True
  484. first_pass = False
  485. elif (ijob == 3):
  486. work[slice2] *= sclr2
  487. work[slice2] += sclr1*matvec(work[slice1])
  488. if resid_ready and callback is not None:
  489. callback(presid / bnrm2)
  490. resid_ready = False
  491. iter_num = iter_num+1
  492. elif (ijob == 4):
  493. if ftflag:
  494. info = -1
  495. ftflag = False
  496. resid, info = _stoptest(work[slice1], atol)
  497. # Inner loop tolerance control
  498. if info or presid > ptol:
  499. ptol_max_factor = min(1.0, 1.5 * ptol_max_factor)
  500. else:
  501. # Inner loop tolerance OK, but outer loop not.
  502. ptol_max_factor = max(1e-16, 0.25 * ptol_max_factor)
  503. if resid != 0:
  504. ptol = presid * min(ptol_max_factor, atol / resid)
  505. else:
  506. ptol = presid * ptol_max_factor
  507. old_ijob = ijob
  508. ijob = 2
  509. if iter_num > maxiter:
  510. info = maxiter
  511. break
  512. if info >= 0 and not (resid <= atol):
  513. # info isn't set appropriately otherwise
  514. info = maxiter
  515. return postprocess(x), info
  516. @non_reentrant()
  517. def qmr(A, b, x0=None, tol=1e-5, maxiter=None, M1=None, M2=None, callback=None,
  518. atol=None):
  519. """Use Quasi-Minimal Residual iteration to solve ``Ax = b``.
  520. Parameters
  521. ----------
  522. A : {sparse matrix, dense matrix, LinearOperator}
  523. The real-valued N-by-N matrix of the linear system.
  524. It is required that the linear operator can produce
  525. ``Ax`` and ``A^T x``.
  526. b : {array, matrix}
  527. Right hand side of the linear system. Has shape (N,) or (N,1).
  528. Returns
  529. -------
  530. x : {array, matrix}
  531. The converged solution.
  532. info : integer
  533. Provides convergence information:
  534. 0 : successful exit
  535. >0 : convergence to tolerance not achieved, number of iterations
  536. <0 : illegal input or breakdown
  537. Other Parameters
  538. ----------------
  539. x0 : {array, matrix}
  540. Starting guess for the solution.
  541. tol, atol : float, optional
  542. Tolerances for convergence, ``norm(residual) <= max(tol*norm(b), atol)``.
  543. The default for ``atol`` is ``'legacy'``, which emulates
  544. a different legacy behavior.
  545. .. warning::
  546. The default value for `atol` will be changed in a future release.
  547. For future compatibility, specify `atol` explicitly.
  548. maxiter : integer
  549. Maximum number of iterations. Iteration will stop after maxiter
  550. steps even if the specified tolerance has not been achieved.
  551. M1 : {sparse matrix, dense matrix, LinearOperator}
  552. Left preconditioner for A.
  553. M2 : {sparse matrix, dense matrix, LinearOperator}
  554. Right preconditioner for A. Used together with the left
  555. preconditioner M1. The matrix M1*A*M2 should have better
  556. conditioned than A alone.
  557. callback : function
  558. User-supplied function to call after each iteration. It is called
  559. as callback(xk), where xk is the current solution vector.
  560. See Also
  561. --------
  562. LinearOperator
  563. Examples
  564. --------
  565. >>> from scipy.sparse import csc_matrix
  566. >>> from scipy.sparse.linalg import qmr
  567. >>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
  568. >>> b = np.array([2, 4, -1], dtype=float)
  569. >>> x, exitCode = qmr(A, b)
  570. >>> print(exitCode) # 0 indicates successful convergence
  571. 0
  572. >>> np.allclose(A.dot(x), b)
  573. True
  574. """
  575. A_ = A
  576. A, M, x, b, postprocess = make_system(A, None, x0, b)
  577. if M1 is None and M2 is None:
  578. if hasattr(A_,'psolve'):
  579. def left_psolve(b):
  580. return A_.psolve(b,'left')
  581. def right_psolve(b):
  582. return A_.psolve(b,'right')
  583. def left_rpsolve(b):
  584. return A_.rpsolve(b,'left')
  585. def right_rpsolve(b):
  586. return A_.rpsolve(b,'right')
  587. M1 = LinearOperator(A.shape, matvec=left_psolve, rmatvec=left_rpsolve)
  588. M2 = LinearOperator(A.shape, matvec=right_psolve, rmatvec=right_rpsolve)
  589. else:
  590. def id(b):
  591. return b
  592. M1 = LinearOperator(A.shape, matvec=id, rmatvec=id)
  593. M2 = LinearOperator(A.shape, matvec=id, rmatvec=id)
  594. n = len(b)
  595. if maxiter is None:
  596. maxiter = n*10
  597. ltr = _type_conv[x.dtype.char]
  598. revcom = getattr(_iterative, ltr + 'qmrrevcom')
  599. get_residual = lambda: np.linalg.norm(A.matvec(x) - b)
  600. atol = _get_atol(tol, atol, np.linalg.norm(b), get_residual, 'qmr')
  601. if atol == 'exit':
  602. return postprocess(x), 0
  603. resid = atol
  604. ndx1 = 1
  605. ndx2 = -1
  606. # Use _aligned_zeros to work around a f2py bug in Numpy 1.9.1
  607. work = _aligned_zeros(11*n,x.dtype)
  608. ijob = 1
  609. info = 0
  610. ftflag = True
  611. iter_ = maxiter
  612. while True:
  613. olditer = iter_
  614. x, iter_, resid, info, ndx1, ndx2, sclr1, sclr2, ijob = \
  615. revcom(b, x, work, iter_, resid, info, ndx1, ndx2, ijob)
  616. if callback is not None and iter_ > olditer:
  617. callback(x)
  618. slice1 = slice(ndx1-1, ndx1-1+n)
  619. slice2 = slice(ndx2-1, ndx2-1+n)
  620. if (ijob == -1):
  621. if callback is not None:
  622. callback(x)
  623. break
  624. elif (ijob == 1):
  625. work[slice2] *= sclr2
  626. work[slice2] += sclr1*A.matvec(work[slice1])
  627. elif (ijob == 2):
  628. work[slice2] *= sclr2
  629. work[slice2] += sclr1*A.rmatvec(work[slice1])
  630. elif (ijob == 3):
  631. work[slice1] = M1.matvec(work[slice2])
  632. elif (ijob == 4):
  633. work[slice1] = M2.matvec(work[slice2])
  634. elif (ijob == 5):
  635. work[slice1] = M1.rmatvec(work[slice2])
  636. elif (ijob == 6):
  637. work[slice1] = M2.rmatvec(work[slice2])
  638. elif (ijob == 7):
  639. work[slice2] *= sclr2
  640. work[slice2] += sclr1*A.matvec(x)
  641. elif (ijob == 8):
  642. if ftflag:
  643. info = -1
  644. ftflag = False
  645. resid, info = _stoptest(work[slice1], atol)
  646. ijob = 2
  647. if info > 0 and iter_ == maxiter and not (resid <= atol):
  648. # info isn't set appropriately otherwise
  649. info = iter_
  650. return postprocess(x), info