Actual source code: cg.c

petsc-3.6.4 2016-04-12
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 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,&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:        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:   PetscBool      diagonalscale;

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

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

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

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

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

137:   switch (ksp->normtype) {
138:   case KSP_NORM_PRECONDITIONED:
139:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
140:     VecNorm(Z,NORM_2,&dp);                /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
141:     break;
142:   case KSP_NORM_UNPRECONDITIONED:
143:     VecNorm(R,NORM_2,&dp);                /*    dp <- r'*r = e'*A'*A*e            */
144:     break;
145:   case KSP_NORM_NATURAL:
146:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
147:     if (cg->singlereduction) {
148:       KSP_MatMult(ksp,Amat,Z,S);
149:       VecXDot(Z,S,&delta);
150:     }
151:     VecXDot(Z,R,&beta);                     /*  beta <- z'*r       */
152:     KSPCheckDot(ksp,beta);
153:     dp = PetscSqrtReal(PetscAbsScalar(beta));                           /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
154:     break;
155:   case KSP_NORM_NONE:
156:     dp = 0.0;
157:     break;
158:   default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
159:   }
160:   KSPLogResidualHistory(ksp,dp);
161:   KSPMonitor(ksp,0,dp);
162:   ksp->rnorm = dp;

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

167:   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
168:     KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
169:   }
170:   if (ksp->normtype != KSP_NORM_NATURAL) {
171:     if (cg->singlereduction) {
172:       KSP_MatMult(ksp,Amat,Z,S);
173:       VecXDot(Z,S,&delta);
174:     }
175:     VecXDot(Z,R,&beta);         /*  beta <- z'*r       */
176:     KSPCheckDot(ksp,beta);
177:   }

179:   i = 0;
180:   do {
181:     ksp->its = i+1;
182:     if (beta == 0.0) {
183:       ksp->reason = KSP_CONVERGED_ATOL;
184:       PetscInfo(ksp,"converged due to beta = 0\n");
185:       break;
186: #if !defined(PETSC_USE_COMPLEX)
187:     } else if ((i > 0) && (beta*betaold < 0.0)) {
188:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
189:       PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
190:       break;
191: #endif
192:     }
193:     if (!i) {
194:       VecCopy(Z,P);         /*     p <- z          */
195:       b    = 0.0;
196:     } else {
197:       b = beta/betaold;
198:       if (eigs) {
199:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
200:         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
201:       }
202:       VecAYPX(P,b,Z);    /*     p <- z + b* p   */
203:     }
204:     dpiold = dpi;
205:     if (!cg->singlereduction || !i) {
206:       KSP_MatMult(ksp,Amat,P,W);          /*     w <- Ap         */
207:       VecXDot(P,W,&dpi);                  /*     dpi <- p'w     */
208:     } else {
209:       VecAYPX(W,beta/betaold,S);                  /*     w <- Ap         */
210:       dpi  = delta - beta*beta*dpiold/(betaold*betaold);             /*     dpi <- p'w     */
211:     }
212:     betaold = beta;
213:     KSPCheckDot(ksp,beta);

215:     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*dpiold) <= 0.0))) {
216:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
217:       PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
218:       break;
219:     }
220:     a = beta/dpi;                                 /*     a = beta/p'w   */
221:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
222:     VecAXPY(X,a,P);          /*     x <- x + ap     */
223:     VecAXPY(R,-a,W);                      /*     r <- r - aw    */
224:     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
225:       KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
226:       if (cg->singlereduction) {
227:         KSP_MatMult(ksp,Amat,Z,S);
228:       }
229:       VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z       */
230:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
231:       VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r       */
232:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
233:       KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
234:       if (cg->singlereduction) {
235:         PetscScalar tmp[2];
236:         Vec         vecs[2];
237:         vecs[0] = S; vecs[1] = R;
238:         KSP_MatMult(ksp,Amat,Z,S);
239:         VecMDot(Z,2,vecs,tmp);
240:         delta = tmp[0]; beta = tmp[1];
241:       } else {
242:         VecXDot(Z,R,&beta);     /*  beta <- r'*z       */
243:       }
244:       KSPCheckDot(ksp,beta);
245:       dp = PetscSqrtReal(PetscAbsScalar(beta));
246:     } else {
247:       dp = 0.0;
248:     }
249:     ksp->rnorm = dp;
250:     KSPLogResidualHistory(ksp,dp);
251:     if (eigs) cg->ned = ksp->its;
252:     KSPMonitor(ksp,i+1,dp);
253:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
254:     if (ksp->reason) break;

256:     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) {
257:       KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
258:       if (cg->singlereduction) {
259:         KSP_MatMult(ksp,Amat,Z,S);
260:       }
261:     }
262:     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
263:       if (cg->singlereduction) {
264:         PetscScalar tmp[2];
265:         Vec         vecs[2];
266:         vecs[0] = S; vecs[1] = R;
267:         VecMDot(Z,2,vecs,tmp);
268:         delta = tmp[0]; beta = tmp[1];
269:       } else {
270:         VecXDot(Z,R,&beta);        /*  beta <- z'*r       */
271:       }
272:       KSPCheckDot(ksp,beta);
273:     }

275:     i++;
276:   } while (i<ksp->max_it);
277:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
278:   return(0);
279: }

283: PetscErrorCode KSPDestroy_CG(KSP ksp)
284: {
285:   KSP_CG         *cg = (KSP_CG*)ksp->data;

289:   /* free space used for singular value calculations */
290:   if (ksp->calc_sings) {
291:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
292:   }
293:   KSPDestroyDefault(ksp);
294:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);
295:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);
296:   return(0);
297: }

299: /*
300:      KSPView_CG - Prints information about the current Krylov method being used

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

306: */
309: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
310: {
311: #if defined(PETSC_USE_COMPLEX)
312:   KSP_CG         *cg = (KSP_CG*)ksp->data;
314:   PetscBool      iascii;

317:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
318:   if (iascii) {
319:     PetscViewerASCIIPrintf(viewer,"  CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
320:   }
321: #endif
322:   return(0);
323: }

325: /*
326:     KSPSetFromOptions_CG - Checks the options database for options related to the
327:                            conjugate gradient method.
328: */
331: PetscErrorCode KSPSetFromOptions_CG(PetscOptions *PetscOptionsObject,KSP ksp)
332: {
334:   KSP_CG         *cg = (KSP_CG*)ksp->data;

337:   PetscOptionsHead(PetscOptionsObject,"KSP CG and CGNE options");
338: #if defined(PETSC_USE_COMPLEX)
339:   PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
340:                           (PetscEnum*)&cg->type,NULL);
341: #endif
342:   PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPI_Allreduce()","KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,NULL);
343:   PetscOptionsTail();
344:   return(0);
345: }

347: /*
348:     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
349:                       This routine is registered below in KSPCreate_CG() and called from the
350:                       routine KSPCGSetType() (see the file cgtype.c).
351: */
354: static PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
355: {
356:   KSP_CG *cg = (KSP_CG*)ksp->data;

359:   cg->type = type;
360:   return(0);
361: }

365: static PetscErrorCode  KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg)
366: {
367:   KSP_CG *cg = (KSP_CG*)ksp->data;

370:   cg->singlereduction = flg;
371:   return(0);
372: }

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

378:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
379: */
380: /*MC
381:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

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

388:    Level: beginner

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

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

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

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

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

407: M*/
410: PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
411: {
413:   KSP_CG         *cg;

416:   PetscNewLog(ksp,&cg);
417: #if !defined(PETSC_USE_COMPLEX)
418:   cg->type = KSP_CG_SYMMETRIC;
419: #else
420:   cg->type = KSP_CG_HERMITIAN;
421: #endif
422:   ksp->data = (void*)cg;

424:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
425:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
426:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);

428:   /*
429:        Sets the functions that are associated with this data structure
430:        (in C++ this is the same as defining virtual functions)
431:   */
432:   ksp->ops->setup          = KSPSetUp_CG;
433:   ksp->ops->solve          = KSPSolve_CG;
434:   ksp->ops->destroy        = KSPDestroy_CG;
435:   ksp->ops->view           = KSPView_CG;
436:   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
437:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
438:   ksp->ops->buildresidual  = KSPBuildResidualDefault;

440:   /*
441:       Attach the function KSPCGSetType_CG() to this object. The routine
442:       KSPCGSetType() checks for this attached function and calls it if it finds
443:       it. (Sort of like a dynamic member function that can be added at run time
444:   */
445:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CG);
446:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",KSPCGUseSingleReduction_CG);
447:   return(0);
448: }