:orphan: # KSPCG The Preconditioned Conjugate Gradient (PCG) iterative method ## 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 PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite. Only left preconditioning is supported; there are several ways to motivate preconditioned CG, but they all produce the same algorithm. One can interpret preconditioning A with B to mean any of the following: (1) Solve a left-preconditioned system BAx = Bb, using inv(B) 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 [2]. 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 ## Developer Notes KSPSolve_CG() should actually query the matrix to determine if it is Hermitian symmetric or not and NOT require the user to indicate it to the `KSP` object. ## References - **** -*** Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems, Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379 - **** -*** Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs, SIAM, 2014. ## See Also [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSetComputeEigenvalues()`, `KSPComputeEigenvalues()` `KSPCGSetType()`, `KSPCGUseSingleReduction()`, `KSPPIPECG`, `KSPGROPPCG` ## Level beginner ## Location src/ksp/ksp/impls/cg/cg.c ## Examples src/ksp/ksp/tutorials/ex59.c
src/ksp/ksp/tutorials/ex71.c
src/ksp/ksp/tutorials/ex78.c
src/ksp/pc/tutorials/ex1.c
src/ksp/pc/tutorials/ex2.c
src/tao/bound/tutorials/jbearing2.c
src/tao/pde_constrained/tutorials/parabolic.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/cg/cg.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)