Actual source code: cg.c


  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.

  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.

 19:     Other basic routines for the KSP objects include
 20:         KSPSetUp_XXX()
 21:         KSPView_XXX()            - Prints details of solver being used.

 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.

 30:     Reference: Hestenes and Steifel, 1952.

 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>
 45: extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP, PetscReal *, PetscReal *);
 46: extern PetscErrorCode KSPComputeEigenvalues_CG(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);

 48: /*
 49:      KSPSetUp_CG - Sets up the workspace needed by the CG method.

 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;

 59:   /* get work vectors needed by CG */
 60:   if (cgP->singlereduction) nwork += 2;
 61:   KSPSetWorkVecs(ksp, nwork);

 63:   /*
 64:      If user requested computations of eigenvalues then allocate
 65:      work space needed
 66:   */
 67:   if (ksp->calc_sings) {
 68:     PetscFree4(cgP->e, cgP->d, cgP->ee, cgP->dd);
 69:     PetscMalloc4(maxit, &cgP->e, maxit, &cgP->d, maxit, &cgP->ee, maxit, &cgP->dd);

 71:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 72:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 73:   }
 74:   return 0;
 75: }

 77: /*
 78:      A macro used in the following KSPSolve_CG and KSPSolve_CG_SingleReduction routines
 79: */
 80: #define VecXDot(x, y, a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x, y, a) : VecTDot(x, y, a))

 82: /*
 83:      KSPSolve_CG - This routine actually applies the conjugate gradient method

 85:      Note : this routine can be replaced with another one (see below) which implements
 86:             another variant of CG.

 88:    Input Parameter:
 89: .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
 90:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
 91: */
 92: static PetscErrorCode KSPSolve_CG(KSP ksp)
 93: {
 94:   PetscInt    i, stored_max_it, eigs;
 95:   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, dpiold;
 96:   PetscReal   dp = 0.0;
 97:   Vec         X, B, Z, R, P, W;
 98:   KSP_CG     *cg;
 99:   Mat         Amat, Pmat;
100:   PetscBool   diagonalscale;

102:   PCGetDiagonalScale(ksp->pc, &diagonalscale);

105:   cg            = (KSP_CG *)ksp->data;
106:   eigs          = ksp->calc_sings;
107:   stored_max_it = ksp->max_it;
108:   X             = ksp->vec_sol;
109:   B             = ksp->vec_rhs;
110:   R             = ksp->work[0];
111:   Z             = ksp->work[1];
112:   P             = ksp->work[2];
113:   W             = Z;

115:   if (eigs) {
116:     e    = cg->e;
117:     d    = cg->d;
118:     e[0] = 0.0;
119:   }
120:   PCGetOperators(ksp->pc, &Amat, &Pmat);

122:   ksp->its = 0;
123:   if (!ksp->guess_zero) {
124:     KSP_MatMult(ksp, Amat, X, R); /*    r <- b - Ax                       */
125:     VecAYPX(R, -1.0, B);
126:   } else {
127:     VecCopy(B, R); /*    r <- b (x is 0)                   */
128:   }
129:   /* 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 */
130:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) VecSetInf(R);

132:   switch (ksp->normtype) {
133:   case KSP_NORM_PRECONDITIONED:
134:     KSP_PCApply(ksp, R, Z);  /*    z <- Br                           */
135:     VecNorm(Z, NORM_2, &dp); /*    dp <- z'*z = e'*A'*B'*B*A*e       */
136:     KSPCheckNorm(ksp, dp);
137:     break;
138:   case KSP_NORM_UNPRECONDITIONED:
139:     VecNorm(R, NORM_2, &dp); /*    dp <- r'*r = e'*A'*A*e            */
140:     KSPCheckNorm(ksp, dp);
141:     break;
142:   case KSP_NORM_NATURAL:
143:     KSP_PCApply(ksp, R, Z); /*    z <- Br                           */
144:     VecXDot(Z, R, &beta);   /*    beta <- z'*r                      */
145:     KSPCheckDot(ksp, beta);
146:     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
147:     break;
148:   case KSP_NORM_NONE:
149:     dp = 0.0;
150:     break;
151:   default:
152:     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
153:   }
154:   KSPLogResidualHistory(ksp, dp);
155:   KSPMonitor(ksp, 0, dp);
156:   ksp->rnorm = dp;

158:   (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP); /* test for convergence */
159:   if (ksp->reason) return 0;

161:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { KSP_PCApply(ksp, R, Z); /*     z <- Br                           */ }
162:   if (ksp->normtype != KSP_NORM_NATURAL) {
163:     VecXDot(Z, R, &beta); /*     beta <- z'*r                      */
164:     KSPCheckDot(ksp, beta);
165:   }

167:   i = 0;
168:   do {
169:     ksp->its = i + 1;
170:     if (beta == 0.0) {
171:       ksp->reason = KSP_CONVERGED_ATOL;
172:       PetscInfo(ksp, "converged due to beta = 0\n");
173:       break;
174: #if !defined(PETSC_USE_COMPLEX)
175:     } else if ((i > 0) && (beta * betaold < 0.0)) {
177:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
178:       PetscInfo(ksp, "diverging due to indefinite preconditioner\n");
179:       break;
180: #endif
181:     }
182:     if (!i) {
183:       VecCopy(Z, P); /*     p <- z                           */
184:       b = 0.0;
185:     } else {
186:       b = beta / betaold;
187:       if (eigs) {
189:         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
190:       }
191:       VecAYPX(P, b, Z); /*     p <- z + b* p                    */
192:     }
193:     dpiold = dpi;
194:     KSP_MatMult(ksp, Amat, P, W); /*     w <- Ap                          */
195:     VecXDot(P, W, &dpi);          /*     dpi <- p'w                       */
196:     KSPCheckDot(ksp, dpi);
197:     betaold = beta;

199:     if ((dpi == 0.0) || ((i > 0) && ((PetscSign(PetscRealPart(dpi)) * PetscSign(PetscRealPart(dpiold))) < 0.0))) {
201:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
202:       PetscInfo(ksp, "diverging due to indefinite or negative definite matrix\n");
203:       break;
204:     }
205:     a = beta / dpi; /*     a = beta/p'w                     */
206:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
207:     VecAXPY(X, a, P);  /*     x <- x + ap                      */
208:     VecAXPY(R, -a, W); /*     r <- r - aw                      */
209:     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
210:       KSP_PCApply(ksp, R, Z);  /*     z <- Br                          */
211:       VecNorm(Z, NORM_2, &dp); /*     dp <- z'*z                       */
212:       KSPCheckNorm(ksp, dp);
213:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
214:       VecNorm(R, NORM_2, &dp); /*     dp <- r'*r                       */
215:       KSPCheckNorm(ksp, dp);
216:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
217:       KSP_PCApply(ksp, R, Z); /*     z <- Br                          */
218:       VecXDot(Z, R, &beta);   /*     beta <- r'*z                     */
219:       KSPCheckDot(ksp, beta);
220:       dp = PetscSqrtReal(PetscAbsScalar(beta));
221:     } else {
222:       dp = 0.0;
223:     }
224:     ksp->rnorm = dp;
225:     KSPLogResidualHistory(ksp, dp);
226:     KSPMonitor(ksp, i + 1, dp);
227:     (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
228:     if (ksp->reason) break;

230:     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i + 2)) { KSP_PCApply(ksp, R, Z); /*     z <- Br                          */ }
231:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i + 2)) {
232:       VecXDot(Z, R, &beta); /*     beta <- z'*r                     */
233:       KSPCheckDot(ksp, beta);
234:     }

236:     i++;
237:   } while (i < ksp->max_it);
238:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
239:   return 0;
240: }

242: /*
243:        KSPSolve_CG_SingleReduction

245:        This variant of CG is identical in exact arithmetic to the standard algorithm,
246:        but is rearranged to use only a single reduction stage per iteration, using additional
247:        intermediate vectors.

249:        See KSPCGUseSingleReduction_CG()

251: */
252: static PetscErrorCode KSPSolve_CG_SingleReduction(KSP ksp)
253: {
254:   PetscInt    i, stored_max_it, eigs;
255:   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, delta, dpiold, tmp[2];
256:   PetscReal   dp = 0.0;
257:   Vec         X, B, Z, R, P, S, W, tmpvecs[2];
258:   KSP_CG     *cg;
259:   Mat         Amat, Pmat;
260:   PetscBool   diagonalscale;

262:   PCGetDiagonalScale(ksp->pc, &diagonalscale);

265:   cg            = (KSP_CG *)ksp->data;
266:   eigs          = ksp->calc_sings;
267:   stored_max_it = ksp->max_it;
268:   X             = ksp->vec_sol;
269:   B             = ksp->vec_rhs;
270:   R             = ksp->work[0];
271:   Z             = ksp->work[1];
272:   P             = ksp->work[2];
273:   S             = ksp->work[3];
274:   W             = ksp->work[4];

276:   if (eigs) {
277:     e    = cg->e;
278:     d    = cg->d;
279:     e[0] = 0.0;
280:   }
281:   PCGetOperators(ksp->pc, &Amat, &Pmat);

283:   ksp->its = 0;
284:   if (!ksp->guess_zero) {
285:     KSP_MatMult(ksp, Amat, X, R); /*    r <- b - Ax                       */
286:     VecAYPX(R, -1.0, B);
287:   } else {
288:     VecCopy(B, R); /*    r <- b (x is 0)                   */
289:   }

291:   switch (ksp->normtype) {
292:   case KSP_NORM_PRECONDITIONED:
293:     KSP_PCApply(ksp, R, Z);  /*    z <- Br                           */
294:     VecNorm(Z, NORM_2, &dp); /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
295:     KSPCheckNorm(ksp, dp);
296:     break;
297:   case KSP_NORM_UNPRECONDITIONED:
298:     VecNorm(R, NORM_2, &dp); /*    dp <- r'*r = e'*A'*A*e            */
299:     KSPCheckNorm(ksp, dp);
300:     break;
301:   case KSP_NORM_NATURAL:
302:     KSP_PCApply(ksp, R, Z); /*    z <- Br                           */
303:     KSP_MatMult(ksp, Amat, Z, S);
304:     VecXDot(Z, S, &delta); /*    delta <- z'*A*z = r'*B*A*B*r      */
305:     VecXDot(Z, R, &beta);  /*    beta <- z'*r                      */
306:     KSPCheckDot(ksp, beta);
307:     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
308:     break;
309:   case KSP_NORM_NONE:
310:     dp = 0.0;
311:     break;
312:   default:
313:     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
314:   }
315:   KSPLogResidualHistory(ksp, dp);
316:   KSPMonitor(ksp, 0, dp);
317:   ksp->rnorm = dp;

319:   (*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP); /* test for convergence */
320:   if (ksp->reason) return 0;

322:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) { KSP_PCApply(ksp, R, Z); /*    z <- Br                           */ }
323:   if (ksp->normtype != KSP_NORM_NATURAL) {
324:     KSP_MatMult(ksp, Amat, Z, S);
325:     VecXDot(Z, S, &delta); /*    delta <- z'*A*z = r'*B*A*B*r      */
326:     VecXDot(Z, R, &beta);  /*    beta <- z'*r                      */
327:     KSPCheckDot(ksp, beta);
328:   }

330:   i = 0;
331:   do {
332:     ksp->its = i + 1;
333:     if (beta == 0.0) {
334:       ksp->reason = KSP_CONVERGED_ATOL;
335:       PetscInfo(ksp, "converged due to beta = 0\n");
336:       break;
337: #if !defined(PETSC_USE_COMPLEX)
338:     } else if ((i > 0) && (beta * betaold < 0.0)) {
340:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
341:       PetscInfo(ksp, "diverging due to indefinite preconditioner\n");
342:       break;
343: #endif
344:     }
345:     if (!i) {
346:       VecCopy(Z, P); /*    p <- z                           */
347:       b = 0.0;
348:     } else {
349:       b = beta / betaold;
350:       if (eigs) {
352:         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
353:       }
354:       VecAYPX(P, b, Z); /*    p <- z + b* p                     */
355:     }
356:     dpiold = dpi;
357:     if (!i) {
358:       KSP_MatMult(ksp, Amat, P, W); /*    w <- Ap                           */
359:       VecXDot(P, W, &dpi);          /*    dpi <- p'w                        */
360:     } else {
361:       VecAYPX(W, beta / betaold, S);                 /*    w <- Ap                           */
362:       dpi = delta - beta * beta * dpiold / (betaold * betaold); /*    dpi <- p'w                        */
363:     }
364:     betaold = beta;
365:     KSPCheckDot(ksp, beta);

367:     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi * dpiold) <= 0.0))) {
369:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
370:       PetscInfo(ksp, "diverging due to indefinite or negative definite matrix\n");
371:       break;
372:     }
373:     a = beta / dpi; /*    a = beta/p'w                      */
374:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
375:     VecAXPY(X, a, P);  /*    x <- x + ap                       */
376:     VecAXPY(R, -a, W); /*    r <- r - aw                       */
377:     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
378:       KSP_PCApply(ksp, R, Z); /*    z <- Br                           */
379:       KSP_MatMult(ksp, Amat, Z, S);
380:       VecNorm(Z, NORM_2, &dp); /*    dp <- z'*z                        */
381:       KSPCheckNorm(ksp, dp);
382:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
383:       VecNorm(R, NORM_2, &dp); /*    dp <- r'*r                        */
384:       KSPCheckNorm(ksp, dp);
385:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
386:       KSP_PCApply(ksp, R, Z); /*    z <- Br                           */
387:       tmpvecs[0] = S;
388:       tmpvecs[1] = R;
389:       KSP_MatMult(ksp, Amat, Z, S);
390:       VecMDot(Z, 2, tmpvecs, tmp); /*    delta <- z'*A*z = r'*B*A*B*r      */
391:       delta = tmp[0];
392:       beta  = tmp[1]; /*    beta <- z'*r                      */
393:       KSPCheckDot(ksp, beta);
394:       dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
395:     } else {
396:       dp = 0.0;
397:     }
398:     ksp->rnorm = dp;
399:     KSPLogResidualHistory(ksp, dp);
400:     KSPMonitor(ksp, i + 1, dp);
401:     (*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP);
402:     if (ksp->reason) break;

404:     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i + 2)) {
405:       KSP_PCApply(ksp, R, Z); /*    z <- Br                           */
406:       KSP_MatMult(ksp, Amat, Z, S);
407:     }
408:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i + 2)) {
409:       tmpvecs[0] = S;
410:       tmpvecs[1] = R;
411:       VecMDot(Z, 2, tmpvecs, tmp);
412:       delta = tmp[0];
413:       beta  = tmp[1];         /*    delta <- z'*A*z = r'*B'*A*B*r     */
414:       KSPCheckDot(ksp, beta); /*    beta <- z'*r                      */
415:     }

417:     i++;
418:   } while (i < ksp->max_it);
419:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
420:   return 0;
421: }

423: /*
424:      KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function
425:                      compositions from KSPCreate_CG. If adding your own KSP implementation,
426:                      you must be sure to free all allocated resources here to prevent
427:                      leaks.
428: */
429: PetscErrorCode KSPDestroy_CG(KSP ksp)
430: {
431:   KSP_CG *cg = (KSP_CG *)ksp->data;

433:   PetscFree4(cg->e, cg->d, cg->ee, cg->dd);
434:   KSPDestroyDefault(ksp);
435:   PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", NULL);
436:   PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", NULL);
437:   return 0;
438: }

440: /*
441:      KSPView_CG - Prints information about the current Krylov method being used.
442:                   If your Krylov method has special options or flags that information
443:                   should be printed here.
444: */
445: PetscErrorCode KSPView_CG(KSP ksp, PetscViewer viewer)
446: {
447:   KSP_CG   *cg = (KSP_CG *)ksp->data;
448:   PetscBool iascii;

450:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
451:   if (iascii) {
452: #if defined(PETSC_USE_COMPLEX)
453:     PetscViewerASCIIPrintf(viewer, "  variant %s\n", KSPCGTypes[cg->type]);
454: #endif
455:     if (cg->singlereduction) PetscViewerASCIIPrintf(viewer, "  using single-reduction variant\n");
456:   }
457:   return 0;
458: }

460: /*
461:     KSPSetFromOptions_CG - Checks the options database for options related to the
462:                            conjugate gradient method.
463: */
464: PetscErrorCode KSPSetFromOptions_CG(KSP ksp, PetscOptionItems *PetscOptionsObject)
465: {
466:   KSP_CG   *cg = (KSP_CG *)ksp->data;
467:   PetscBool flg;

469:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP CG and CGNE options");
470: #if defined(PETSC_USE_COMPLEX)
471:   PetscOptionsEnum("-ksp_cg_type", "Matrix is Hermitian or complex symmetric", "KSPCGSetType", KSPCGTypes, (PetscEnum)cg->type, (PetscEnum *)&cg->type, NULL);
472: #endif
473:   PetscOptionsBool("-ksp_cg_single_reduction", "Merge inner products into single MPI_Allreduce()", "KSPCGUseSingleReduction", cg->singlereduction, &cg->singlereduction, &flg);
474:   if (flg) KSPCGUseSingleReduction(ksp, cg->singlereduction);
475:   PetscOptionsHeadEnd();
476:   return 0;
477: }

