Actual source code: cg.c

petsc-3.12.5 2020-03-29
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:       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: SETERRQ1(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)) {
176:       if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite preconditioner, beta %g, betaold %g",(double)beta,(double)betaold);
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) {
188:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
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))) {
200:       if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite matrix, dpi %g, dpiold %g",(double)PetscRealPart(dpi),(double)PetscRealPart(dpiold));
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:     if (eigs) cg->ned = ksp->its;
227:     KSPMonitor(ksp,i+1,dp);
228:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
229:     if (ksp->reason) break;

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

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

245: /*    
246:        KSPSolve_CG_SingleReduction

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

252:        See KSPCGUseSingleReduction_CG()

254: */
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:       KSPCheckNorm(ksp,dp);
297:       break;
298:     case KSP_NORM_UNPRECONDITIONED:
299:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r = e'*A'*A*e            */
300:       KSPCheckNorm(ksp,dp);
301:       break;
302:     case KSP_NORM_NATURAL:
303:       KSP_PCApply(ksp,R,Z);               /*    z <- Br                           */
304:       KSP_MatMult(ksp,Amat,Z,S);
305:       VecXDot(Z,S,&delta);                /*    delta <- z'*A*z = r'*B*A*B*r      */
306:       VecXDot(Z,R,&beta);                 /*    beta <- z'*r                      */
307:       KSPCheckDot(ksp,beta);
308:       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
309:       break;
310:     case KSP_NORM_NONE:
311:       dp = 0.0;
312:       break;
313:     default: SETERRQ1(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)) {
323:     KSP_PCApply(ksp,R,Z);                 /*    z <- Br                           */
324:   }
325:   if (ksp->normtype != KSP_NORM_NATURAL) {
326:     KSP_MatMult(ksp,Amat,Z,S);
327:     VecXDot(Z,S,&delta);                  /*    delta <- z'*A*z = r'*B*A*B*r      */
328:     VecXDot(Z,R,&beta);                   /*    beta <- z'*r                      */
329:     KSPCheckDot(ksp,beta);
330:   }

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

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

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

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

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

434:   /* free space used for singular value calculations */
435:   if (ksp->calc_sings) {
436:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
437:   }
438:   KSPDestroyDefault(ksp);
439:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);
440:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);
441:   return(0);
442: }

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

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

468: /*
469:     KSPSetFromOptions_CG - Checks the options database for options related to the
470:                            conjugate gradient method.
471: */
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: */
493: PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
494: {
495:   KSP_CG *cg = (KSP_CG*)ksp->data;

498:   cg->type = type;
499:   return(0);
500: }

502: /*
503:     KSPCGUseSingleReduction_CG

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

514:   cg->singlereduction = flg;
515:   if (cg->singlereduction) {
516:     ksp->ops->solve = KSPSolve_CG_SingleReduction;
517:   } else {
518:     ksp->ops->solve = KSPSolve_CG;
519:   }
520:   return(0);
521: }

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

527:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
528: */
529: /*MC
530:      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method

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

537:    Level: beginner

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

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

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

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

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

566: M*/
567: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
568: {
570:   KSP_CG         *cg;

573:   PetscNewLog(ksp,&cg);
574: #if !defined(PETSC_USE_COMPLEX)
575:   cg->type = KSP_CG_SYMMETRIC;
576: #else
577:   cg->type = KSP_CG_HERMITIAN;
578: #endif
579:   ksp->data = (void*)cg;

581:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
582:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
583:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
584:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

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

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