1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPBCGSL.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>KSPBCGSL</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.4.2 2013-07-02</b></div>
<A NAME="KSPBCGSL"><H1>KSPBCGSL</H1></A>
Implements a slight variant of the Enhanced BiCGStab(L) algorithm in (3) and (2). The variation concerns cases when either kappa0**2 or kappa1**2 is negative due to round-off. Kappa0 has also been pulled out of the denominator in the formula for ghat.
<H3><FONT COLOR="#CC3333">References</FONT></H3>
1. 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), pp 235-54, 1996.
2. G.L.G. Sleijpen, H.A. van der Vorst, D.R. Fokkema,
"BiCGStab(L) and other hybrid Bi-CG methods",
Numerical Algorithms, 7, pp 75-109, 1994.
3. D.R. Fokkema, "Enhanced implementation of BiCGStab(L)
for solving linear systems of equations", preprint
from www.citeseer.com.
<P>
Contributed by: Joel M. Malard, email jm.malard@pnl.gov
<P>
<H3><FONT COLOR="#CC3333">Options Database Keys</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_bcgsl_ell <ell> Number of Krylov search directions, defaults to 2 </B></TD><TD>- - <A HREF="../KSP/KSPBCGSLSetEll.html#KSPBCGSLSetEll">KSPBCGSLSetEll</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_bcgsl_cxpol </B></TD><TD>- Use a convex function of the MinRes and OR polynomials after the BiCG step instead of default MinRes -- <A HREF="../KSP/KSPBCGSLSetPol.html#KSPBCGSLSetPol">KSPBCGSLSetPol</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_bcgsl_mrpoly </B></TD><TD>- Use the default MinRes polynomial after the BiCG step -- <A HREF="../KSP/KSPBCGSLSetPol.html#KSPBCGSLSetPol">KSPBCGSLSetPol</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_bcgsl_xres <res> Threshold used to decide when to refresh computed residuals </B></TD><TD>- - <A HREF="../KSP/KSPBCGSLSetXRes.html#KSPBCGSLSetXRes">KSPBCGSLSetXRes</A>()
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>-ksp_bcgsl_pinv <true/false> </B></TD><TD>- (de)activate use of pseudoinverse -- <A HREF="../KSP/KSPBCGSLSetUsePseudoinverse.html#KSPBCGSLSetUsePseudoinverse">KSPBCGSLSetUsePseudoinverse</A>()
</TD></TR></TABLE>
<P>
<P>
<H3><FONT COLOR="#CC3333">See Also</FONT></H3>
<A HREF="../KSP/KSPCreate.html#KSPCreate">KSPCreate</A>(), <A HREF="../KSP/KSPSetType.html#KSPSetType">KSPSetType</A>(), <A HREF="../KSP/KSPType.html#KSPType">KSPType</A> (for list of available types), <A HREF="../KSP/KSP.html#KSP">KSP</A>, <A HREF="../KSP/KSPFGMRES.html#KSPFGMRES">KSPFGMRES</A>, <A HREF="../KSP/KSPBCGS.html#KSPBCGS">KSPBCGS</A>, <A HREF="../KSP/KSPSetPCSide.html#KSPSetPCSide">KSPSetPCSide</A>(), <A HREF="../KSP/KSPBCGSLSetEll.html#KSPBCGSLSetEll">KSPBCGSLSetEll</A>(), <A HREF="../KSP/KSPBCGSLSetXRes.html#KSPBCGSLSetXRes">KSPBCGSLSetXRes</A>()
<BR>
<P>
<P><B><P><B><FONT COLOR="#CC3333">Level:</FONT></B>beginner
<BR><FONT COLOR="#CC3333">Location:</FONT></B><A HREF="../../../src/ksp/ksp/impls/bcgsl/bcgsl.c.html#KSPBCGSL">src/ksp/ksp/impls/bcgsl/bcgsl.c</A>
<BR><A HREF="./index.html">Index of all KSP routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|