Fork me on GitHub

7. The Transpose-Free Quasi-Minimal Residual Method

7.1. The tfqmr Module

class tfqmr.TFQMR(op, **kwargs)

Bases: pykrylov.generic.generic.KrylovMethod

A pure Python implementation of the transpose-free quasi-minimal residual (TFQMR) algorithm. TFQMR may be used to solve unsymmetric systems of linear equations, i.e., systems of the form

A x = b

where the operator A may be unsymmetric.

TFQMR requires 2 operator-vector products with A, 4 dot products and 10 daxpys per iteration. It does not require products with the adjoint of A.

If a preconditioner is supplied, TFQMR needs to solve 2 preconditioning systems per iteration. Our implementation is inspired by the original description in [Freund] and that of [Kelley].

References:

[Freund]R. W. Freund, A Transpose-Free Quasi-Minimal Residual Method for Non-Hermitian Linear Systems, SIAM Journal on Scientific Computing, 14 (2), pp. 470–482, 1993.
solve(rhs, **kwargs)

Solve a linear system with rhs as right-hand side by the TFQMR method. The vector rhs should be a Numpy array.

Keywords:
guess:Initial guess (Numpy array, default: 0)
matvec_max:Max. number of matrix-vector produts (2n)

7.2. Example

Here is an example using TFQMR on a linear system. The coefficient matrix is read from file in Matrix Market format:

import numpy as np
from pykrylov.tfqmr import TFQMR as KSolver
from pykrylov.linop import PysparseLinearOperator
from pysparse import spmatrix
from pysparse.pysparseMatrix import PysparseMatrix as sp

A = sp(matrix=spmatrix.ll_mat_from_mtx('jpwh_991.mtx'))
n = A.shape[0]
e = np.ones(n)
rhs = A*e

ks = KSolver( PysparseLinearOperator(A),
              matvec_max=2*n,
              verbose=False,
              reltol = 1.0e-5 )
ks.solve(rhs, guess = 1+np.arange(n, dtype=np.float))

print 'Number of matvecs: ', ks.nMatvec
print 'Initial/final res: %8.2e/%8.2e' % (ks.residNorm0, ks.residNorm)
print 'Direct error: %8.2e' % (np.linalg.norm(ks.bestSolution-e)/sqrt(n))

Running this script produces the following output:

Number of matvecs:  70
Initial/final res: 8.64e+03/6.23e-04
Direct error: 2.77e-05