KSPBCGSL#
Implements a slight variant of the Enhanced BiCGStab(L) algorithm in (3) and (2). 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.
References#
**** -*** G.L.G. Sleijpen, H.A. van der Vorst, “An overview of approaches for the stable computation of hybrid BiCG methods”, Applied Numerical Mathematics: Transactions f IMACS, 19(3), 1996.
**** -*** G.L.G. Sleijpen, H.A. van der Vorst, D.R. Fokkema, “BiCGStab(L) and other hybrid BiCG methods”, Numerical Algorithms, 7, 1994.
**** -*** D.R. Fokkema, “Enhanced implementation of BiCGStab(L) for solving linear systems of equations”, preprint from www.citeseer.com.
Contributed by: Joel M. Malard, email jm.malard@pnl.gov
Options Database Keys#
-ksp_bcgsl_ell
Number of Krylov search directions, defaults to 2 - - KSPBCGSLSetEll()-ksp_bcgsl_cxpol - Use a convex function of the MinRes and OR polynomials after the BiCG step instead of default MinRes – KSPBCGSLSetPol()
-ksp_bcgsl_mrpoly - Use the default MinRes polynomial after the BiCG step – KSPBCGSLSetPol()
-ksp_bcgsl_xres
Threshold used to decide when to refresh computed residuals - - KSPBCGSLSetXRes()-ksp_bcgsl_pinv <true/false> - (de)activate use of pseudoinverse – KSPBCGSLSetUsePseudoinverse()
See Also#
KSPCreate()
, KSPSetType()
, KSPType
, KSP
, KSPFGMRES
, KSPBCGS
, KSPSetPCSide()
, KSPBCGSLSetEll()
, KSPBCGSLSetXRes()
Level#
beginner
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