xref: /petsc/src/ksp/ksp/impls/cg/cg.c (revision d7547e516efde0cd36ffdeebcfafd4768debadcc)
1 
2 /*
3     This file implements the conjugate gradient method in PETSc as part of
4     KSP. You can use this as a starting point for implementing your own
5     Krylov method that is not provided with PETSc.
6 
7     The following basic routines are required for each Krylov method.
8         KSPCreate_XXX()          - Creates the Krylov context
9         KSPSetFromOptions_XXX()  - Sets runtime options
10         KSPSolve_XXX()           - Runs the Krylov method
11         KSPDestroy_XXX()         - Destroys the Krylov context, freeing all
12                                    memory it needed
13     Here the "_XXX" denotes a particular implementation, in this case
14     we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines
15     are actually called via the common user interface routines
16     KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the
17     application code interface remains identical for all preconditioners.
18 
19     Other basic routines for the KSP objects include
20         KSPSetUp_XXX()
21         KSPView_XXX()            - Prints details of solver being used.
22 
23     Detailed Notes:
24     By default, this code implements the CG (Conjugate Gradient) method,
25     which is valid for real symmetric (and complex Hermitian) positive
26     definite matrices. Note that for the complex Hermitian case, the
27     VecDot() arguments within the code MUST remain in the order given
28     for correct computation of inner products.
29 
30     Reference: Hestenes and Steifel, 1952.
31 
32     By switching to the indefinite vector inner product, VecTDot(), the
33     same code is used for the complex symmetric case as well.  The user
34     must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option
35     -ksp_cg_type symmetric to invoke this variant for the complex case.
36     Note, however, that the complex symmetric code is NOT valid for
37     all such matrices ... and thus we don't recommend using this method.
38 */
39 /*
40     cgimpl.h defines the simple data structured used to store information
41     related to the type of matrix (e.g. complex symmetric) being solved and
42     data used during the optional Lanczo process used to compute eigenvalues
43 */
44 #include <../src/ksp/ksp/impls/cg/cgimpl.h> /*I "petscksp.h" I*/
45 extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP, PetscReal *, PetscReal *);
46 extern PetscErrorCode KSPComputeEigenvalues_CG(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);
47 
48 /*
49      KSPSetUp_CG - Sets up the workspace needed by the CG method.
50 
51       This is called once, usually automatically by KSPSolve() or KSPSetUp()
52      but can be called directly by KSPSetUp()
53 */
54 static PetscErrorCode KSPSetUp_CG(KSP ksp)
55 {
56   KSP_CG  *cgP   = (KSP_CG *)ksp->data;
57   PetscInt maxit = ksp->max_it, nwork = 3;
58 
59   PetscFunctionBegin;
60   /* get work vectors needed by CG */
61   if (cgP->singlereduction) nwork += 2;
62   PetscCall(KSPSetWorkVecs(ksp, nwork));
63 
64   /*
65      If user requested computations of eigenvalues then allocate
66      work space needed
67   */
68   if (ksp->calc_sings) {
69     PetscCall(PetscFree4(cgP->e, cgP->d, cgP->ee, cgP->dd));
70     PetscCall(PetscMalloc4(maxit, &cgP->e, maxit, &cgP->d, maxit, &cgP->ee, maxit, &cgP->dd));
71 
72     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
73     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 /*
79      A macro used in the following KSPSolve_CG and KSPSolve_CG_SingleReduction routines
80 */
81 #define VecXDot(x, y, a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x, y, a) : VecTDot(x, y, a))
82 
83 /*
84      KSPSolve_CG - This routine actually applies the conjugate gradient method
85 
86      Note : this routine can be replaced with another one (see below) which implements
87             another variant of CG.
88 
89    Input Parameter:
90 .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
91             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
92 */
93 static PetscErrorCode KSPSolve_CG(KSP ksp)
94 {
95   PetscInt    i, stored_max_it, eigs;
96   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, dpiold;
97   PetscReal   dp = 0.0;
98   Vec         X, B, Z, R, P, W;
99   KSP_CG     *cg;
100   Mat         Amat, Pmat;
101   PetscBool   diagonalscale;
102 
103   PetscFunctionBegin;
104   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
105   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
106 
107   cg            = (KSP_CG *)ksp->data;
108   eigs          = ksp->calc_sings;
109   stored_max_it = ksp->max_it;
110   X             = ksp->vec_sol;
111   B             = ksp->vec_rhs;
112   R             = ksp->work[0];
113   Z             = ksp->work[1];
114   P             = ksp->work[2];
115   W             = Z;
116 
117   if (eigs) {
118     e    = cg->e;
119     d    = cg->d;
120     e[0] = 0.0;
121   }
122   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
123 
124   ksp->its = 0;
125   if (!ksp->guess_zero) {
126     PetscCall(KSP_MatMult(ksp, Amat, X, R)); /*    r <- b - Ax                       */
127     PetscCall(VecAYPX(R, -1.0, B));
128   } else {
129     PetscCall(VecCopy(B, R)); /*    r <- b (x is 0)                   */
130   }
131   /* This may be true only on a subset of MPI ranks; setting it here so it will be detected by the first norm computation below */
132   if (ksp->reason == KSP_DIVERGED_PC_FAILED) PetscCall(VecSetInf(R));
133 
134   switch (ksp->normtype) {
135   case KSP_NORM_PRECONDITIONED:
136     PetscCall(KSP_PCApply(ksp, R, Z));  /*    z <- Br                           */
137     PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z = e'*A'*B'*B*A*e       */
138     KSPCheckNorm(ksp, dp);
139     break;
140   case KSP_NORM_UNPRECONDITIONED:
141     PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r = e'*A'*A*e            */
142     KSPCheckNorm(ksp, dp);
143     break;
144   case KSP_NORM_NATURAL:
145     PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
146     PetscCall(VecXDot(Z, R, &beta));   /*    beta <- z'*r                      */
147     KSPCheckDot(ksp, beta);
148     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
149     break;
150   case KSP_NORM_NONE:
151     dp = 0.0;
152     break;
153   default:
154     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
155   }
156   PetscCall(KSPLogResidualHistory(ksp, dp));
157   PetscCall(KSPMonitor(ksp, 0, dp));
158   ksp->rnorm = dp;
159 
160   PetscCall((*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP)); /* test for convergence */
161   if (ksp->reason) PetscFunctionReturn(0);
162 
163   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                           */ }
164   if (ksp->normtype != KSP_NORM_NATURAL) {
165     PetscCall(VecXDot(Z, R, &beta)); /*     beta <- z'*r                      */
166     KSPCheckDot(ksp, beta);
167   }
168 
169   i = 0;
170   do {
171     ksp->its = i + 1;
172     if (beta == 0.0) {
173       ksp->reason = KSP_CONVERGED_ATOL;
174       PetscCall(PetscInfo(ksp, "converged due to beta = 0\n"));
175       break;
176 #if !defined(PETSC_USE_COMPLEX)
177     } else if ((i > 0) && (beta * betaold < 0.0)) {
178       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite preconditioner, beta %g, betaold %g", (double)beta, (double)betaold);
179       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
180       PetscCall(PetscInfo(ksp, "diverging due to indefinite preconditioner\n"));
181       break;
182 #endif
183     }
184     if (!i) {
185       PetscCall(VecCopy(Z, P)); /*     p <- z                           */
186       b = 0.0;
187     } else {
188       b = beta / betaold;
189       if (eigs) {
190         PetscCheck(ksp->max_it == stored_max_it, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Can not change maxit AND calculate eigenvalues");
191         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
192       }
193       PetscCall(VecAYPX(P, b, Z)); /*     p <- z + b* p                    */
194     }
195     dpiold = dpi;
196     PetscCall(KSP_MatMult(ksp, Amat, P, W)); /*     w <- Ap                          */
197     PetscCall(VecXDot(P, W, &dpi));          /*     dpi <- p'w                       */
198     KSPCheckDot(ksp, dpi);
199     betaold = beta;
200 
201     if ((dpi == 0.0) || ((i > 0) && ((PetscSign(PetscRealPart(dpi)) * PetscSign(PetscRealPart(dpiold))) < 0.0))) {
202       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite matrix, dpi %g, dpiold %g", (double)PetscRealPart(dpi), (double)PetscRealPart(dpiold));
203       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
204       PetscCall(PetscInfo(ksp, "diverging due to indefinite or negative definite matrix\n"));
205       break;
206     }
207     a = beta / dpi; /*     a = beta/p'w                     */
208     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
209     PetscCall(VecAXPY(X, a, P));  /*     x <- x + ap                      */
210     PetscCall(VecAXPY(R, -a, W)); /*     r <- r - aw                      */
211     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
212       PetscCall(KSP_PCApply(ksp, R, Z));  /*     z <- Br                          */
213       PetscCall(VecNorm(Z, NORM_2, &dp)); /*     dp <- z'*z                       */
214       KSPCheckNorm(ksp, dp);
215     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
216       PetscCall(VecNorm(R, NORM_2, &dp)); /*     dp <- r'*r                       */
217       KSPCheckNorm(ksp, dp);
218     } else if (ksp->normtype == KSP_NORM_NATURAL) {
219       PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                          */
220       PetscCall(VecXDot(Z, R, &beta));   /*     beta <- r'*z                     */
221       KSPCheckDot(ksp, beta);
222       dp = PetscSqrtReal(PetscAbsScalar(beta));
223     } else {
224       dp = 0.0;
225     }
226     ksp->rnorm = dp;
227     PetscCall(KSPLogResidualHistory(ksp, dp));
228     PetscCall(KSPMonitor(ksp, i + 1, dp));
229     PetscCall((*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP));
230     if (ksp->reason) break;
231 
232     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i + 2)) { PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                          */ }
233     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i + 2)) {
234       PetscCall(VecXDot(Z, R, &beta)); /*     beta <- z'*r                     */
235       KSPCheckDot(ksp, beta);
236     }
237 
238     i++;
239   } while (i < ksp->max_it);
240   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
241   PetscFunctionReturn(0);
242 }
243 
244 /*
245        KSPSolve_CG_SingleReduction
246 
247        This variant of CG is identical in exact arithmetic to the standard algorithm,
248        but is rearranged to use only a single reduction stage per iteration, using additional
249        intermediate vectors.
250 
251        See KSPCGUseSingleReduction_CG()
252 
253 */
254 static PetscErrorCode KSPSolve_CG_SingleReduction(KSP ksp)
255 {
256   PetscInt    i, stored_max_it, eigs;
257   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, delta, dpiold, tmp[2];
258   PetscReal   dp = 0.0;
259   Vec         X, B, Z, R, P, S, W, tmpvecs[2];
260   KSP_CG     *cg;
261   Mat         Amat, Pmat;
262   PetscBool   diagonalscale;
263 
264   PetscFunctionBegin;
265   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
266   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
267 
268   cg            = (KSP_CG *)ksp->data;
269   eigs          = ksp->calc_sings;
270   stored_max_it = ksp->max_it;
271   X             = ksp->vec_sol;
272   B             = ksp->vec_rhs;
273   R             = ksp->work[0];
274   Z             = ksp->work[1];
275   P             = ksp->work[2];
276   S             = ksp->work[3];
277   W             = ksp->work[4];
278 
279   if (eigs) {
280     e    = cg->e;
281     d    = cg->d;
282     e[0] = 0.0;
283   }
284   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
285 
286   ksp->its = 0;
287   if (!ksp->guess_zero) {
288     PetscCall(KSP_MatMult(ksp, Amat, X, R)); /*    r <- b - Ax                       */
289     PetscCall(VecAYPX(R, -1.0, B));
290   } else {
291     PetscCall(VecCopy(B, R)); /*    r <- b (x is 0)                   */
292   }
293 
294   switch (ksp->normtype) {
295   case KSP_NORM_PRECONDITIONED:
296     PetscCall(KSP_PCApply(ksp, R, Z));  /*    z <- Br                           */
297     PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
298     KSPCheckNorm(ksp, dp);
299     break;
300   case KSP_NORM_UNPRECONDITIONED:
301     PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r = e'*A'*A*e            */
302     KSPCheckNorm(ksp, dp);
303     break;
304   case KSP_NORM_NATURAL:
305     PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
306     PetscCall(KSP_MatMult(ksp, Amat, Z, S));
307     PetscCall(VecXDot(Z, S, &delta)); /*    delta <- z'*A*z = r'*B*A*B*r      */
308     PetscCall(VecXDot(Z, R, &beta));  /*    beta <- z'*r                      */
309     KSPCheckDot(ksp, beta);
310     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
311     break;
312   case KSP_NORM_NONE:
313     dp = 0.0;
314     break;
315   default:
316     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
317   }
318   PetscCall(KSPLogResidualHistory(ksp, dp));
319   PetscCall(KSPMonitor(ksp, 0, dp));
320   ksp->rnorm = dp;
321 
322   PetscCall((*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP)); /* test for convergence */
323   if (ksp->reason) PetscFunctionReturn(0);
324 
325   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */ }
326   if (ksp->normtype != KSP_NORM_NATURAL) {
327     PetscCall(KSP_MatMult(ksp, Amat, Z, S));
328     PetscCall(VecXDot(Z, S, &delta)); /*    delta <- z'*A*z = r'*B*A*B*r      */
329     PetscCall(VecXDot(Z, R, &beta));  /*    beta <- z'*r                      */
330     KSPCheckDot(ksp, beta);
331   }
332 
333   i = 0;
334   do {
335     ksp->its = i + 1;
336     if (beta == 0.0) {
337       ksp->reason = KSP_CONVERGED_ATOL;
338       PetscCall(PetscInfo(ksp, "converged due to beta = 0\n"));
339       break;
340 #if !defined(PETSC_USE_COMPLEX)
341     } else if ((i > 0) && (beta * betaold < 0.0)) {
342       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite preconditioner");
343       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
344       PetscCall(PetscInfo(ksp, "diverging due to indefinite preconditioner\n"));
345       break;
346 #endif
347     }
348     if (!i) {
349       PetscCall(VecCopy(Z, P)); /*    p <- z                           */
350       b = 0.0;
351     } else {
352       b = beta / betaold;
353       if (eigs) {
354         PetscCheck(ksp->max_it == stored_max_it, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Can not change maxit AND calculate eigenvalues");
355         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
356       }
357       PetscCall(VecAYPX(P, b, Z)); /*    p <- z + b* p                     */
358     }
359     dpiold = dpi;
360     if (!i) {
361       PetscCall(KSP_MatMult(ksp, Amat, P, W)); /*    w <- Ap                           */
362       PetscCall(VecXDot(P, W, &dpi));          /*    dpi <- p'w                        */
363     } else {
364       PetscCall(VecAYPX(W, beta / betaold, S));                 /*    w <- Ap                           */
365       dpi = delta - beta * beta * dpiold / (betaold * betaold); /*    dpi <- p'w                        */
366     }
367     betaold = beta;
368     KSPCheckDot(ksp, beta);
369 
370     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi * dpiold) <= 0.0))) {
371       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite matrix");
372       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
373       PetscCall(PetscInfo(ksp, "diverging due to indefinite or negative definite matrix\n"));
374       break;
375     }
376     a = beta / dpi; /*    a = beta/p'w                      */
377     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
378     PetscCall(VecAXPY(X, a, P));  /*    x <- x + ap                       */
379     PetscCall(VecAXPY(R, -a, W)); /*    r <- r - aw                       */
380     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
381       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
382       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
383       PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z                        */
384       KSPCheckNorm(ksp, dp);
385     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
386       PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r                        */
387       KSPCheckNorm(ksp, dp);
388     } else if (ksp->normtype == KSP_NORM_NATURAL) {
389       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
390       tmpvecs[0] = S;
391       tmpvecs[1] = R;
392       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
393       PetscCall(VecMDot(Z, 2, tmpvecs, tmp)); /*    delta <- z'*A*z = r'*B*A*B*r      */
394       delta = tmp[0];
395       beta  = tmp[1]; /*    beta <- z'*r                      */
396       KSPCheckDot(ksp, beta);
397       dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
398     } else {
399       dp = 0.0;
400     }
401     ksp->rnorm = dp;
402     PetscCall(KSPLogResidualHistory(ksp, dp));
403     PetscCall(KSPMonitor(ksp, i + 1, dp));
404     PetscCall((*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP));
405     if (ksp->reason) break;
406 
407     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i + 2)) {
408       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
409       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
410     }
411     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i + 2)) {
412       tmpvecs[0] = S;
413       tmpvecs[1] = R;
414       PetscCall(VecMDot(Z, 2, tmpvecs, tmp));
415       delta = tmp[0];
416       beta  = tmp[1];         /*    delta <- z'*A*z = r'*B'*A*B*r     */
417       KSPCheckDot(ksp, beta); /*    beta <- z'*r                      */
418     }
419 
420     i++;
421   } while (i < ksp->max_it);
422   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
423   PetscFunctionReturn(0);
424 }
425 
426 /*
427      KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function
428                      compositions from KSPCreate_CG. If adding your own KSP implementation,
429                      you must be sure to free all allocated resources here to prevent
430                      leaks.
431 */
432 PetscErrorCode KSPDestroy_CG(KSP ksp)
433 {
434   KSP_CG *cg = (KSP_CG *)ksp->data;
435 
436   PetscFunctionBegin;
437   PetscCall(PetscFree4(cg->e, cg->d, cg->ee, cg->dd));
438   PetscCall(KSPDestroyDefault(ksp));
439   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", NULL));
440   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", NULL));
441   PetscFunctionReturn(0);
442 }
443 
444 /*
445      KSPView_CG - Prints information about the current Krylov method being used.
446                   If your Krylov method has special options or flags that information
447                   should be printed here.
448 */
449 PetscErrorCode KSPView_CG(KSP ksp, PetscViewer viewer)
450 {
451   KSP_CG   *cg = (KSP_CG *)ksp->data;
452   PetscBool iascii;
453 
454   PetscFunctionBegin;
455   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
456   if (iascii) {
457 #if defined(PETSC_USE_COMPLEX)
458     PetscCall(PetscViewerASCIIPrintf(viewer, "  variant %s\n", KSPCGTypes[cg->type]));
459 #endif
460     if (cg->singlereduction) PetscCall(PetscViewerASCIIPrintf(viewer, "  using single-reduction variant\n"));
461   }
462   PetscFunctionReturn(0);
463 }
464 
465 /*
466     KSPSetFromOptions_CG - Checks the options database for options related to the
467                            conjugate gradient method.
468 */
469 PetscErrorCode KSPSetFromOptions_CG(KSP ksp, PetscOptionItems *PetscOptionsObject)
470 {
471   KSP_CG   *cg = (KSP_CG *)ksp->data;
472   PetscBool flg;
473 
474   PetscFunctionBegin;
475   PetscOptionsHeadBegin(PetscOptionsObject, "KSP CG and CGNE options");
476 #if defined(PETSC_USE_COMPLEX)
477   PetscCall(PetscOptionsEnum("-ksp_cg_type", "Matrix is Hermitian or complex symmetric", "KSPCGSetType", KSPCGTypes, (PetscEnum)cg->type, (PetscEnum *)&cg->type, NULL));
478 #endif
479   PetscCall(PetscOptionsBool("-ksp_cg_single_reduction", "Merge inner products into single MPI_Allreduce()", "KSPCGUseSingleReduction", cg->singlereduction, &cg->singlereduction, &flg));
480   if (flg) PetscCall(KSPCGUseSingleReduction(ksp, cg->singlereduction));
481   PetscOptionsHeadEnd();
482   PetscFunctionReturn(0);
483 }
484 
485 /*
486     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
487                       This routine is registered below in KSPCreate_CG() and called from the
488                       routine KSPCGSetType() (see the file cgtype.c).
489 */
490 PetscErrorCode KSPCGSetType_CG(KSP ksp, KSPCGType type)
491 {
492   KSP_CG *cg = (KSP_CG *)ksp->data;
493 
494   PetscFunctionBegin;
495   cg->type = type;
496   PetscFunctionReturn(0);
497 }
498 
499 /*
500     KSPCGUseSingleReduction_CG
501 
502     This routine sets a flag to use a variant of CG. Note that (in somewhat
503     atypical fashion) it also swaps out the routine called when KSPSolve()
504     is invoked.
505 */
506 static PetscErrorCode KSPCGUseSingleReduction_CG(KSP ksp, PetscBool flg)
507 {
508   KSP_CG *cg = (KSP_CG *)ksp->data;
509 
510   PetscFunctionBegin;
511   cg->singlereduction = flg;
512   if (cg->singlereduction) {
513     ksp->ops->solve = KSPSolve_CG_SingleReduction;
514   } else {
515     ksp->ops->solve = KSPSolve_CG;
516   }
517   PetscFunctionReturn(0);
518 }
519 
520 PETSC_INTERN PetscErrorCode KSPBuildResidual_CG(KSP ksp, Vec t, Vec v, Vec *V)
521 {
522   PetscFunctionBegin;
523   PetscCall(VecCopy(ksp->work[0], v));
524   *V = v;
525   PetscFunctionReturn(0);
526 }
527 
528 /*MC
529      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method
530 
531    Options Database Keys:
532 +   -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see `KSPCGSetType()`
533 .   -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
534 -   -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single `MPI_Allreduce()` call, see `KSPCGUseSingleReduction()`
535 
536    Level: beginner
537 
538    Notes:
539     The PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite.
540 
541    Only left preconditioning is supported; there are several ways to motivate preconditioned CG, but they all produce the same algorithm.
542    One can interpret preconditioning A with B to mean any of the following\:
543 .n  (1) Solve a left-preconditioned system BAx = Bb, using inv(B) to define an inner product in the algorithm.
544 .n  (2) Solve a right-preconditioned system ABy = b, x = By, using B to define an inner product in the algorithm.
545 .n  (3) Solve a symmetrically-preconditioned system, E^TAEy = E^Tb, x = Ey, where B = EE^T.
546 .n  (4) Solve Ax=b with CG, but use the inner product defined by B to define the method [2].
547 .n  In all cases, the resulting algorithm only requires application of B to vectors.
548 
549    For complex numbers there are two different CG methods, one for Hermitian symmetric matrices and one for non-Hermitian symmetric matrices. Use
550    `KSPCGSetType()` to indicate which type you are using.
551 
552    One can use `KSPSetComputeEigenvalues()` and `KSPComputeEigenvalues()` to compute the eigenvalues of the (preconditioned) operator
553 
554    Developer Notes:
555     KSPSolve_CG() should actually query the matrix to determine if it is Hermitian symmetric or not and NOT require the user to
556    indicate it to the `KSP` object.
557 
558    References:
559 +  * - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
560    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
561 -  * - Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs,
562     SIAM, 2014.
563 
564 .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSetComputeEigenvalues()`, `KSPComputeEigenvalues()`
565           `KSPCGSetType()`, `KSPCGUseSingleReduction()`, `KSPPIPECG`, `KSPGROPPCG`
566 M*/
567 
568 /*
569     KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
570        function pointers for all the routines it needs to call (KSPSolve_CG() etc)
571 
572     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
573 */
574 PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
575 {
576   KSP_CG *cg;
577 
578   PetscFunctionBegin;
579   PetscCall(PetscNew(&cg));
580 #if !defined(PETSC_USE_COMPLEX)
581   cg->type = KSP_CG_SYMMETRIC;
582 #else
583   cg->type = KSP_CG_HERMITIAN;
584 #endif
585   ksp->data = (void *)cg;
586 
587   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
588   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
589   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2));
590   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
591 
592   /*
593        Sets the functions that are associated with this data structure
594        (in C++ this is the same as defining virtual functions)
595   */
596   ksp->ops->setup          = KSPSetUp_CG;
597   ksp->ops->solve          = KSPSolve_CG;
598   ksp->ops->destroy        = KSPDestroy_CG;
599   ksp->ops->view           = KSPView_CG;
600   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
601   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
602   ksp->ops->buildresidual  = KSPBuildResidual_CG;
603 
604   /*
605       Attach the function KSPCGSetType_CG() to this object. The routine
606       KSPCGSetType() checks for this attached function and calls it if it finds
607       it. (Sort of like a dynamic member function that can be added at run time
608   */
609   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", KSPCGSetType_CG));
610   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", KSPCGUseSingleReduction_CG));
611   PetscFunctionReturn(0);
612 }
613