KSPCG#

The Preconditioned Conjugate Gradient (PCG) iterative method [HS52] and [MalekStrakovs14] for solving linear systems using KSP.

Options Database Keys#

  • -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see KSPCGSetType()

  • -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric

  • -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single MPI_Allreduce() call, see KSPCGUseSingleReduction()

Notes#

The KSPCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite.

KSPCG is the best Krylov method, KSPType, when the matrix and preconditioner are symmetric positive definite (SPD).

Only left preconditioning is supported with KSPCG; there are several ways to motivate preconditioned CG, but they all produce the same algorithm. One can interpret preconditioning AA with BB to mean any of the following:

   (1) Solve a left-preconditioned system $BAx = Bb $, using $ B^{-1}$ to define an inner product in the algorithm.
   (2) Solve a right-preconditioned system $ABy = b, x = By,$ using $B$ to define an inner product in the algorithm.
   (3) Solve a symmetrically-preconditioned system, $ E^TAEy = E^Tb, x = Ey, $ where $B = EE^T.$
   (4) Solve $Ax=b$ with CG, but use the inner product defined by $B$ to define the method.
   In all cases, the resulting algorithm only requires application of $B$ to vectors.

For complex numbers there are two different CG methods, one for Hermitian symmetric matrices and one for non-Hermitian symmetric matrices. Use KSPCGSetType() to indicate which type you are using.

One can use KSPSetComputeEigenvalues() and KSPComputeEigenvalues() to compute the eigenvalues of the (preconditioned) operator

There are two pipelined implementations of CG in PETSc KSPPIPECG and KSPGROPPCG. These may perform better for very large numbers of MPI processes since they overlap communication and computation so the reduction operations in CG, that is inner products and norms, do not dominate the compute time.

Developer Note#

KSPSolve_CG() should actually query the matrix to determine if it is Hermitian or symmetric and NOT require the user to indicate it to the KSP object.

References#

[HS52]

Magnus R. Hestenes and Eduard Steifel. Methods of conjugate gradients for solving linear systems. J. Research of the National Bureau of Standards, 49:409–436, 1952.

[MalekStrakovs14]

Josef Málek and Zdeněk Strakoš. Preconditioning and the conjugate gradient method in the context of solving PDEs. SIAM, 2014.

See Also#

KSP: Linear System Solvers, KSPCreate(), KSPSetType(), KSPType, KSP, KSPSetComputeEigenvalues(), KSPComputeEigenvalues() KSPCGSetType(), KSPCGUseSingleReduction(), KSPPIPECG, KSPGROPPCG

Level#

beginner

Location#

src/ksp/ksp/impls/cg/cg.c

Examples#

src/snes/tutorials/ex11.c
src/ksp/ksp/tutorials/ex71.c
src/tao/pde_constrained/tutorials/parabolic.c
src/ksp/pc/tutorials/ex1.c
src/ksp/ksp/tutorials/ex78.c
src/ksp/ksp/tutorials/ex59.c
src/tao/bound/tutorials/jbearing2.c
src/ksp/pc/tutorials/ex2.c


Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages