def gcrotmk(A, b, x0=None, tol=1e-5, maxiter=1000, M=None, callback=None, m=20, k=None, CU=None, discard_C=False, truncate='oldest', atol=None): """ Solve a matrix equation using flexible GCROT(m,k) algorithm. Parameters ---------- A : {sparse matrix, ndarray, LinearOperator} The real or complex N-by-N matrix of the linear system. Alternatively, ``A`` can be a linear operator which can produce ``Ax`` using, e.g., ``scipy.sparse.linalg.LinearOperator``. b : ndarray Right hand side of the linear system. Has shape (N,) or (N,1). x0 : ndarray Starting guess for the solution. tol, atol : float, optional Tolerances for convergence, ``norm(residual) <= max(tol*norm(b), atol)``. The default for ``atol`` is `tol`. .. warning:: The default value for `atol` will be changed in a future release. For future compatibility, specify `atol` explicitly. maxiter : int, optional Maximum number of iterations. Iteration will stop after maxiter steps even if the specified tolerance has not been achieved. M : {sparse matrix, ndarray, LinearOperator}, optional Preconditioner for A. The preconditioner should approximate the inverse of A. gcrotmk is a 'flexible' algorithm and the preconditioner can vary from iteration to iteration. Effective preconditioning dramatically improves the rate of convergence, which implies that fewer iterations are needed to reach a given error tolerance. callback : function, optional User-supplied function to call after each iteration. It is called as callback(xk), where xk is the current solution vector. m : int, optional Number of inner FGMRES iterations per each outer iteration. Default: 20 k : int, optional Number of vectors to carry between inner FGMRES iterations. According to [2]_, good values are around m. Default: m CU : list of tuples, optional List of tuples ``(c, u)`` which contain the columns of the matrices C and U in the GCROT(m,k) algorithm. For details, see [2]_. The list given and vectors contained in it are modified in-place. If not given, start from empty matrices. The ``c`` elements in the tuples can be ``None``, in which case the vectors are recomputed via ``c = A u`` on start and orthogonalized as described in [3]_. discard_C : bool, optional Discard the C-vectors at the end. Useful if recycling Krylov subspaces for different linear systems. truncate : {'oldest', 'smallest'}, optional Truncation scheme to use. Drop: oldest vectors, or vectors with smallest singular values using the scheme discussed in [1,2]. See [2]_ for detailed comparison. Default: 'oldest' Returns ------- x : ndarray The solution found. info : int Provides convergence information: * 0 : successful exit * >0 : convergence to tolerance not achieved, number of iterations References ---------- .. [1] E. de Sturler, ''Truncation strategies for optimal Krylov subspace methods'', SIAM J. Numer. Anal. 36, 864 (1999). .. [2] J.E. Hicken and D.W. Zingg, ''A simplified and flexible variant of GCROT for solving nonsymmetric linear systems'', SIAM J. Sci. Comput. 32, 172 (2010). .. [3] M.L. Parks, E. de Sturler, G. Mackey, D.D. Johnson, S. Maiti, ''Recycling Krylov subspaces for sequences of linear systems'', SIAM J. Sci. Comput. 28, 1651 (2006). """ A, M, x, b, postprocess = make_system(A, M, x0, b) if not np.isfinite(b).all(): raise ValueError("RHS must contain only finite numbers") if truncate not in ('oldest', 'smallest'): raise ValueError("Invalid value for 'truncate': %r" % (truncate, )) if atol is None: warnings.warn( "scipy.sparse.linalg.gcrotmk called without specifying `atol`. " "The default value will change in the future. To preserve " "current behavior, set ``atol=tol``.", category=DeprecationWarning, stacklevel=2) atol = tol matvec = A.matvec psolve = M.matvec if CU is None: CU = [] if k is None: k = m axpy, dot, scal = None, None, None r = b - matvec(x) axpy, dot, scal, nrm2 = get_blas_funcs(['axpy', 'dot', 'scal', 'nrm2'], (x, r)) b_norm = nrm2(b) if b_norm == 0: x = b return (postprocess(x), 0) if discard_C: CU[:] = [(None, u) for c, u in CU] # Reorthogonalize old vectors if CU: # Sort already existing vectors to the front CU.sort(key=lambda cu: cu[0] is not None) # Fill-in missing ones C = np.empty((A.shape[0], len(CU)), dtype=r.dtype, order='F') us = [] j = 0 while CU: # More memory-efficient: throw away old vectors as we go c, u = CU.pop(0) if c is None: c = matvec(u) C[:, j] = c j += 1 us.append(u) # Orthogonalize Q, R, P = qr(C, overwrite_a=True, mode='economic', pivoting=True) del C # C := Q cs = list(Q.T) # U := U P R^-1, back-substitution new_us = [] for j in range(len(cs)): u = us[P[j]] for i in range(j): u = axpy(us[P[i]], u, u.shape[0], -R[i, j]) if abs(R[j, j]) < 1e-12 * abs(R[0, 0]): # discard rest of the vectors break u = scal(1.0 / R[j, j], u) new_us.append(u) # Form the new CU lists CU[:] = list(zip(cs, new_us))[::-1] if CU: axpy, dot = get_blas_funcs(['axpy', 'dot'], (r, )) # Solve first the projection operation with respect to the CU # vectors. This corresponds to modifying the initial guess to # be # # x' = x + U y # y = argmin_y || b - A (x + U y) ||^2 # # The solution is y = C^H (b - A x) for c, u in CU: yc = dot(c, r) x = axpy(u, x, x.shape[0], yc) r = axpy(c, r, r.shape[0], -yc) # GCROT main iteration for j_outer in range(maxiter): # -- callback if callback is not None: callback(x) beta = nrm2(r) # -- check stopping condition beta_tol = max(atol, tol * b_norm) if beta <= beta_tol and (j_outer > 0 or CU): # recompute residual to avoid rounding error r = b - matvec(x) beta = nrm2(r) if beta <= beta_tol: j_outer = -1 break ml = m + max(k - len(CU), 0) cs = [c for c, u in CU] try: Q, R, B, vs, zs, y, pres = _fgmres(matvec, r / beta, ml, rpsolve=psolve, atol=max(atol, tol * b_norm) / beta, cs=cs) y *= beta except LinAlgError: # Floating point over/underflow, non-finite result from # matmul etc. -- report failure. break # # At this point, # # [A U, A Z] = [C, V] G; G = [ I B ] # [ 0 H ] # # where [C, V] has orthonormal columns, and r = beta v_0. Moreover, # # || b - A (x + Z y + U q) ||_2 = || r - C B y - V H y - C q ||_2 = min! # # from which y = argmin_y || beta e_1 - H y ||_2, and q = -B y # # # GCROT(m,k) update # # Define new outer vectors # ux := (Z - U B) y ux = zs[0] * y[0] for z, yc in zip(zs[1:], y[1:]): ux = axpy(z, ux, ux.shape[0], yc) # ux += z*yc by = B.dot(y) for cu, byc in zip(CU, by): c, u = cu ux = axpy(u, ux, ux.shape[0], -byc) # ux -= u*byc # cx := V H y hy = Q.dot(R.dot(y)) cx = vs[0] * hy[0] for v, hyc in zip(vs[1:], hy[1:]): cx = axpy(v, cx, cx.shape[0], hyc) # cx += v*hyc # Normalize cx, maintaining cx = A ux # This new cx is orthogonal to the previous C, by construction try: alpha = 1 / nrm2(cx) if not np.isfinite(alpha): raise FloatingPointError() except (FloatingPointError, ZeroDivisionError): # Cannot update, so skip it continue cx = scal(alpha, cx) ux = scal(alpha, ux) # Update residual and solution gamma = dot(cx, r) r = axpy(cx, r, r.shape[0], -gamma) # r -= gamma*cx x = axpy(ux, x, x.shape[0], gamma) # x += gamma*ux # Truncate CU if truncate == 'oldest': while len(CU) >= k and CU: del CU[0] elif truncate == 'smallest': if len(CU) >= k and CU: # cf. [1,2] D = solve(R[:-1, :].T, B.T).T W, sigma, V = svd(D) # C := C W[:,:k-1], U := U W[:,:k-1] new_CU = [] for j, w in enumerate(W[:, :k - 1].T): c, u = CU[0] c = c * w[0] u = u * w[0] for cup, wp in zip(CU[1:], w[1:]): cp, up = cup c = axpy(cp, c, c.shape[0], wp) u = axpy(up, u, u.shape[0], wp) # Reorthogonalize at the same time; not necessary # in exact arithmetic, but floating point error # tends to accumulate here for cp, up in new_CU: alpha = dot(cp, c) c = axpy(cp, c, c.shape[0], -alpha) u = axpy(up, u, u.shape[0], -alpha) alpha = nrm2(c) c = scal(1.0 / alpha, c) u = scal(1.0 / alpha, u) new_CU.append((c, u)) CU[:] = new_CU # Add new vector to CU CU.append((cx, ux)) # Include the solution vector to the span CU.append((None, x.copy())) if discard_C: CU[:] = [(None, uz) for cz, uz in CU] return postprocess(x), j_outer + 1
def bicgstab(A, b, x0=None, tol=1e-5, criteria='rr', maxiter=None, M=None, callback=None, residuals=None): """Biconjugate Gradient Algorithm with Stabilization. Solves the linear system Ax = b. Left preconditioning is supported. Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria criteria : string Stopping criteria, let r=r_k, x=x_k 'rr': ||r|| < tol ||b|| 'rr+': ||r|| < tol (||b|| + ||A||_F ||x||) if ||b||=0, then set ||b||=1 for these tests. maxiter : int maximum number of iterations allowed M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list residual history in the 2-norm, including the initial residual Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. Examples -------- >>> from pyamg.krylov import bicgstab >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = bicgstab(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 4.68163 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 231-234, 2003 http://www-users.cs.umn.edu/~saad/books.html """ # Convert inputs to linear system, with error checking A, M, x, b, postprocess = make_system(A, M, x0, b) # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._bicgstab') # Check iteration numbers if maxiter is None: maxiter = len(x) + 5 elif maxiter < 1: raise ValueError('Number of iterations must be positive') # Prep for method r = b - A @ x normr = norm(r) if residuals is not None: residuals[:] = [normr] # Check initial guess, if b != 0, normb = norm(b) normb = norm(b) if normb == 0.0: normb = 1.0 # reset so that tol is unscaled # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': if sparse.issparse(A.A): normA = norm(A.A.data) elif isinstance(A.A, np.ndarray): normA = norm(np.ravel(A.A)) else: raise ValueError( 'Unable to use ||A||_F with the current matrix format.') rtol = tol * (normA * np.linalg.norm(x) + normb) else: raise ValueError('Invalid stopping criteria.') if normr < rtol: return (postprocess(x), 0) # Is thisAa one dimensional matrix? # Use a matvec to access A[0,0] if A.shape[0] == 1: entry = np.ravel(A @ np.array([1.0], dtype=x.dtype)) return (postprocess(b / entry), 0) rstar = r.copy() p = r.copy() rrstarOld = np.inner(rstar.conjugate(), r) it = 0 # Begin BiCGStab while True: Mp = M @ p AMp = A @ Mp # alpha = (r_j, rstar) / (A*M*p_j, rstar) alpha = rrstarOld / np.inner(rstar.conjugate(), AMp) # s_j = r_j - alpha*A*M*p_j s = r - alpha * AMp Ms = M @ s AMs = A @ Ms # omega = (A*M*s_j, s_j)/(A*M*s_j, A*M*s_j) omega = np.inner(AMs.conjugate(), s) / np.inner(AMs.conjugate(), AMs) # x_{j+1} = x_j + alpha*M*p_j + omega*M*s_j x = x + alpha * Mp + omega * Ms # r_{j+1} = s_j - omega*A*M*s r = s - omega * AMs # beta_j = (r_{j+1}, rstar)/(r_j, rstar) * (alpha/omega) rrstarNew = np.inner(rstar.conjugate(), r) beta = (rrstarNew / rrstarOld) * (alpha / omega) rrstarOld = rrstarNew # p_{j+1} = r_{j+1} + beta*(p_j - omega*A*M*p) p = r + beta * (p - omega * AMp) it += 1 normr = norm(r) if residuals is not None: residuals.append(normr) if callback is not None: callback(x) # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': rtol = tol * (normA * np.linalg.norm(x) + normb) if normr < rtol: return (postprocess(x), 0) if it == maxiter: return (postprocess(x), it)
def gmres_mgs(A, b, x0=None, tol=1e-5, restrt=None, maxiter=None, M=None, callback=None, residuals=None, reorth=False): """Generalized Minimum Residual Method (GMRES) based on MGS. GMRES iteratively refines the initial solution guess to the system Ax = b. Modified Gram-Schmidt version. Left preconditioning, leading to preconditioned residuals. Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria, let r=r_k ||M r|| < tol ||M b|| if ||b||=0, then set ||M b||=1 for these tests. restrt : None, int - if int, restrt is max number of inner iterations and maxiter is the max number of outer iterations - if None, do not restart GMRES, and max number of inner iterations is maxiter maxiter : None, int - if restrt is None, maxiter is the max number of inner iterations and GMRES does not restart - if restrt is int, maxiter is the max number of outer iterations, and restrt is the max number of inner iterations - defaults to min(n,40) if restart=None M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list preconditioned residual history in the 2-norm, including the initial preconditioned residual reorth : boolean If True, then a check is made whether to re-orthogonalize the Krylov space each GMRES iteration Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. For robustness, modified Gram-Schmidt is used to orthogonalize the Krylov Space Givens Rotations are used to provide the residual norm each iteration The residual is the *preconditioned* residual. Examples -------- >>> from pyamg.krylov import gmres >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = gmres(A,b, maxiter=2, tol=1e-8, orthog='mgs') >>> print(f'{norm(b - A*x):.6}') 6.54282 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 151-172, pp. 272-275, 2003 http://www-users.cs.umn.edu/~saad/books.html .. [2] C. T. Kelley, http://www4.ncsu.edu/~ctk/matlab_roots.html """ # Convert inputs to linear system, with error checking A, M, x, b, postprocess = make_system(A, M, x0, b) n = A.shape[0] # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._gmres_mgs') # Get fast access to underlying BLAS routines # dotc is the conjugate dot, dotu does no conjugation [lartg] = get_lapack_funcs(['lartg'], [x]) if np.iscomplexobj(np.zeros((1, ), dtype=x.dtype)): [axpy, dotu, dotc, scal] =\ get_blas_funcs(['axpy', 'dotu', 'dotc', 'scal'], [x]) else: # real type [axpy, dotu, dotc, scal] =\ get_blas_funcs(['axpy', 'dot', 'dot', 'scal'], [x]) # Set number of outer and inner iterations # If no restarts, # then set max_inner=maxiter and max_outer=n # If restarts are set, # then set max_inner=restart and max_outer=maxiter if restrt: if maxiter: max_outer = maxiter else: max_outer = 1 if restrt > n: warn('Setting restrt to maximum allowed, n.') restrt = n max_inner = restrt else: max_outer = 1 if maxiter > n: warn('Setting maxiter to maximum allowed, n.') maxiter = n elif maxiter is None: maxiter = min(n, 40) max_inner = maxiter # Is this a one dimensional matrix? if n == 1: entry = np.ravel(A @ np.array([1.0], dtype=x.dtype)) return (postprocess(b / entry), 0) # Prep for method r = b - A @ x # Apply preconditioner r = M @ r normr = norm(r) if residuals is not None: residuals[:] = [normr] # initial residual # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normMb = 1.0 # reset so that tol is unscaled else: normMb = norm(M @ b) # set the stopping criteria (see the docstring) if normr < tol * normMb: return (postprocess(x), 0) # Use separate variable to track iterations. If convergence fails, we # cannot simply report niter = (outer-1)*max_outer + inner. Numerical # error could cause the inner loop to halt while the actual ||r|| > tolerance. niter = 0 # Begin GMRES for _outer in range(max_outer): # Preallocate for Givens Rotations, Hessenberg matrix and Krylov Space # Space required is O(n*max_inner). # NOTE: We are dealing with row-major matrices, so we traverse in a # row-major fashion, # i.e., H and V's transpose is what we store. Q = [] # Givens Rotations # Upper Hessenberg matrix, which is then # converted to upper tri with Givens Rots H = np.zeros((max_inner + 1, max_inner + 1), dtype=x.dtype) V = np.zeros((max_inner + 1, n), dtype=x.dtype) # Krylov Space # vs store the pointers to each column of V. # This saves a considerable amount of time. vs = [] # v = r/normr V[0, :] = scal(1.0 / normr, r) vs.append(V[0, :]) # This is the RHS vector for the problem in the Krylov Space g = np.zeros((n, ), dtype=x.dtype) g[0] = normr for inner in range(max_inner): # New Search Direction v = V[inner + 1, :] v[:] = np.ravel(M @ (A @ vs[-1])) vs.append(v) normv_old = norm(v) # Modified Gram Schmidt for k in range(inner + 1): vk = vs[k] alpha = dotc(vk, v) H[inner, k] = alpha v[:] = axpy(vk, v, n, -alpha) normv = norm(v) H[inner, inner + 1] = normv # Re-orthogonalize if (reorth is True) and (normv_old == normv_old + 0.001 * normv): for k in range(inner + 1): vk = vs[k] alpha = dotc(vk, v) H[inner, k] = H[inner, k] + alpha v[:] = axpy(vk, v, n, -alpha) # Check for breakdown if H[inner, inner + 1] != 0.0: v[:] = scal(1.0 / H[inner, inner + 1], v) # Apply previous Givens rotations to H if inner > 0: apply_givens(Q, H[inner, :], inner) # Calculate and apply next complex-valued Givens Rotation # for the last inner iteration, when inner = n-1. # ==> Note that if max_inner = n, then this is unnecessary if inner != n - 1: if H[inner, inner + 1] != 0: [c, s, r] = lartg(H[inner, inner], H[inner, inner + 1]) Qblock = np.array([[c, s], [-np.conjugate(s), c]], dtype=x.dtype) Q.append(Qblock) # Apply Givens Rotation to g, # the RHS for the linear system in the Krylov Subspace. g[inner:inner + 2] = np.dot(Qblock, g[inner:inner + 2]) # Apply effect of Givens Rotation to H H[inner, inner] = dotu(Qblock[0, :], H[inner, inner:inner + 2]) H[inner, inner + 1] = 0.0 niter += 1 # Do not update normr if last inner iteration, because # normr is calculated directly after this loop ends. if inner < max_inner - 1: normr = np.abs(g[inner + 1]) if normr < tol * normMb: break if residuals is not None: residuals.append(normr) if callback is not None: y = sp.linalg.solve(H[0:inner + 1, 0:inner + 1].T, g[0:inner + 1]) update = np.ravel(V[:inner + 1, :].T.dot(y.reshape(-1, 1))) callback(x + update) # end inner loop, back to outer loop # Find best update to x in Krylov Space V. Solve inner x inner system. y = sp.linalg.solve(H[0:inner + 1, 0:inner + 1].T, g[0:inner + 1]) update = np.ravel(V[:inner + 1, :].T.dot(y.reshape(-1, 1))) x = x + update r = b - A @ x # Apply preconditioner r = M @ r normr = norm(r) # Allow user access to the iterates if callback is not None: callback(x) if residuals is not None: residuals.append(normr) # Has GMRES stagnated? indices = (x != 0) if indices.any(): change = np.max(np.abs(update[indices] / x[indices])) if change < 1e-12: # No change, halt return (postprocess(x), -1) # test for convergence if normr < tol * normMb: return (postprocess(x), 0) # end outer loop return (postprocess(x), niter)
def minimal_residual(A, b, x0=None, tol=1e-5, maxiter=None, M=None, callback=None, residuals=None): """Minimal residual (MR) algorithm. 1D projection method. Solves the linear system Ax = b. Left preconditioning is supported. Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria, let r=r_k ||M r|| < tol ||M b|| if ||b||=0, then set ||M b||=1 for these tests. maxiter : int maximum number of iterations allowed M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list preconditioned residual history in the 2-norm, including the initial preconditioned residual Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. minimal residual algorithm: Preconditioned version: r = b - A x r = b - A x, z = M r while not converged: while not converged: p = A r p = M A z alpha = (p,r) / (p,p) alpha = (p, z) / (p, p) x = x + alpha r x = x + alpha z r = r - alpha p z = z - alpha p See Also -------- _steepest_descent Examples -------- >>> from pyamg.krylov import minimal_residual >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = minimal_residual(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 7.26369 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 137--142, 2003 http://www-users.cs.umn.edu/~saad/books.html """ A, M, x, b, postprocess = make_system(A, M, x0, b) # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._minimal_residual') # determine maxiter if maxiter is None: maxiter = int(1.3 * len(b)) + 2 elif maxiter < 1: raise ValueError('Number of iterations must be positive') # setup method r = b - A @ x z = M @ r normr = norm(z) # store initial residual if residuals is not None: residuals[:] = [normr] # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normMb = 1.0 # reset so that tol is unscaled else: normMb = norm(M @ b) # set the stopping criteria (see the docstring) if normr < tol * normMb: return (postprocess(x), 0) # How often should r be recomputed recompute_r = 50 it = 0 while True: p = M @ (A @ z) # (p, z) = (M A M r, M r) = (M A z, z) pz = np.inner(p.conjugate(), z) # check curvature of M^-1 A if pz < 0.0: warn( '\nIndefinite matrix detected in minimal residual, stopping.\n' ) return (postprocess(x), -1) alpha = pz / np.inner(p.conjugate(), p) x = x + alpha * z it += 1 if np.mod(it, recompute_r) and it > 0: r = b - A @ x z = M @ r else: z = z - alpha * p normr = norm(z) if residuals is not None: residuals.append(normr) if callback is not None: callback(x) # set the stopping criteria (see the docstring) if normr < tol * normMb: return (postprocess(x), 0) if it == maxiter: return (postprocess(x), it)
def fgmres(A, b, x0=None, tol=1e-5, restrt=None, maxiter=None, M=None, callback=None, residuals=None): """Flexible Generalized Minimum Residual Method (fGMRES). fGMRES iteratively refines the initial solution guess to the system Ax = b. fGMRES is flexible in the sense that the right preconditioner (M) can vary from iteration to iteration. Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria, let r=r_k ||r|| < tol ||b|| if ||b||=0, then set ||b||=1 for these tests. restrt : None, int - if int, restrt is max number of inner iterations and maxiter is the max number of outer iterations - if None, do not restart GMRES, and max number of inner iterations is maxiter maxiter : None, int - if restrt is None, maxiter is the max number of inner iterations and GMRES does not restart - if restrt is int, maxiter is the max number of outer iterations, and restrt is the max number of inner iterations - defaults to min(n,40) if restart=None M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. M need not be stationary for fgmres callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list residual history in the 2-norm, including the initial residual reorth : boolean If True, then a check is made whether to re-orthogonalize the Krylov space each GMRES iteration Returns ------- xk, info xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. fGMRES allows for non-stationary preconditioners, as opposed to GMRES For robustness, Householder reflections are used to orthonormalize the Krylov Space Givens Rotations are used to provide the residual norm each iteration Flexibility implies that the right preconditioner, M, can vary from iteration to iteration Examples -------- >>> from pyamg.krylov import fgmres >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = fgmres(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 6.54282 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 151-172, pp. 272-275, 2003 http://www-users.cs.umn.edu/~saad/books.html """ # Convert inputs to linear system, with error checking A, M, x, b, postprocess = make_system(A, M, x0, b) n = A.shape[0] # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._fgmres') # Get fast access to underlying BLAS routines [lartg] = get_lapack_funcs(['lartg'], [x]) # Set number of outer and inner iterations # If no restarts, # then set max_inner=maxiter and max_outer=n # If restarts are set, # then set max_inner=restart and max_outer=maxiter if restrt: if maxiter: max_outer = maxiter else: max_outer = 1 if restrt > n: warn('Setting restrt to maximum allowed, n.') restrt = n max_inner = restrt else: max_outer = 1 if maxiter > n: warn('Setting maxiter to maximum allowed, n.') maxiter = n elif maxiter is None: maxiter = min(n, 40) max_inner = maxiter # Is this a one dimensional matrix? if n == 1: entry = np.ravel(A @ np.array([1.0], dtype=x.dtype)) return (postprocess(b / entry), 0) # Prep for method r = b - A @ x normr = norm(r) if residuals is not None: residuals[:] = [normr] # initial residual # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normb = 1.0 # reset so that tol is unscaled if normr < tol * normb: return (postprocess(x), 0) # Use separate variable to track iterations. If convergence fails, # we cannot simply report niter = (outer-1)*max_outer + inner. Numerical # error could cause the inner loop to halt while the actual ||r|| > tol. niter = 0 # Begin fGMRES for _outer in range(max_outer): # Calculate vector w, which defines the Householder reflector # Take shortcut in calculating, # w = r + sign(r[1])*||r||_2*e_1 w = r beta = _mysign(w[0]) * normr w[0] += beta w /= norm(w) # Preallocate for Krylov vectors, Householder reflectors and Hessenberg # matrix # Space required is O(n*max_inner) # Givens Rotations Q = np.zeros((4 * max_inner, ), dtype=x.dtype) # upper Hessenberg matrix (made upper tri with Givens Rotations) H = np.zeros((max_inner, max_inner), dtype=x.dtype) W = np.zeros((max_inner, n), dtype=x.dtype) # Householder reflectors # For fGMRES, preconditioned vectors must be stored # No Horner-like scheme exists that allow us to avoid this Z = np.zeros((n, max_inner), dtype=x.dtype) W[0, :] = w # Multiply r with (I - 2*w*w.T), i.e. apply the Householder reflector # This is the RHS vector for the problem in the Krylov Space g = np.zeros((n, ), dtype=x.dtype) g[0] = -beta for inner in range(max_inner): # Calculate Krylov vector in two steps # (1) Calculate v = P_j = (I - 2*w*w.T)v, where k = inner v = -2.0 * np.conjugate(w[inner]) * w v[inner] += 1.0 # (2) Calculate the rest, v = P_1*P_2*P_3...P_{j-1}*ej. # for j in range(inner-1,-1,-1): # v = v - 2.0*dot(conjugate(W[j,:]), v)*W[j,:] amg_core.apply_householders(v, np.ravel(W), n, inner - 1, -1, -1) # Apply preconditioner v = M @ v # Check for nan, inf # if isnan(v).any() or isinf(v).any(): # warn('inf or nan after application of preconditioner') # return(postprocess(x), -1) Z[:, inner] = v # Calculate new search direction v = A @ v # Factor in all Householder orthogonal reflections on new search # direction # for j in range(inner+1): # v = v - 2.0*dot(conjugate(W[j,:]), v)*W[j,:] amg_core.apply_householders(v, np.ravel(W), n, 0, inner + 1, 1) # Calculate next Householder reflector, w # w = v[inner+1:] + sign(v[inner+1])*||v[inner+1:]||_2*e_{inner+1) # Note that if max_inner = n, then this is unnecessary for # the last inner iteration, when inner = n-1. Here we do # not need to calculate a Householder reflector or Givens # rotation because nnz(v) is already the desired length, # i.e. we do not need to zero anything out. if inner != n - 1: if inner < (max_inner - 1): w = W[inner + 1, :] vslice = v[inner + 1:] alpha = norm(vslice) if alpha != 0: alpha = _mysign(vslice[0]) * alpha # do not need the final reflector for future calculations if inner < (max_inner - 1): w[inner + 1:] = vslice w[inner + 1] += alpha w /= norm(w) # Apply new reflector to v # v = v - 2.0*w*(w.T*v) v[inner + 1] = -alpha v[inner + 2:] = 0.0 if inner > 0: # Apply all previous Givens Rotations to v amg_core.apply_givens(Q, v, n, inner) # Calculate the next Givens rotation, where j = inner Note that if # max_inner = n, then this is unnecessary for the last inner # iteration, when inner = n-1. Here we do not need to # calculate a Householder reflector or Givens rotation because # nnz(v) is already the desired length, i.e. we do not need to zero # anything out. if inner != n - 1: if v[inner + 1] != 0: [c, s, r] = lartg(v[inner], v[inner + 1]) Qblock = np.array([[c, s], [-np.conjugate(s), c]], dtype=x.dtype) Q[(inner * 4):((inner + 1) * 4)] = np.ravel(Qblock).copy() # Apply Givens Rotation to g, the RHS for the linear system # in the Krylov Subspace. Note that this dot does a matrix # multiply, not an actual dot product where a conjugate # transpose is taken g[inner:inner + 2] = np.dot(Qblock, g[inner:inner + 2]) # Apply effect of Givens Rotation to v v[inner] = np.dot(Qblock[0, :], v[inner:inner + 2]) v[inner + 1] = 0.0 # Write to upper Hessenberg Matrix, # the LHS for the linear system in the Krylov Subspace H[:, inner] = v[0:max_inner] # Don't update normr if last inner iteration, because # normr is calculated directly after this loop ends. if inner < max_inner - 1: normr = np.abs(g[inner + 1]) if normr < tol * normb: break if residuals is not None: residuals.append(normr) if callback is not None: y = sp.linalg.solve(H[0:(inner + 1), 0:(inner + 1)], g[0:(inner + 1)]) update = np.dot(Z[:, 0:inner + 1], y) callback(x + update) niter += 1 # end inner loop, back to outer loop # Find best update to x in Krylov Space, V. Solve inner+1 x inner+1 # system. Apparently this is the best way to solve a triangular system # in the magical world of scipy # piv = arange(inner+1) # y = lu_solve((H[0:(inner+1),0:(inner+1)], piv), # g[0:(inner+1)], trans=0) y = sp.linalg.solve(H[0:(inner + 1), 0:(inner + 1)], g[0:(inner + 1)]) # No Horner like scheme exists because the preconditioner can change # each iteration # Hence, we must store each preconditioned vector update = np.dot(Z[:, 0:inner + 1], y) x = x + update r = b - A @ x # Allow user access to the iterates if callback is not None: callback(x) normr = norm(r) if residuals is not None: residuals.append(normr) # Has fGMRES stagnated? indices = (x != 0) if indices.any(): change = np.max(np.abs(update[indices] / x[indices])) if change < 1e-12: # No change, halt return (postprocess(x), -1) # test for convergence if normr < tol * normb: return (postprocess(x), 0) # end outer loop return (postprocess(x), niter)
def steepest_descent(A, b, x0=None, tol=1e-5, criteria='rr', maxiter=None, M=None, callback=None, residuals=None): """Steepest descent algorithm. Solves the linear system Ax = b. Left preconditioning is supported. Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria criteria : string Stopping criteria, let r=r_k, x=x_k 'rr': ||r|| < tol ||b|| 'rr+': ||r|| < tol (||b|| + ||A||_F ||x||) 'MrMr': ||M r|| < tol ||M b|| 'rMr': <r, Mr>^1/2 < tol if ||b||=0, then set ||b||=1 for these tests. maxiter : int maximum number of iterations allowed M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list residual history in the 2-norm, including the initial residual Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status of cg == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. See Also -------- _minimal_residual Examples -------- >>> from pyamg.krylov import steepest_descent >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = steepest_descent(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 7.89436 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 137--142, 2003 http://www-users.cs.umn.edu/~saad/books.html """ A, M, x, b, postprocess = make_system(A, M, x0, b) # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._steepest_descent') # determine maxiter if maxiter is None: maxiter = int(len(b)) elif maxiter < 1: raise ValueError('Number of iterations must be positive') # setup method r = b - A @ x z = M @ r rz = np.inner(r.conjugate(), z) normr = np.linalg.norm(r) if residuals is not None: residuals[:] = [normr] # initial residual # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normb = 1.0 # reset so that tol is unscaled # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': if sparse.issparse(A.A): normA = norm(A.A.data) elif isinstance(A.A, np.ndarray): normA = norm(np.ravel(A.A)) else: raise ValueError( 'Unable to use ||A||_F with the current matrix format.') rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = norm(z) normMb = norm(M @ b) rtol = tol * normMb elif criteria == 'rMr': normr = np.sqrt(rz) rtol = tol else: raise ValueError('Invalid stopping criteria.') # How often should r be recomputed recompute_r = 50 it = 0 while True: q = A @ z zAz = np.inner(z.conjugate(), q) # check curvature of A if zAz < 0.0: warn( '\nIndefinite matrix detected in steepest descent, aborting\n') return (postprocess(x), -1) alpha = rz / zAz # step size x = x + alpha * z it += 1 if np.mod(it, recompute_r) and it > 0: r = b - A @ x else: r = r - alpha * q z = M @ r rz = np.inner(r.conjugate(), z) if rz < 0.0: # check curvature of M warn( '\nIndefinite preconditioner detected in steepest descent, stopping.\n' ) return (postprocess(x), -1) normr = norm(r) if residuals is not None: residuals.append(normr) if callback is not None: callback(x) # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = norm(z) rtol = tol * normMb elif criteria == 'rMr': normr = np.sqrt(rz) rtol = tol if normr < rtol: return (postprocess(x), 0) if rz == 0.0: # important to test after testing normr < tol. rz == 0.0 is an # indicator of convergence when r = 0.0 warn( '\nSingular preconditioner detected in steepest descent, stopping.\n' ) return (postprocess(x), -1) if it == maxiter: return (postprocess(x), it)
def cgne(A, b, x0=None, tol=1e-5, criteria='rr', maxiter=None, M=None, callback=None, residuals=None): """Conjugate Gradient, Normal Error algorithm. Applies CG to the normal equations, A A.H x = b. Left preconditioning is supported. Note that unless A is well-conditioned, the use of CGNE is inadvisable Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria criteria : string Stopping criteria, let r=r_k, x=x_k 'rr': ||r|| < tol ||b|| 'rr+': ||r|| < tol (||b|| + ||A||_F ||x||) 'MrMr': ||M r|| < tol ||M b|| 'rMr': <r, Mr>^1/2 < tol if ||b||=0, then set ||b||=1 for these tests. maxiter : int maximum number of iterations allowed M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A A.H x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list residual history in the 2-norm, including the initial residual Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. Examples -------- >>> from pyamg.krylov import cgne >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = cgne(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 46.1547 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 276-7, 2003 http://www-users.cs.umn.edu/~saad/books.html """ # Store the conjugate transpose explicitly as it will be used much later on if sparse.isspmatrix(A): AH = A.H else: # avoid doing this since A may be a different sparse type AH = aslinearoperator(np.asarray(A).conj().T) # Convert inputs to linear system, with error checking A, M, x, b, postprocess = make_system(A, M, x0, b) n = A.shape[0] # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._cgne') # How often should r be recomputed recompute_r = 8 # Check iteration numbers. CGNE suffers from loss of orthogonality quite # easily, so we arbitrarily let the method go up to 130% over the # theoretically necessary limit of maxiter=n if maxiter is None: maxiter = int(np.ceil(1.3 * n)) + 2 elif maxiter < 1: raise ValueError('Number of iterations must be positive') elif maxiter > (1.3 * n): warn('maximum allowed inner iterations (maxiter) are the 130% times' 'the number of dofs') maxiter = int(np.ceil(1.3 * n)) + 2 # Prep for method r = b - A @ x normr = norm(r) # Apply preconditioner and calculate initial search direction z = M @ r p = AH @ z old_zr = np.inner(z.conjugate(), r) if residuals is not None: residuals[:] = [normr] # initial residual # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normb = 1.0 # reset so that tol is unscaled # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': if sparse.issparse(A.A): normA = norm(A.A.data) elif isinstance(A.A, np.ndarray): normA = norm(np.ravel(A.A)) else: raise ValueError( 'Unable to use ||A||_F with the current matrix format.') rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = norm(z) normMb = norm(M @ b) rtol = tol * normMb elif criteria == 'rMr': normr = np.sqrt(old_zr) rtol = tol else: raise ValueError('Invalid stopping criteria.') if normr < rtol: return (postprocess(x), 0) # Begin CGNE it = 0 while True: # Step number in Saad's pseudocode # alpha = (z_j, r_j) / (p_j, p_j) alpha = old_zr / np.inner(p.conjugate(), p) # x_{j+1} = x_j + alpha*p_j x += alpha * p # r_{j+1} = r_j - alpha*w_j, where w_j = A*p_j if np.mod(it, recompute_r) and it > 0: r -= alpha * (A @ p) else: r = b - A @ x # z_{j+1} = M*r_{j+1} z = M @ r # beta = (z_{j+1}, r_{j+1}) / (z_j, r_j) new_zr = np.inner(z.conjugate(), r) beta = new_zr / old_zr old_zr = new_zr # p_{j+1} = A.H*z_{j+1} + beta*p_j p *= beta p += AH @ z it += 1 normr = np.linalg.norm(r) if residuals is not None: residuals.append(normr) if callback is not None: callback(x) # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = norm(z) rtol = tol * normMb elif criteria == 'rMr': normr = np.sqrt(new_zr) rtol = tol if normr < rtol: return (postprocess(x), 0) if it == maxiter: return (postprocess(x), it)
def cr(A, b, x0=None, tol=1e-5, criteria='rr', maxiter=None, M=None, callback=None, residuals=None): """Conjugate Residual algorithm. Solves the linear system Ax = b. Left preconditioning is supported. The matrix A must be Hermitian symmetric (but not necessarily definite). Parameters ---------- A : array, matrix, sparse matrix, LinearOperator n x n, linear system to solve b : array, matrix right hand side, shape is (n,) or (n,1) x0 : array, matrix initial guess, default is a vector of zeros tol : float Tolerance for stopping criteria criteria : string Stopping criteria, let r=r_k, x=x_k 'rr': ||r|| < tol ||b|| 'rr+': ||r|| < tol (||b|| + ||A||_F ||x||) 'MrMr': ||M r|| < tol ||M b|| if ||b||=0, then set ||b||=1 for these tests. maxiter : int maximum number of iterations allowed M : array, matrix, sparse matrix, LinearOperator n x n, inverted preconditioner, i.e. solve M A x = M b. callback : function User-supplied function is called after each iteration as callback(xk), where xk is the current solution vector residuals : list residual history in the 2-norm, including the initial residual Returns ------- (xk, info) xk : an updated guess after k iterations to the solution of Ax = b info : halting status == ======================================= 0 successful exit >0 convergence to tolerance not achieved, return iteration count instead. <0 numerical breakdown, or illegal input == ======================================= Notes ----- The LinearOperator class is in scipy.sparse.linalg.interface. Use this class if you prefer to define A or M as a mat-vec routine as opposed to explicitly constructing the matrix. Examples -------- >>> from pyamg.krylov import cr >>> from pyamg.util.linalg import norm >>> import numpy as np >>> from pyamg.gallery import poisson >>> A = poisson((10,10)) >>> b = np.ones((A.shape[0],)) >>> (x,flag) = cr(A,b, maxiter=2, tol=1e-8) >>> print(f'{norm(b - A*x):.6}') 6.54282 References ---------- .. [1] Yousef Saad, "Iterative Methods for Sparse Linear Systems, Second Edition", SIAM, pp. 262-67, 2003 http://www-users.cs.umn.edu/~saad/books.html """ A, M, x, b, postprocess = make_system(A, M, x0, b) # Ensure that warnings are always reissued from this function warnings.filterwarnings('always', module='pyamg.krylov._cr') # determine maxiter if maxiter is None: maxiter = int(1.3*len(b)) + 2 elif maxiter < 1: raise ValueError('Number of iterations must be positive') # setup method r = b - A @ x z = M @ r p = z.copy() zz = np.inner(z.conjugate(), z) normr = np.linalg.norm(r) if residuals is not None: residuals[:] = [normr] # initial residual # Check initial guess if b != 0, normb = norm(b) if normb == 0.0: normb = 1.0 # reset so that tol is unscaled # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': if sparse.issparse(A.A): normA = norm(A.A.data) elif isinstance(A.A, np.ndarray): normA = norm(np.ravel(A.A)) else: raise ValueError('Unable to use ||A||_F with the current matrix format.') rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = np.sqrt(zz) normMb = norm(M @ b) rtol = tol * normMb else: raise ValueError('Invalid stopping criteria.') if normr < rtol: return (postprocess(x), 0) # How often should r be recomputed recompute_r = 8 Az = A @ z rAz = np.inner(r.conjugate(), Az) Ap = A @ p it = 0 while True: rAz_old = rAz alpha = rAz / np.inner(Ap.conjugate(), Ap) # 3 x += alpha * p # 4 if np.mod(it, recompute_r) and it > 0: # 5 r -= alpha * Ap else: r = b - A @ x z = M @ r Az = A @ z rAz = np.inner(r.conjugate(), Az) beta = rAz/rAz_old # 6 p *= beta # 7 p += z Ap *= beta # 8 Ap += Az it += 1 zz = np.inner(z.conjugate(), z) normr = np.linalg.norm(r) if residuals is not None: residuals.append(normr) if callback is not None: callback(x) # set the stopping criteria (see the docstring) if criteria == 'rr': rtol = tol * normb elif criteria == 'rr+': rtol = tol * (normA * np.linalg.norm(x) + normb) elif criteria == 'MrMr': normr = norm(z) rtol = tol * normMb if normr < rtol: return (postprocess(x), 0) if zz == 0.0: # rz == 0.0 is an indicator of convergence when r = 0.0 warn('\nSingular preconditioner detected in CR, ceasing iterations\n') return (postprocess(x), -1) if it == maxiter: return (postprocess(x), it)