petsc-3.13.6 2020-09-29
KSPPIPECGRR
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 Section 1.5 Writing Application Codes with PETSc.
KSPPIPECGRR improves the robustness of KSPPIPECG 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.
See also KSPPIPECG, which is identical to KSPPIPECGRR without residual replacements.
See also KSPPIPECR, where the reduction is only overlapped with the matrix-vector product.
Notes
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.
Contributed by
Siegfried Cools, Universiteit Antwerpen, Dept. Mathematics & Computer Science,
European FP7 Project on EXascale Algorithms and Advanced Computational Techniques (EXA2CT) / Research Foundation Flanders (FWO)
Reference
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.
See Also
KSPCreate(), KSPSetType(), KSPPIPECR, KSPGROPPCG, KSPPIPECG, KSPPGMRES, KSPCG, KSPPIPEBCGS, KSPCGUseSingleReduction()
Level
intermediate
Location
src/ksp/ksp/impls/cg/pipecgrr/pipecgrr.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages