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
|
<!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/KSPPIPECGRR.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<TITLE>KSPPIPECGRR</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>petsc-3.10.3 2018-12-18</b></div>
<div id="bugreport" align=right><a href="mailto:petsc-maint@mcs.anl.gov?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: petsc-3.10.3 v3.10.3 docs/manualpages/KSP/KSPPIPECGRR.html "><small>Report Typos and Errors</small></a></div>
<A NAME="KSPPIPECGRR"><H1>KSPPIPECGRR</H1></A>
Pipelined conjugate gradient method with automated residual replacements. This method has only a single non-blocking reduction per iteration, compared to 2 blocking for standard CG. The
non-blocking reduction is overlapped by the matrix-vector product and preconditioner application.
<P>
<A HREF="../KSP/KSPPIPECGRR.html#KSPPIPECGRR">KSPPIPECGRR</A> improves the robustness of <A HREF="../KSP/KSPPIPECG.html#KSPPIPECG">KSPPIPECG</A> by adding an automated residual replacement strategy.
True residual and other auxiliary variables are computed explicitly in a number of dynamically determined
iterations to counteract the accumulation of rounding errors and thus attain a higher maximal final accuracy.
<P>
See also <A HREF="../KSP/KSPPIPECG.html#KSPPIPECG">KSPPIPECG</A>, which is identical to <A HREF="../KSP/KSPPIPECGRR.html#KSPPIPECGRR">KSPPIPECGRR</A> without residual replacements.
See also <A HREF="../KSP/KSPPIPECR.html#KSPPIPECR">KSPPIPECR</A>, where the reduction is only overlapped with the matrix-vector product.
<P>
<P>
<H3><FONT COLOR="#CC3333">Notes</FONT></H3>
MPI configuration may be necessary for reductions to make asynchronous progress, which is important for
performance of pipelined methods. See the FAQ on the PETSc website for details.
<P>
<H3><FONT COLOR="#CC3333">Contributed by</FONT></H3>
Siegfried Cools, Universiteit Antwerpen, Dept. Mathematics & Computer Science,
European FP7 Project on EXascale Algorithms and Advanced Computational Techniques (EXA2CT) / Research Foundation Flanders (FWO)
<P>
<H3><FONT COLOR="#CC3333">Reference</FONT></H3>
S. Cools, E.F. Yetkin, E. Agullo, L. Giraud, W. Vanroose, "Analyzing the effect of local rounding error
propagation on the maximal attainable accuracy of the pipelined Conjugate Gradients method",
SIAM Journal on Matrix Analysis and Applications (SIMAX), 39(1):426–450, 2018.
<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/KSPPIPECR.html#KSPPIPECR">KSPPIPECR</A>, <A HREF="../KSP/KSPGROPPCG.html#KSPGROPPCG">KSPGROPPCG</A>, <A HREF="../KSP/KSPPIPECG.html#KSPPIPECG">KSPPIPECG</A>, <A HREF="../KSP/KSPPGMRES.html#KSPPGMRES">KSPPGMRES</A>, <A HREF="../KSP/KSPCG.html#KSPCG">KSPCG</A>, <A HREF="../KSP/KSPPIPEBCGS.html#KSPPIPEBCGS">KSPPIPEBCGS</A>, <A HREF="../KSP/KSPCGUseSingleReduction.html#KSPCGUseSingleReduction">KSPCGUseSingleReduction</A>()
<BR><P><B></B><H3><FONT COLOR="#CC3333">Level</FONT></H3>intermediate<BR>
<H3><FONT COLOR="#CC3333">Location</FONT></H3>
</B><A HREF="../../../src/ksp/ksp/impls/cg/pipecgrr/pipecgrr.c.html#KSPPIPECGRR">src/ksp/ksp/impls/cg/pipecgrr/pipecgrr.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>
|