:orphan: # KSPCGUseSingleReduction Merge the two inner products needed in `KSPCG` into a single `MPI_Allreduce()` call. ## Synopsis ``` #include "petscksp.h" PetscErrorCode KSPCGUseSingleReduction(KSP ksp, PetscBool flg) ``` Logically Collective ## Input Parameters - ***ksp -*** the iterative context - ***flg -*** turn on or off the single reduction ## Options Database Key - ***-ksp_cg_single_reduction -*** Merge inner products into single `MPI_Allreduce()` ## Notes The algorithm used in this case is described as Method 1 in [1]. V. Eijkhout credits the algorithm initially to Chronopoulos and Gear. It requires two extra work vectors than the conventional implementation in PETSc. See also `KSPPIPECG`, `KSPPIPECR`, and `KSPGROPPCG` that use non-blocking reductions. [](sec_pipelineksp), ## References - ***[1] -*** Lapack Working Note 56, "Conjugate Gradient Algorithms with Reduced Synchronization Overhead Distributed Memory Multiprocessors", by E. F. D'Azevedo, V. L. Eijkhout, and C. H. Romine, December 3, 1999. ## See Also [](ch_ksp), [](sec_pipelineksp), `KSP`, `KSPCG`, `KSPGMRES`, `KSPPIPECG`, `KSPPIPECR`, and `KSPGROPPCG` ## Level intermediate ## Location src/ksp/ksp/impls/cg/cgtype.c ## Implementations KSPCGUseSingleReduction_CG in src/ksp/ksp/impls/cg/cg.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/cg/cgtype.c) [Index of all KSP routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)