BiCGStab(l)

BiCGStab(l) solves the problem $Ax = b$ approximately for $x$ where $A$ is a general, linear operator and $b$ the right-hand side vector. The methods combines BiCG with $l$ GMRES iterations, resulting in a short-reccurence iteration. As a result the memory is fixed as well as the computational costs per iteration.

Usage

IterativeSolvers.bicgstabl!Function
bicgstabl!(x, A, b, l; kwargs...) -> x, [history]

Arguments

  • A: linear operator;
  • b: right hand side (vector);
  • l::Int = 2: Number of GMRES steps.

Keywords

  • max_mv_products::Int = size(A, 2): maximum number of matrix vector products.

For BiCGStab(l) this is a less dubious term than "number of iterations";

  • Pl = Identity(): left preconditioner of the method;
  • 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 approximate residual in the kth iteration;
    Note
    1. The true residual norm is never computed during the iterations, only an approximation;
    2. If a left preconditioner is given, the stopping condition is based on the preconditioned residual.

Return values

if log is false

  • x: approximate solution.

if log is true

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

Implementation details

The method is based on the original article [Sleijpen1993], but does not implement later improvements. The normal equations arising from the GMRES steps are solved without orthogonalization. Hence the method should only be reliable for relatively small values of $l$.

The r and u factors are pre-allocated as matrices of size $n \times (l + 1)$, so that BLAS2 methods can be used. Also the random shadow residual is pre-allocated as a vector. Hence the storage costs are approximately $2l + 3$ vectors.

Tip

BiCGStabl(l) can be used as an iterator.

  • Sleijpen1993Sleijpen, Gerard LG, and Diederik R. Fokkema. "BiCGstab(l) for linear equations involving unsymmetric matrices with complex spectrum." Electronic Transactions on Numerical Analysis 1.11 (1993): 2000.