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 singleMPI_Allreduce()
Notes#
The algorithm used in this case is described as Method 1 in [DAzevedoER93]. 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. Pipelined Krylov Methods,
References#
EF D’Azevedo, VL Eijkhout, and CH Romine. Conjugate gradient algorithms with reduced synchronization overhead on distributed memory multiprocessors. Computer Science Technical Report CS-93-185, University of Tennessee, Knoxville, 1993.
See Also#
KSP: Linear System Solvers, Pipelined Krylov Methods, KSP
, KSPCG
, KSPGMRES
, KSPPIPECG
, KSPPIPECR
, and KSPGROPPCG
Level#
intermediate
Location#
Implementations#
KSPCGUseSingleReduction_CG() in src/ksp/ksp/impls/cg/cg.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages