Actual source code: cg.c

  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*/

 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:   KSPDefaultGetWork(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)));
 75:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 76:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 77:   }
 78:   return(0);
 79: }

 81: /*
 82:        KSPSolve_CG - This routine actually applies the conjugate gradient  method

 84:    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.

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

105:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
106:   if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,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,&pflag);

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:     if (PetscIsInfOrNanScalar(beta)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
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(((PetscObject)ksp)->comm,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:     if (PetscIsInfOrNanScalar(beta)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
177:   }
178: 
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(((PetscObject)ksp)->comm,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:      if (PetscIsInfOrNanScalar(dpi)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");

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) {
222:        d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
223:      }
224:      VecAXPY(X,a,P);          /*     x <- x + ap     */
225:      VecAXPY(R,-a,W);                      /*     r <- r - aw    */
226:      if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
227:        KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
228:        if (cg->singlereduction) {
229:          KSP_MatMult(ksp,Amat,Z,S);
230:        }
231:        VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z       */
232:      } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
233:        VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r       */
234:      } else if (ksp->normtype == KSP_NORM_NATURAL) {
235:        KSP_PCApply(ksp,R,Z);               /*     z <- Br         */
236:        if (cg->singlereduction) {
237:          PetscScalar tmp[2];
238:          Vec         vecs[2];
239:          vecs[0] = S; vecs[1] = R;
240:          KSP_MatMult(ksp,Amat,Z,S);
241:          /*VecXDot(Z,S,&delta);
242:            VecXDot(Z,R,&beta); */    /*  beta <- r'*z       */
243:          VecMDot(Z,2,vecs,tmp);
244:          delta = tmp[0]; beta = tmp[1];
245:        } else {
246:          VecXDot(Z,R,&beta);     /*  beta <- r'*z       */
247:        }
248:        if (PetscIsInfOrNanScalar(beta)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
249:        dp = PetscSqrtReal(PetscAbsScalar(beta));
250:      } else {
251:        dp = 0.0;
252:      }
253:      ksp->rnorm = dp;
254:      KSPLogResidualHistory(ksp,dp);
255:      KSPMonitor(ksp,i+1,dp);
256:      (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
257:      if (ksp->reason) break;

259:      if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)){
260:        KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
261:        if (cg->singlereduction) {
262:          KSP_MatMult(ksp,Amat,Z,S);
263:        }
264:      }
265:      if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)){
266:        if (cg->singlereduction) {
267:          PetscScalar tmp[2];
268:          Vec         vecs[2];
269:          vecs[0] = S; vecs[1] = R;
270:          /* VecXDot(Z,R,&beta);   */     /*  beta <- z'*r       */
271:          /* VecXDot(Z,S,&delta);*/
272:          VecMDot(Z,2,vecs,tmp);
273:          delta = tmp[0]; beta = tmp[1];
274:        } else {
275:          VecXDot(Z,R,&beta);        /*  beta <- z'*r       */
276:        }
277:        if (PetscIsInfOrNanScalar(beta)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
278:      }

280:      i++;
281:   } while (i<ksp->max_it);
282:   if (i >= ksp->max_it) {
283:     ksp->reason = KSP_DIVERGED_ITS;
284:   }
285:   return(0);
286: }

290: PetscErrorCode KSPDestroy_CG(KSP ksp)
291: {
292:   KSP_CG         *cg = (KSP_CG*)ksp->data;

296:   /* free space used for singular value calculations */
297:   if (ksp->calc_sings) {
298:     PetscFree4(cg->e,cg->d,cg->ee,cg->dd);
299:   }
300:   KSPDefaultDestroy(ksp);
301:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","",PETSC_NULL);
302:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGUseSingleReduction_C","",PETSC_NULL);
303:   return(0);
304: }

306: /*
307:      KSPView_CG - Prints information about the current Krylov method being used

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

313: */
316: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
317: {
318: #if defined(PETSC_USE_COMPLEX)
319:   KSP_CG         *cg = (KSP_CG *)ksp->data;
321:   PetscBool      iascii;

324:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
325:   if (iascii) {
326:     PetscViewerASCIIPrintf(viewer,"  CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
327:   } else {
328:     SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
329:   }
330: #endif
331:   return(0);
332: }

334: /*
335:     KSPSetFromOptions_CG - Checks the options database for options related to the 
336:                            conjugate gradient method.
337: */
340: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
341: {
343:   KSP_CG         *cg = (KSP_CG *)ksp->data;

346:   PetscOptionsHead("KSP CG and CGNE options");
347: #if defined(PETSC_USE_COMPLEX)
348:   PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
349:                           (PetscEnum*)&cg->type,PETSC_NULL);
350: #endif
351:   PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPI_Allreduce()",
352:                            "KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,PETSC_NULL);
353:   PetscOptionsTail();
354:   return(0);
355: }

357: /*
358:     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
359:                       This routine is registered below in KSPCreate_CG() and called from the 
360:                       routine KSPCGSetType() (see the file cgtype.c).

363: */
367: PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
368: {
369:   KSP_CG *cg = (KSP_CG *)ksp->data;

372:   cg->type = type;
373:   return(0);
374: }

380: PetscErrorCode  KSPCGUseSingleReduction_CG(KSP ksp,PetscBool  flg)
381: {
382:   KSP_CG *cg  = (KSP_CG *)ksp->data;

385:   cg->singlereduction = flg;
386:   return(0);
387: }

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

395: */
396: /*MC
397:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

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

404:    Level: beginner

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

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

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

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

420: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
421:            KSPCGSetType(), KSPCGUseSingleReduction()

423: M*/
427: PetscErrorCode  KSPCreate_CG(KSP ksp)
428: {
430:   KSP_CG         *cg;

433:   PetscNewLog(ksp,KSP_CG,&cg);
434: #if !defined(PETSC_USE_COMPLEX)
435:   cg->type                       = KSP_CG_SYMMETRIC;
436: #else
437:   cg->type                       = KSP_CG_HERMITIAN;
438: #endif
439:   ksp->data                      = (void*)cg;

441:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
442:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
443:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
444:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

446:   /*
447:        Sets the functions that are associated with this data structure 
448:        (in C++ this is the same as defining virtual functions)
449:   */
450:   ksp->ops->setup                = KSPSetUp_CG;
451:   ksp->ops->solve                = KSPSolve_CG;
452:   ksp->ops->destroy              = KSPDestroy_CG;
453:   ksp->ops->view                 = KSPView_CG;
454:   ksp->ops->setfromoptions       = KSPSetFromOptions_CG;
455:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
456:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

458:   /*
459:       Attach the function KSPCGSetType_CG() to this object. The routine 
460:       KSPCGSetType() checks for this attached function and calls it if it finds
461:       it. (Sort of like a dynamic member function that can be added at run time
462:   */
463:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG", KSPCGSetType_CG);
464:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGUseSingleReduction_C","KSPCGUseSingleReduction_CG", KSPCGUseSingleReduction_CG);
465:   return(0);
466: }