Actual source code: cg.c

petsc-3.4.5 2014-06-29
  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 vai 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,PetscScalar,&cgP->e,maxit+1,PetscScalar,&cgP->d,maxit+1,PetscReal,&cgP->ee,maxit+1,PetscReal,&cgP->dd);
 74:     PetscLogObjectMemory(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:        KSPSolve_CG - This routine actually applies the conjugate gradient  method

 85:    This routine is MUCH too messy. I has too many options (norm type and single reduction) embedded making the code confusing and likely to be buggy.

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

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

109:   cg            = (KSP_CG*)ksp->data;
110:   eigs          = ksp->calc_sings;
111:   stored_max_it = ksp->max_it;
112:   X             = ksp->vec_sol;
113:   B             = ksp->vec_rhs;
114:   R             = ksp->work[0];
115:   Z             = ksp->work[1];
116:   P             = ksp->work[2];
117:   if (cg->singlereduction) {
118:     S = ksp->work[3];
119:     W = ksp->work[4];
120:   } else {
121:     S = 0;                      /* unused */
122:     W = Z;
123:   }

125: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))

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

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

138:   switch (ksp->normtype) {
139:   case KSP_NORM_PRECONDITIONED:
140:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
141:     VecNorm(Z,NORM_2,&dp);                /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
142:     break;
143:   case KSP_NORM_UNPRECONDITIONED:
144:     VecNorm(R,NORM_2,&dp);                /*    dp <- r'*r = e'*A'*A*e            */
145:     break;
146:   case KSP_NORM_NATURAL:
147:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
148:     if (cg->singlereduction) {
149:       KSP_MatMult(ksp,Amat,Z,S);
150:       VecXDot(Z,S,&delta);
151:     }
152:     VecXDot(Z,R,&beta);                     /*  beta <- z'*r       */
153:     if (PetscIsInfOrNanScalar(beta)) {
154:       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
155:       else {
156:         ksp->reason = KSP_DIVERGED_NANORINF;
157:         return(0);
158:       }
159:     }
160:     dp = PetscSqrtReal(PetscAbsScalar(beta));                           /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
161:     break;
162:   case KSP_NORM_NONE:
163:     dp = 0.0;
164:     break;
165:   default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
166:   }
167:   KSPLogResidualHistory(ksp,dp);
168:   KSPMonitor(ksp,0,dp);
169:   ksp->rnorm = dp;

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

174:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
175:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
176:   }
177:   if (ksp->normtype != KSP_NORM_NATURAL) {
178:     if (cg->singlereduction) {
179:       KSP_MatMult(ksp,Amat,Z,S);
180:       VecXDot(Z,S,&delta);
181:     }
182:     VecXDot(Z,R,&beta);         /*  beta <- z'*r       */
183:     if (PetscIsInfOrNanScalar(beta)) {
184:       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
185:       else {
186:         ksp->reason = KSP_DIVERGED_NANORINF;
187:         return(0);
188:       }
189:     }
190:   }

192:   i = 0;
193:   do {
194:     ksp->its = i+1;
195:     if (beta == 0.0) {
196:       ksp->reason = KSP_CONVERGED_ATOL;
197:       PetscInfo(ksp,"converged due to beta = 0\n");
198:       break;
199: #if !defined(PETSC_USE_COMPLEX)
200:     } else if ((i > 0) && (beta*betaold < 0.0)) {
201:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
202:       PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
203:       break;
204: #endif
205:     }
206:     if (!i) {
207:       VecCopy(Z,P);         /*     p <- z          */
208:       b    = 0.0;
209:     } else {
210:       b = beta/betaold;
211:       if (eigs) {
212:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
213:         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
214:       }
215:       VecAYPX(P,b,Z);    /*     p <- z + b* p   */
216:     }
217:     dpiold = dpi;
218:     if (!cg->singlereduction || !i) {
219:       KSP_MatMult(ksp,Amat,P,W);          /*     w <- Ap         */
220:       VecXDot(P,W,&dpi);                  /*     dpi <- p'w     */
221:     } else {
222:       VecAYPX(W,beta/betaold,S);                  /*     w <- Ap         */
223:       dpi  = delta - beta*beta*dpiold/(betaold*betaold);             /*     dpi <- p'w     */
224:     }
225:     betaold = beta;
226:     if (PetscIsInfOrNanScalar(dpi)) {
227:       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
228:       else {
229:         ksp->reason = KSP_DIVERGED_NANORINF;
230:         return(0);
231:       }
232:     }

234:     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*dpiold) <= 0.0))) {
235:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
236:       PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
237:       break;
238:     }
239:     a = beta/dpi;                                 /*     a = beta/p'w   */
240:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
241:     VecAXPY(X,a,P);          /*     x <- x + ap     */
242:     VecAXPY(R,-a,W);                      /*     r <- r - aw    */
243:     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
244:       KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
245:       if (cg->singlereduction) {
246:         KSP_MatMult(ksp,Amat,Z,S);
247:       }
248:       VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z       */
249:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
250:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r       */
251:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
252:       KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
253:       if (cg->singlereduction) {
254:         PetscScalar tmp[2];
255:         Vec         vecs[2];
256:         vecs[0] = S; vecs[1] = R;
257:         KSP_MatMult(ksp,Amat,Z,S);
258:         VecMDot(Z,2,vecs,tmp);
259:         delta = tmp[0]; beta = tmp[1];
260:       } else {
261:         VecXDot(Z,R,&beta);     /*  beta <- r'*z       */
262:       }
263:       if (PetscIsInfOrNanScalar(beta)) {
264:         if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
265:         else {
266:           ksp->reason = KSP_DIVERGED_NANORINF;
267:           return(0);
268:         }
269:       }
270:       dp = PetscSqrtReal(PetscAbsScalar(beta));
271:     } else {
272:       dp = 0.0;
273:     }
274:     ksp->rnorm = dp;
275:     KSPLogResidualHistory(ksp,dp);
276:     KSPMonitor(ksp,i+1,dp);
277:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
278:     if (ksp->reason) break;

