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);
 70:     PetscLogObjectMemory((PetscObject)ksp,2*maxit*(sizeof(PetscScalar)+sizeof(PetscReal)));

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

 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))

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

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

 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;

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

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

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

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

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

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

159:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
160:     KSP_PCApply(ksp,R,Z);                /*     z <- Br                           */
161:   }
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)) {
231:       KSP_PCApply(ksp,R,Z);               /*     z <- Br                          */
232:     }
233:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
234:       VecXDot(Z,R,&beta);                 /*     beta <- z'*r                     */
235:       KSPCheckDot(ksp,beta);
236:     }

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

244: /*
245:        KSPSolve_CG_SingleReduction

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.

251:        See KSPCGUseSingleReduction_CG()

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;

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

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

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

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

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

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

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

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

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

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

412:     i++;
413:   } while (i<ksp->max_it);
414:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
415:   return 0;
416: }

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

428:   PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
429:   KSPDestroyDefault(ksp);
430:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);
431:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);
432:   return 0;
433: }

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

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

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

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

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

486:   cg->type = type;
487:   return 0;
488: }

490: /*
491:     KSPCGUseSingleReduction_CG

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

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

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

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

521:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
522: */
523: /*MC
524:      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method

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

531:    Level: beginner

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

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

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

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

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

557: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
558:            KSPCGSetType(), KSPCGUseSingleReduction(), KSPPIPECG, KSPGROPPCG

560: M*/
561: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
562: {
563:   KSP_CG         *cg;

565:   PetscNewLog(ksp,&cg);
566: #if !defined(PETSC_USE_COMPLEX)
567:   cg->type = KSP_CG_SYMMETRIC;
568: #else
569:   cg->type = KSP_CG_HERMITIAN;
570: #endif
571:   ksp->data = (void*)cg;

573:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
574:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
575:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
576:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

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

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