tfqmr.py 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179
  1. import numpy as np
  2. from .iterative import _get_atol_rtol
  3. from .utils import make_system
  4. __all__ = ['tfqmr']
  5. def tfqmr(A, b, x0=None, *, rtol=1e-5, atol=0., maxiter=None, M=None,
  6. callback=None, show=False):
  7. """
  8. Solve ``Ax = b`` with the Transpose-Free Quasi-Minimal Residual method.
  9. Parameters
  10. ----------
  11. A : {sparse array, ndarray, LinearOperator}
  12. The real or complex N-by-N matrix of the linear system.
  13. Alternatively, `A` can be a linear operator which can
  14. produce ``Ax`` using, e.g.,
  15. `scipy.sparse.linalg.LinearOperator`.
  16. b : {ndarray}
  17. Right hand side of the linear system. Has shape (N,) or (N,1).
  18. x0 : {ndarray}
  19. Starting guess for the solution.
  20. rtol, atol : float, optional
  21. Parameters for the convergence test. For convergence,
  22. ``norm(b - A @ x) <= max(rtol*norm(b), atol)`` should be satisfied.
  23. The default is ``rtol=1e-5``, the default for ``atol`` is ``0.0``.
  24. maxiter : int, optional
  25. Maximum number of iterations. Iteration will stop after maxiter
  26. steps even if the specified tolerance has not been achieved.
  27. Default is ``min(10000, ndofs * 10)``, where ``ndofs = A.shape[0]``.
  28. M : {sparse array, ndarray, LinearOperator}
  29. Inverse of the preconditioner of A. M should approximate the
  30. inverse of A and be easy to solve for (see Notes). Effective
  31. preconditioning dramatically improves the rate of convergence,
  32. which implies that fewer iterations are needed to reach a given
  33. error tolerance. By default, no preconditioner is used.
  34. callback : function, optional
  35. User-supplied function to call after each iteration. It is called
  36. as ``callback(xk)``, where ``xk`` is the current solution vector.
  37. show : bool, optional
  38. Specify ``show = True`` to show the convergence, ``show = False`` is
  39. to close the output of the convergence.
  40. Default is `False`.
  41. Returns
  42. -------
  43. x : ndarray
  44. The converged solution.
  45. info : int
  46. Provides convergence information:
  47. - 0 : successful exit
  48. - >0 : convergence to tolerance not achieved, number of iterations
  49. - <0 : illegal input or breakdown
  50. Notes
  51. -----
  52. The Transpose-Free QMR algorithm is derived from the CGS algorithm.
  53. However, unlike CGS, the convergence curves for the TFQMR method is
  54. smoothed by computing a quasi minimization of the residual norm. The
  55. implementation supports left preconditioner, and the "residual norm"
  56. to compute in convergence criterion is actually an upper bound on the
  57. actual residual norm ``||b - Axk||``.
  58. References
  59. ----------
  60. .. [1] R. W. Freund, A Transpose-Free Quasi-Minimal Residual Algorithm for
  61. Non-Hermitian Linear Systems, SIAM J. Sci. Comput., 14(2), 470-482,
  62. 1993.
  63. .. [2] Y. Saad, Iterative Methods for Sparse Linear Systems, 2nd edition,
  64. SIAM, Philadelphia, 2003.
  65. .. [3] C. T. Kelley, Iterative Methods for Linear and Nonlinear Equations,
  66. number 16 in Frontiers in Applied Mathematics, SIAM, Philadelphia,
  67. 1995.
  68. Examples
  69. --------
  70. >>> import numpy as np
  71. >>> from scipy.sparse import csc_array
  72. >>> from scipy.sparse.linalg import tfqmr
  73. >>> A = csc_array([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=float)
  74. >>> b = np.array([2, 4, -1], dtype=float)
  75. >>> x, exitCode = tfqmr(A, b, atol=0.0)
  76. >>> print(exitCode) # 0 indicates successful convergence
  77. 0
  78. >>> np.allclose(A.dot(x), b)
  79. True
  80. """
  81. # Check data type
  82. dtype = A.dtype
  83. if np.issubdtype(dtype, np.int64):
  84. dtype = float
  85. A = A.astype(dtype)
  86. if np.issubdtype(b.dtype, np.int64):
  87. b = b.astype(dtype)
  88. A, M, x, b = make_system(A, M, x0, b)
  89. # Check if the R.H.S is a zero vector
  90. if np.linalg.norm(b) == 0.:
  91. x = b.copy()
  92. return (x, 0)
  93. ndofs = A.shape[0]
  94. if maxiter is None:
  95. maxiter = min(10000, ndofs * 10)
  96. if x0 is None:
  97. r = b.copy()
  98. else:
  99. r = b - A.matvec(x)
  100. u = r
  101. w = r.copy()
  102. # Take rstar as b - Ax0, that is rstar := r = b - Ax0 mathematically
  103. rstar = r
  104. v = M.matvec(A.matvec(r))
  105. uhat = v
  106. d = theta = eta = 0.
  107. # at this point we know rstar == r, so rho is always real
  108. rho = np.inner(rstar.conjugate(), r).real
  109. rhoLast = rho
  110. r0norm = np.sqrt(rho)
  111. tau = r0norm
  112. if r0norm == 0:
  113. return (x, 0)
  114. # we call this to get the right atol and raise errors as necessary
  115. atol, _ = _get_atol_rtol('tfqmr', r0norm, atol, rtol)
  116. for iter in range(maxiter):
  117. even = iter % 2 == 0
  118. if (even):
  119. vtrstar = np.inner(rstar.conjugate(), v)
  120. # Check breakdown
  121. if vtrstar == 0.:
  122. return (x, -1)
  123. alpha = rho / vtrstar
  124. uNext = u - alpha * v # [1]-(5.6)
  125. w -= alpha * uhat # [1]-(5.8)
  126. d = u + (theta**2 / alpha) * eta * d # [1]-(5.5)
  127. # [1]-(5.2)
  128. theta = np.linalg.norm(w) / tau
  129. c = np.sqrt(1. / (1 + theta**2))
  130. tau *= theta * c
  131. # Calculate step and direction [1]-(5.4)
  132. eta = (c**2) * alpha
  133. z = M.matvec(d)
  134. x += eta * z
  135. if callback is not None:
  136. callback(x)
  137. # Convergence criterion
  138. if tau * np.sqrt(iter+1) < atol:
  139. if (show):
  140. print("TFQMR: Linear solve converged due to reach TOL "
  141. f"iterations {iter+1}")
  142. return (x, 0)
  143. if (not even):
  144. # [1]-(5.7)
  145. rho = np.inner(rstar.conjugate(), w)
  146. beta = rho / rhoLast
  147. u = w + beta * u
  148. v = beta * uhat + (beta**2) * v
  149. uhat = M.matvec(A.matvec(u))
  150. v += uhat
  151. else:
  152. uhat = M.matvec(A.matvec(uNext))
  153. u = uNext
  154. rhoLast = rho
  155. if (show):
  156. print("TFQMR: Linear solve not converged due to reach MAXIT "
  157. f"iterations {iter+1}")
  158. return (x, maxiter)