Actual source code: cg.c

petsc-3.7.5 2017-01-01
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>       /*I "petscksp.h" I*/
 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: */
 56: PetscErrorCode KSPSetUp_CG(KSP ksp)
 57: {
 58:   KSP_CG         *cgP = (KSP_CG*)ksp->data;
 60:   PetscInt       maxit = ksp->max_it,nwork = 3;

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

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

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

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

 87: /*
 88:      KSPSolve_CG - This routine actually applies the conjugate gradient method

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

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

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

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

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

127:   ksp->its = 0;
128:   if (!ksp->guess_zero) {
129:     KSP_MatMult(ksp,Amat,X,R);            /*    r <- b - Ax                       */
130:     VecAYPX(R,-1.0,B);
131:   } else {
132:     VecCopy(B,R);                         /*    r <- b (x is 0)                   */
133:   }

135:   switch (ksp->normtype) {
136:     case KSP_NORM_PRECONDITIONED:
137:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
138:       VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
139:       break;
140:     case KSP_NORM_UNPRECONDITIONED:
141:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r = e'*A'*A*e            */
142:       break;
143:     case KSP_NORM_NATURAL:
144:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
145:       VecXDot(Z,R,&beta);                 /*    beta <- z'*r                      */
146:       KSPCheckDot(ksp,beta);
147:       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
148:       break;
149:     case KSP_NORM_NONE:
150:       dp = 0.0;
151:       break;
152:     default: SETERRQ1(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)) {
162:     KSP_PCApply(ksp,R,Z);                /*     z <- Br                           */
163:   }
164:   if (ksp->normtype != KSP_NORM_NATURAL) {
165:     VecXDot(Z,R,&beta);                  /*     beta <- z'*r                      */
166:     KSPCheckDot(ksp,beta);
167:   }

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

200:     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*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:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
213:       VecNorm(R,NORM_2,&dp);              /*     dp <- r'*r                       */
214:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
215:       KSP_PCApply(ksp,R,Z);               /*     z <- Br                          */
216:       VecXDot(Z,R,&beta);                 /*     beta <- r'*z                     */
217:       KSPCheckDot(ksp,beta);
218:       dp = PetscSqrtReal(PetscAbsScalar(beta));
219:     } else {
220:       dp = 0.0;
221:     }
222:     ksp->rnorm = dp;
223:     KSPLogResidualHistory(ksp,dp);
224:     if (eigs) cg->ned = ksp->its;
225:     KSPMonitor(ksp,i+1,dp);
226:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
227:     if (ksp->reason) break;

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

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

243: /*    
244:        KSPSolve_CG_SingleReduction

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

250:        See KSPCGUseSingleReduction_CG()

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

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

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

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

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

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

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

320:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
321:     KSP_PCApply(ksp,R,Z);                 /*    z <- Br                           */
322:   }
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)) {
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) {
350:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
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))) {
367:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
368:       PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
369:       break;
370:     }
371:     a = beta/dpi;                                              /*    a = beta/p'w                      */
372:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
373:     VecAXPY(X,a,P);                       /*    x <- x + ap                       */
374:     VecAXPY(R,-a,W);                      /*    r <- r - aw                       */
375:     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
376:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
377:       KSP_MatMult(ksp,Amat,Z,S);
378:       VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z                        */
379:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
380:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r                        */
381:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
382:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
383:       tmpvecs[0] = S; tmpvecs[1] = R;
384:       KSP_MatMult(ksp,Amat,Z,S);
385:       VecMDot(Z,2,tmpvecs,tmp);          /*    delta <- z'*A*z = r'*B*A*B*r      */
386:       delta = tmp[0]; beta = tmp[1];                           /*    beta <- z'*r                      */
387:       KSPCheckDot(ksp,beta);
388:       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
389:     } else {
390:       dp = 0.0;
391:     }
392:     ksp->rnorm = dp;
393:     KSPLogResidualHistory(ksp,dp);
394:     if (eigs) cg->ned = ksp->its;
395:     KSPMonitor(ksp,i+1,dp);
396:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
397:     if (ksp->reason) break;

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

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

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

430:   /* free space used for singular value calculations */
431:   if (ksp->calc_sings) {
432:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
433:   }
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: */
447: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
448: {
449:   KSP_CG         *cg = (KSP_CG*)ksp->data;
451:   PetscBool      iascii;

454:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
455:   if (iascii) {
456: #if defined(PETSC_USE_COMPLEX)
457:     PetscViewerASCIIPrintf(viewer,"  CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
458: #endif
459:     if (cg->singlereduction) {
460:       PetscViewerASCIIPrintf(viewer,"  CG: using single-reduction variant\n");
461:     }
462:   }
463:   return(0);
464: }

466: /*
467:     KSPSetFromOptions_CG - Checks the options database for options related to the
468:                            conjugate gradient method.
469: */
472: PetscErrorCode KSPSetFromOptions_CG(PetscOptionItems *PetscOptionsObject,KSP ksp)
473: {
475:   KSP_CG         *cg = (KSP_CG*)ksp->data;

478:   PetscOptionsHead(PetscOptionsObject,"KSP CG and CGNE options");
479: #if defined(PETSC_USE_COMPLEX)
480:   PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
481:                           (PetscEnum*)&cg->type,NULL);
482: #endif
483:   PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPIU_Allreduce()","KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,NULL);
484:   PetscOptionsTail();
485:   return(0);
486: }

488: /*
489:     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
490:                       This routine is registered below in KSPCreate_CG() and called from the
491:                       routine KSPCGSetType() (see the file cgtype.c).
492: */
495: static PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
496: {
497:   KSP_CG *cg = (KSP_CG*)ksp->data;

500:   cg->type = type;
501:   return(0);
502: }

504: /*
505:     KSPCGUseSingleReduction_CG

507:     This routine sets a flag to use a variant of CG. Note that (in somewhat
508:     atypical fashion) it also swaps out the routine called when KSPSolve()
509:     is invoked.
510: */
513: static PetscErrorCode  KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg)
514: {
515:   KSP_CG *cg = (KSP_CG*)ksp->data;

518:   cg->singlereduction = flg;
519:   if (cg->singlereduction) {
520:     ksp->ops->solve = KSPSolve_CG_SingleReduction;
521:   } else {
522:     ksp->ops->solve = KSPSolve_CG;
523:   }
524:   return(0);
525: }

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

531:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
532: */
533: /*MC
534:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

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

541:    Level: beginner

543:    Notes: The PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite
544:           Only left preconditioning is supported.

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

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

552:    References:
553: .   1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
554:    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379

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

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

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

576:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
577:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
578:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);

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

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