Actual source code: cg.c

petsc-3.9.4 2018-09-11
Report Typos and Errors

  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 are
 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;
 58:   PetscInt       maxit = ksp->max_it,nwork = 3;

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

 65:   /*
 66:      If user requested computations of eigenvalues then allocate work
 67:      work space needed
 68:   */
 69:   if (ksp->calc_sings) {
 70:     /* get space to store tridiagonal matrix for Lanczos */
 71:     PetscMalloc4(maxit+1,&cgP->e,maxit+1,&cgP->d,maxit+1,&cgP->ee,maxit+1,&cgP->dd);
 72:     PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));

 74:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 75:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 76:   }
 77:   return(0);
 78: }

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

 85: /*
 86:      KSPSolve_CG - This routine actually applies the conjugate gradient method

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

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

107:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
108:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

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

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

123:   ksp->its = 0;
124:   if (!ksp->guess_zero) {
125:     KSP_MatMult(ksp,Amat,X,R);            /*    r <- b - Ax                       */
126:     VecAYPX(R,-1.0,B);
127:   } else {
128:     VecCopy(B,R);                         /*    r <- b (x is 0)                   */
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:       break;
136:     case KSP_NORM_UNPRECONDITIONED:
137:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r = e'*A'*A*e            */
138:       break;
139:     case KSP_NORM_NATURAL:
140:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
141:       VecXDot(Z,R,&beta);                 /*    beta <- z'*r                      */
142:       KSPCheckDot(ksp,beta);
143:       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
144:       break;
145:     case KSP_NORM_NONE:
146:       dp = 0.0;
147:       break;
148:     default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
149:   }
150:   KSPLogResidualHistory(ksp,dp);
151:   KSPMonitor(ksp,0,dp);
152:   ksp->rnorm = dp;

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

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

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

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

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

233:     i++;
234:   } while (i<ksp->max_it);
235:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
236:   return(0);
237: }

239: /*    
240:        KSPSolve_CG_SingleReduction

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

246:        See KSPCGUseSingleReduction_CG()

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

261:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
262:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

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

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

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

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

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

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

324:   i = 0;
325:   do {
326:     ksp->its = i+1;
327:     if (beta == 0.0) {
328:       ksp->reason = KSP_CONVERGED_ATOL;
329:       PetscInfo(ksp,"converged due to beta = 0\n");
330:       break;
331: #if !defined(PETSC_USE_COMPLEX)
332:     } else if ((i > 0) && (beta*betaold < 0.0)) {
333:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
334:       PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
335:       break;
336: #endif
337:     }
338:     if (!i) {
339:       VecCopy(Z,P);                       /*    p <- z                           */
340:       b    = 0.0;
341:     } else {
342:       b = beta/betaold;
343:       if (eigs) {
344:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
345:         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
346:       }
347:       VecAYPX(P,b,Z);                     /*    p <- z + b* p                     */
348:     }
349:     dpiold = dpi;
350:     if (!i) {
351:       KSP_MatMult(ksp,Amat,P,W);          /*    w <- Ap                           */
352:       VecXDot(P,W,&dpi);                  /*    dpi <- p'w                        */
353:     } else {
354:       VecAYPX(W,beta/betaold,S);          /*    w <- Ap                           */
355:       dpi  = delta - beta*beta*dpiold/(betaold*betaold);       /*    dpi <- p'w                        */
356:     }
357:     betaold = beta;
358:     KSPCheckDot(ksp,beta);

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

393:     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) {
394:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
395:       KSP_MatMult(ksp,Amat,Z,S);
396:     }
397:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
398:       tmpvecs[0] = S; tmpvecs[1] = R;
399:       VecMDot(Z,2,tmpvecs,tmp);
400:       delta = tmp[0]; beta = tmp[1];                           /*    delta <- z'*A*z = r'*B'*A*B*r     */
401:       KSPCheckDot(ksp,beta);                                   /*    beta <- z'*r                      */
402:     }

404:     i++;
405:   } while (i<ksp->max_it);
406:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
407:   return(0);
408: }

410: /*
411:      KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function
412:                      compositions from KSPCreate_CG. If adding your own KSP implementation,
413:                      you must be sure to free all allocated resources here to prevent
414:                      leaks.
415: */
416: PetscErrorCode KSPDestroy_CG(KSP ksp)
417: {
418:   KSP_CG         *cg = (KSP_CG*)ksp->data;

422:   /* free space used for singular value calculations */
423:   if (ksp->calc_sings) {
424:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
425:   }
426:   KSPDestroyDefault(ksp);
427:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);
428:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);
429:   return(0);
430: }

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

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

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

466:   PetscOptionsHead(PetscOptionsObject,"KSP CG and CGNE options");
467: #if defined(PETSC_USE_COMPLEX)
468:   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 MPIU_Allreduce()","KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,NULL);
472:   PetscOptionsTail();
473:   return(0);
474: }

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

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

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

515:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
516: */
517: /*MC
518:      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method

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

525:    Level: beginner

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

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

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

543:    References:
544: .   1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
545:    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
546: .   2. - Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs, 
547:     SIAM, 2014.

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

552: M*/
553: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
554: {
556:   KSP_CG         *cg;

559:   PetscNewLog(ksp,&cg);
560: #if !defined(PETSC_USE_COMPLEX)
561:   cg->type = KSP_CG_SYMMETRIC;
562: #else
563:   cg->type = KSP_CG_HERMITIAN;
564: #endif
565:   ksp->data = (void*)cg;

567:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
568:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
569:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
570:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

572:   /*
573:        Sets the functions that are associated with this data structure
574:        (in C++ this is the same as defining virtual functions)
575:   */
576:   ksp->ops->setup          = KSPSetUp_CG;
577:   ksp->ops->solve          = KSPSolve_CG;
578:   ksp->ops->destroy        = KSPDestroy_CG;
579:   ksp->ops->view           = KSPView_CG;
580:   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
581:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
582:   ksp->ops->buildresidual  = KSPBuildResidualDefault;

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