280:     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) {
281:       KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
282:       if (cg->singlereduction) {
283:         KSP_MatMult(ksp,Amat,Z,S);
284:       }
285:     }
286:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
287:       if (cg->singlereduction) {
288:         PetscScalar tmp[2];
289:         Vec         vecs[2];
290:         vecs[0] = S; vecs[1] = R;
291:         VecMDot(Z,2,vecs,tmp);
292:         delta = tmp[0]; beta = tmp[1];
293:       } else {
294:         VecXDot(Z,R,&beta);        /*  beta <- z'*r       */
295:       }
296:       if (PetscIsInfOrNanScalar(beta)) {
297:         if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
298:         else {
299:           ksp->reason = KSP_DIVERGED_NANORINF;
300:           return(0);
301:         }
302:       }
303:     }

305:     i++;
306:   } while (i<ksp->max_it);
307:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
308:   return(0);
309: }

313: PetscErrorCode KSPDestroy_CG(KSP ksp)
314: {
315:   KSP_CG         *cg = (KSP_CG*)ksp->data;

319:   /* free space used for singular value calculations */
320:   if (ksp->calc_sings) {
321:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
322:   }
323:   KSPDestroyDefault(ksp);
324:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);
325:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);
326:   return(0);
327: }

329: /*
330:      KSPView_CG - Prints information about the current Krylov method being used

332:       Currently this only prints information to a file (or stdout) about the
333:       symmetry of the problem. If your Krylov method has special options or
334:       flags that information should be printed here.

336: */
339: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
340: {
341: #if defined(PETSC_USE_COMPLEX)
342:   KSP_CG         *cg = (KSP_CG*)ksp->data;
344:   PetscBool      iascii;

347:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
348:   if (iascii) {
349:     PetscViewerASCIIPrintf(viewer,"  CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
350:   }
351: #endif
352:   return(0);
353: }

355: /*
356:     KSPSetFromOptions_CG - Checks the options database for options related to the
357:                            conjugate gradient method.
358: */
361: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
362: {
364:   KSP_CG         *cg = (KSP_CG*)ksp->data;

367:   PetscOptionsHead("KSP CG and CGNE options");
368: #if defined(PETSC_USE_COMPLEX)
369:   PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
370:                           (PetscEnum*)&cg->type,NULL);
371: #endif
372:   PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPI_Allreduce()",
373:                           "KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,NULL);
374:   PetscOptionsTail();
375:   return(0);
376: }

378: /*
379:     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
380:                       This routine is registered below in KSPCreate_CG() and called from the
381:                       routine KSPCGSetType() (see the file cgtype.c).
382: */
385: static PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
386: {
387:   KSP_CG *cg = (KSP_CG*)ksp->data;

390:   cg->type = type;
391:   return(0);
392: }

396: static PetscErrorCode  KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg)
397: {
398:   KSP_CG *cg = (KSP_CG*)ksp->data;

401:   cg->singlereduction = flg;
402:   return(0);
403: }

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

409:     It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
410: */
411: /*MC
412:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

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

419:    Level: beginner

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

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

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

430:    References:
431:    Methods of Conjugate Gradients for Solving Linear Systems, Magnus R. Hestenes and Eduard Stiefel,
432:    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
433:    pp. 409--436.

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

438: M*/
441: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
442: {
444:   KSP_CG         *cg;

447:   PetscNewLog(ksp,KSP_CG,&cg);
448: #if !defined(PETSC_USE_COMPLEX)
449:   cg->type = KSP_CG_SYMMETRIC;
450: #else
451:   cg->type = KSP_CG_HERMITIAN;
452: #endif
453:   ksp->data = (void*)cg;

455:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
456:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
457:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
458:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

460:   /*
461:        Sets the functions that are associated with this data structure
462:        (in C++ this is the same as defining virtual functions)
463:   */
464:   ksp->ops->setup          = KSPSetUp_CG;
465:   ksp->ops->solve          = KSPSolve_CG;
466:   ksp->ops->destroy        = KSPDestroy_CG;
467:   ksp->ops->view           = KSPView_CG;
468:   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
469:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
470:   ksp->ops->buildresidual  = KSPBuildResidualDefault;

472:   /*
473:       Attach the function KSPCGSetType_CG() to this object. The routine
474:       KSPCGSetType() checks for this attached function and calls it if it finds
475:       it. (Sort of like a dynamic member function that can be added at run time
476:   */
477:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CG);
478:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",KSPCGUseSingleReduction_CG);
479:   return(0);
480: }