lgmres.py 8.6 KB

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