479: /*
480:     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
481:                       This routine is registered below in KSPCreate_CG() and called from the
482:                       routine KSPCGSetType() (see the file cgtype.c).
483: */
484: PetscErrorCode KSPCGSetType_CG(KSP ksp, KSPCGType type)
485: {
486:   KSP_CG *cg = (KSP_CG *)ksp->data;

488:   cg->type = type;
489:   return 0;
490: }

492: /*
493:     KSPCGUseSingleReduction_CG

495:     This routine sets a flag to use a variant of CG. Note that (in somewhat
496:     atypical fashion) it also swaps out the routine called when KSPSolve()
497:     is invoked.
498: */
499: static PetscErrorCode KSPCGUseSingleReduction_CG(KSP ksp, PetscBool flg)
500: {
501:   KSP_CG *cg = (KSP_CG *)ksp->data;

503:   cg->singlereduction = flg;
504:   if (cg->singlereduction) {
505:     ksp->ops->solve = KSPSolve_CG_SingleReduction;
506:   } else {
507:     ksp->ops->solve = KSPSolve_CG;
508:   }
509:   return 0;
510: }

512: PETSC_INTERN PetscErrorCode KSPBuildResidual_CG(KSP ksp, Vec t, Vec v, Vec *V)
513: {
514:   VecCopy(ksp->work[0], v);
515:   *V = v;
516:   return 0;
517: }

519: /*MC
520:      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method

522:    Options Database Keys:
523: +   -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see `KSPCGSetType()`
524: .   -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
525: -   -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single `MPI_Allreduce()` call, see `KSPCGUseSingleReduction()`

527:    Level: beginner

529:    Notes:
530:     The PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite.

532:    Only left preconditioning is supported; there are several ways to motivate preconditioned CG, but they all produce the same algorithm.
533:    One can interpret preconditioning A with B to mean any of the following\:
534: .n  (1) Solve a left-preconditioned system BAx = Bb, using inv(B) to define an inner product in the algorithm.
535: .n  (2) Solve a right-preconditioned system ABy = b, x = By, using B to define an inner product in the algorithm.
536: .n  (3) Solve a symmetrically-preconditioned system, E^TAEy = E^Tb, x = Ey, where B = EE^T.
537: .n  (4) Solve Ax=b with CG, but use the inner product defined by B to define the method [2].
538: .n  In all cases, the resulting algorithm only requires application of B to vectors.

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

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

545:    Developer Notes:
546:     KSPSolve_CG() should actually query the matrix to determine if it is Hermitian symmetric or not and NOT require the user to
547:    indicate it to the `KSP` object.

549:    References:
550: +  * - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
551:    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
552: -  * - Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs,
553:     SIAM, 2014.

555: .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSetComputeEigenvalues()`, `KSPComputeEigenvalues()`
556:           `KSPCGSetType()`, `KSPCGUseSingleReduction()`, `KSPPIPECG`, `KSPGROPPCG`
557: M*/

559: /*
560:     KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
561:        function pointers for all the routines it needs to call (KSPSolve_CG() etc)

563:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
564: */
565: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
566: {
567:   KSP_CG *cg;

569:   PetscNew(&cg);
570: #if !defined(PETSC_USE_COMPLEX)
571:   cg->type = KSP_CG_SYMMETRIC;
572: #else
573:   cg->type = KSP_CG_HERMITIAN;
574: #endif
575:   ksp->data = (void *)cg;

577:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3);
578:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2);
579:   KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2);
580:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1);

582:   /*
583:        Sets the functions that are associated with this data structure
584:        (in C++ this is the same as defining virtual functions)
585:   */
586:   ksp->ops->setup          = KSPSetUp_CG;
587:   ksp->ops->solve          = KSPSolve_CG;
588:   ksp->ops->destroy        = KSPDestroy_CG;
589:   ksp->ops->view           = KSPView_CG;
590:   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
591:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
592:   ksp->ops->buildresidual  = KSPBuildResidual_CG;

594:   /*
595:       Attach the function KSPCGSetType_CG() to this object. The routine
596:       KSPCGSetType() checks for this attached function and calls it if it finds
597:       it. (Sort of like a dynamic member function that can be added at run time
598:   */
599:   PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", KSPCGSetType_CG);
600:   PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", KSPCGUseSingleReduction_CG);
601:   return 0;
602: }