hifir4py.ksp.gmres_hif

gmres_hif(A, b, M: Optional[hifir4py.hif.HIF] = None, **kw) Tuple[numpy.ndarray, int, dict][source]

HIF-preconditioned (right) restarted GMRES

This function implements the right-preconditioned GMRES(restart) with HIF preconditioners for consistent systems. This function can be either called as a “one-time” solver in that one can factorize HIF and compute the solution at the same time, or can be called repeatedly where the user can maintain the preconditioner and workspace buffer outside of this function. The implementation follows the standard right-preconditioned GMRES with restart; see Iterative Methods for Sparse Linear Systems by Saad.

Parameters:
  • A (csr_matrix) – Input CRS matrix
  • b (ndarray) – Input RHS vector (1D)
  • M (HIF, optional) – HIF preconditioner
  • x0 (ndarray, optional) – Initial guess, default is zeros; if passed in, this will be overwritten by the solution.
  • verbose ({0,1,2}, optional) – Verbose level, 0 disables the verbose printing, while 1 enables HIF verbose logging; default is 1.
  • restart (int, optional) – Restart dimension, default is 30
  • rtol (float, optional) – Relative tolerance threshold, default is 1e-6. The solver will terminate if \(\Vert b-Ax\Vert\le\text{rtol}\Vert b\Vert\).
  • work (GMRES_WorkSpace, optional) – Workspace buffer, by default, the solver creates all buffer space.
Returns:

  • x (ndarray) – The computed solution, which overwrites x0 if passed in.
  • flag (int) – Solver termination flags: The solver either successfully finishes (0), reaches the maximum iteration bound (1), breaks down due to NaN/Inf residuals (2), or stagnates (3).
  • info (dict) – Statistics of the solver, which are stored in a dictionary with fields: total iterations (key "iters"), final relative residual (key "relres"), relative residual history array (key "resids"), and factorization and solve runtimes (key "times").

Notes

For repeated calls of GMRES solver, one should at least maintain the buffer outside of this function, so that it won’t waste time in memory allocation and deallocation. In addition, one can maintain the preconditioner object M and update it whenever necessary; this is useful in solving dynamic and/or transient problems.

Notice that one can pass the key-value pairs **kw in HIF.factorize to this function call if he/she plans to let this solver build preconditioner.

Examples

There are many ways one can use this GMRES solver.

>>> import numpy as np
>>> from scipy.sparse import rand
>>> from hifir4py.ksp import *
>>> A = rand(10, 10, 0.5)
>>> b = A.dot(np.ones(10))

The simplest way is

>>> x, flag, _ = gmres_hif(A, b)

which is equivalent to the following call

>>> x, flag, _ = gmres_hif(A, b, M=None, verbose=1, restart=30, maxit=500,
... rtol=1e-6, x0=None, work=None)

In addition, if M is None, then one can pass any key-value pairs to this function that serve as control parameters of the preconditioner.

The following illustrates how to maintain the preconditioner and workspace buffer outside of the solver.

>>> from hifir4py import *
>>> M = HIF(A)
>>> work = GMRES_WorkSpace(b.size)
>>> x, flag, info = gmres_hif(A, b, M=M, work=work)
>>> # One can then update M and repeated call gmres_hif as above
>>> print("total #iters = {}".format(info["iters"]))