Chebyshev iteration

Chebyshev iteration solves the problem $Ax=b$ approximately for $x$ where $A$ is a symmetric, definite linear operator and $b$ the right-hand side vector. The methods assumes the interval $[\lambda_{min}, \lambda_{max}]$ containing all eigenvalues of $A$ is known, so that $x$ can be iteratively constructed via a Chebyshev polynomial with zeros in this interval. This polynomial ultimately acts as a filter that removes components in the direction of the eigenvectors from the initial residual.

The main advantage with respect to Conjugate Gradients is that BLAS1 operations such as inner products are avoided.

Usage

IterativeSolvers.chebyshev!Function
chebyshev!(x, A, b, λmin::Real, λmax::Real; kwargs...) -> x, [history]

Solve Ax = b for symmetric, definite matrices A using Chebyshev iteration.

Arguments

  • x: initial guess, will be updated in-place;
  • A: linear operator;
  • b: right-hand side;
  • λmin::Real: lower bound for the real eigenvalues
  • λmax::Real: upper bound for the real eigenvalues

Keywords

  • initially_zero::Bool = false: if true assumes that iszero(x) so that one matrix-vector product can be saved when computing the initial residual vector;
  • abstol::Real = zero(real(eltype(b))), reltol::Real = sqrt(eps(real(eltype(b)))): absolute and relative tolerance for the stopping condition |r_k| ≤ max(reltol * |r_0|, abstol), where r_k = A * x_k - b is the residual in the kth iteration;
  • maxiter::Int = size(A, 2): maximum number of inner iterations of GMRES;
  • Pl = Identity(): left preconditioner;
  • log::Bool = false: keep track of the residual norm in each iteration;
  • verbose::Bool = false: print convergence information during the iterations.

Return values

if log is false

  • x: approximate solution.

if log is true

  • x: approximate solution;
  • history: convergence history.
source

Implementation details

BLAS1 operations

Although the method is often used to avoid computation of inner products, the stopping criterion is still based on the residual norm. Hence the current implementation is not free of BLAS1 operations.

Tip

Chebyshev iteration can be used as an iterator.