lgmres.py 8.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237
  1. # Copyright (C) 2009, Pauli Virtanen <pav@iki.fi>
  2. # Distributed under the same license as SciPy.
  3. import warnings
  4. import numpy as np
  5. from numpy.linalg import LinAlgError
  6. from scipy.linalg import get_blas_funcs
  7. from .utils import make_system
  8. from ._gcrotmk import _fgmres
  9. __all__ = ['lgmres']
  10. def lgmres(A, b, x0=None, tol=1e-5, maxiter=1000, M=None, callback=None,
  11. inner_m=30, outer_k=3, outer_v=None, store_outer_Av=True,
  12. prepend_outer_v=False, atol=None):
  13. """
  14. Solve a matrix equation using the LGMRES algorithm.
  15. The LGMRES algorithm [1]_ [2]_ is designed to avoid some problems
  16. in the convergence in restarted GMRES, and often converges in fewer
  17. iterations.
  18. Parameters
  19. ----------
  20. A : {sparse matrix, ndarray, LinearOperator}
  21. The real or complex N-by-N matrix of the linear system.
  22. Alternatively, ``A`` can be a linear operator which can
  23. produce ``Ax`` using, e.g.,
  24. ``scipy.sparse.linalg.LinearOperator``.
  25. b : ndarray
  26. Right hand side of the linear system. Has shape (N,) or (N,1).
  27. x0 : ndarray
  28. Starting guess for the solution.
  29. tol, atol : float, optional
  30. Tolerances for convergence, ``norm(residual) <= max(tol*norm(b), atol)``.
  31. The default for ``atol`` is `tol`.
  32. .. warning::
  33. The default value for `atol` will be changed in a future release.
  34. For future compatibility, specify `atol` explicitly.
  35. maxiter : int, optional
  36. Maximum number of iterations. Iteration will stop after maxiter
  37. steps even if the specified tolerance has not been achieved.
  38. M : {sparse matrix, ndarray, LinearOperator}, optional
  39. Preconditioner for A. The preconditioner should approximate the
  40. inverse of A. 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, optional
  44. User-supplied function to call after each iteration. It is called
  45. as callback(xk), where xk is the current solution vector.
  46. inner_m : int, optional
  47. Number of inner GMRES iterations per each outer iteration.
  48. outer_k : int, optional
  49. Number of vectors to carry between inner GMRES iterations.
  50. According to [1]_, good values are in the range of 1...3.
  51. However, note that if you want to use the additional vectors to
  52. accelerate solving multiple similar problems, larger values may
  53. be beneficial.
  54. outer_v : list of tuples, optional
  55. List containing tuples ``(v, Av)`` of vectors and corresponding
  56. matrix-vector products, used to augment the Krylov subspace, and
  57. carried between inner GMRES iterations. The element ``Av`` can
  58. be `None` if the matrix-vector product should be re-evaluated.
  59. This parameter is modified in-place by `lgmres`, and can be used
  60. to pass "guess" vectors in and out of the algorithm when solving
  61. similar problems.
  62. store_outer_Av : bool, optional
  63. Whether LGMRES should store also A@v in addition to vectors `v`
  64. in the `outer_v` list. Default is True.
  65. prepend_outer_v : bool, optional
  66. Whether to put outer_v augmentation vectors before Krylov iterates.
  67. In standard LGMRES, prepend_outer_v=False.
  68. Returns
  69. -------
  70. x : ndarray
  71. The converged solution.
  72. info : int
  73. Provides convergence information:
  74. - 0 : successful exit
  75. - >0 : convergence to tolerance not achieved, number of iterations
  76. - <0 : illegal input or breakdown
  77. Notes
  78. -----
  79. The LGMRES algorithm [1]_ [2]_ is designed to avoid the
  80. slowing of convergence in restarted GMRES, due to alternating
  81. residual vectors. Typically, it often outperforms GMRES(m) of
  82. comparable memory requirements by some measure, or at least is not
  83. much worse.
  84. Another advantage in this algorithm is that you can supply it with
  85. 'guess' vectors in the `outer_v` argument that augment the Krylov
  86. subspace. If the solution lies close to the span of these vectors,
  87. the algorithm converges faster. This can be useful if several very
  88. similar matrices need to be inverted one after another, such as in
  89. Newton-Krylov iteration where the Jacobian matrix often changes
  90. little in the nonlinear steps.
  91. References
  92. ----------
  93. .. [1] A.H. Baker and E.R. Jessup and T. Manteuffel, "A Technique for
  94. Accelerating the Convergence of Restarted GMRES", SIAM J. Matrix
  95. Anal. Appl. 26, 962 (2005).
  96. .. [2] A.H. Baker, "On Improving the Performance of the Linear Solver
  97. restarted GMRES", PhD thesis, University of Colorado (2003).
  98. Examples
  99. --------
  100. >>> import numpy as np
  101. >>> from scipy.sparse import csc_matrix
  102. >>> from scipy.sparse.linalg import lgmres
  103. >>> A = csc_matrix([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
  104. >>> b = np.array([2, 4, -1], dtype=float)
  105. >>> x, exitCode = lgmres(A, b)
  106. >>> print(exitCode) # 0 indicates successful convergence
  107. 0
  108. >>> np.allclose(A.dot(x), b)
  109. True
  110. """
  111. A,M,x,b,postprocess = make_system(A,M,x0,b)
  112. if not np.isfinite(b).all():
  113. raise ValueError("RHS must contain only finite numbers")
  114. if atol is None:
  115. warnings.warn("scipy.sparse.linalg.lgmres called without specifying `atol`. "
  116. "The default value will change in the future. To preserve "
  117. "current behavior, set ``atol=tol``.",
  118. category=DeprecationWarning, stacklevel=2)
  119. atol = tol
  120. matvec = A.matvec
  121. psolve = M.matvec
  122. if outer_v is None:
  123. outer_v = []
  124. axpy, dot, scal = None, None, None
  125. nrm2 = get_blas_funcs('nrm2', [b])
  126. b_norm = nrm2(b)
  127. if b_norm == 0:
  128. x = b
  129. return (postprocess(x), 0)
  130. ptol_max_factor = 1.0
  131. for k_outer in range(maxiter):
  132. r_outer = matvec(x) - b
  133. # -- callback
  134. if callback is not None:
  135. callback(x)
  136. # -- determine input type routines
  137. if axpy is None:
  138. if np.iscomplexobj(r_outer) and not np.iscomplexobj(x):
  139. x = x.astype(r_outer.dtype)
  140. axpy, dot, scal, nrm2 = get_blas_funcs(['axpy', 'dot', 'scal', 'nrm2'],
  141. (x, r_outer))
  142. # -- check stopping condition
  143. r_norm = nrm2(r_outer)
  144. if r_norm <= max(atol, tol * b_norm):
  145. break
  146. # -- inner LGMRES iteration
  147. v0 = -psolve(r_outer)
  148. inner_res_0 = nrm2(v0)
  149. if inner_res_0 == 0:
  150. rnorm = nrm2(r_outer)
  151. raise RuntimeError("Preconditioner returned a zero vector; "
  152. "|v| ~ %.1g, |M v| = 0" % rnorm)
  153. v0 = scal(1.0/inner_res_0, v0)
  154. ptol = min(ptol_max_factor, max(atol, tol*b_norm)/r_norm)
  155. try:
  156. Q, R, B, vs, zs, y, pres = _fgmres(matvec,
  157. v0,
  158. inner_m,
  159. lpsolve=psolve,
  160. atol=ptol,
  161. outer_v=outer_v,
  162. prepend_outer_v=prepend_outer_v)
  163. y *= inner_res_0
  164. if not np.isfinite(y).all():
  165. # Overflow etc. in computation. There's no way to
  166. # recover from this, so we have to bail out.
  167. raise LinAlgError()
  168. except LinAlgError:
  169. # Floating point over/underflow, non-finite result from
  170. # matmul etc. -- report failure.
  171. return postprocess(x), k_outer + 1
  172. # Inner loop tolerance control
  173. if pres > ptol:
  174. ptol_max_factor = min(1.0, 1.5 * ptol_max_factor)
  175. else:
  176. ptol_max_factor = max(1e-16, 0.25 * ptol_max_factor)
  177. # -- GMRES terminated: eval solution
  178. dx = zs[0]*y[0]
  179. for w, yc in zip(zs[1:], y[1:]):
  180. dx = axpy(w, dx, dx.shape[0], yc) # dx += w*yc
  181. # -- Store LGMRES augmentation vectors
  182. nx = nrm2(dx)
  183. if nx > 0:
  184. if store_outer_Av:
  185. q = Q.dot(R.dot(y))
  186. ax = vs[0]*q[0]
  187. for v, qc in zip(vs[1:], q[1:]):
  188. ax = axpy(v, ax, ax.shape[0], qc)
  189. outer_v.append((dx/nx, ax/nx))
  190. else:
  191. outer_v.append((dx/nx, None))
  192. # -- Retain only a finite number of augmentation vectors
  193. while len(outer_v) > outer_k:
  194. del outer_v[0]
  195. # -- Apply step
  196. x += dx
  197. else:
  198. # didn't converge ...
  199. return postprocess(x), maxiter
  200. return postprocess(x), 0