KSPBCGSL#

Implements a slight variant of the Enhanced BiCGStab(L) algorithm in [Fok96] and [SVdVF94], see also [SVdV95]. The variation concerns cases when either kappa02 or kappa12 is negative due to round-off. Kappa0 has also been pulled out of the denominator in the formula for ghat.

Options Database Keys#

  • -ksp_bcgsl_ell - Number of Krylov search directions, defaults to 2, cf. KSPBCGSLSetEll()

  • -ksp_bcgsl_cxpol - Use a convex function of the MinRes and OR polynomials after the BiCG step instead of default MinRes, cf. KSPBCGSLSetPol()

  • -ksp_bcgsl_mrpoly - Use the default MinRes polynomial after the BiCG step, cf. KSPBCGSLSetPol()

  • -ksp_bcgsl_xres - Threshold used to decide when to refresh computed residuals, cf. KSPBCGSLSetXRes()

  • -ksp_bcgsl_pinv <true/false> - (de)activate use of pseudoinverse, cf. KSPBCGSLSetUsePseudoinverse()

Contributed by#

Joel M. Malard, email jm.malard@pnl.gov

Note#

The “sub-iterations” of this solver are not reported by -ksp_monitor or recorded in KSPSetResidualHistory() since the solution is not directly computed for these sub-iterations.

Developer Notes#

This has not been completely cleaned up into PETSc style.

All the BLAS and LAPACK calls in the source should be removed and replaced with loops and the macros for block solvers converted from LINPACK.

References#

Fok96

Diederik R Fokkema. Enhanced implementation of BiCGstab(l) for solving linear systems of equations. Citeseer, 1996.

SVdV95

Gerard LG Sleijpen and Henk A Van der Vorst. An overview of approaches for the stable computation of hybrid BiCG methods. Applied numerical mathematics, 19(3):235–254, 1995.

SVdVF94

Gerard LG Sleijpen, Henk A Van der Vorst, and Diederik R Fokkema. BiCGstab(l) and other hybrid Bi-CG methods. Numerical Algorithms, 7:75–109, 1994.

See Also#

KSP: Linear System Solvers, KSPFBCGS, KSPFBCGSR, KSPBCGS, KSPPIPEBCGS, KSPCreate(), KSPSetType(), KSPType, KSP, KSPFGMRES, KSPBCGS, KSPSetPCSide(), KSPBCGSLSetEll(), KSPBCGSLSetXRes(), KSPBCGSLSetUsePseudoinverse(), KSPBCGSLSetPol()

Level#

intermediate

Location#

src/ksp/ksp/impls/bcgsl/bcgsl.c


Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages