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_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:        cgctx.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/cgctx.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: */
 56: PetscErrorCode KSPSetUp_CG(KSP ksp)
 57: {
 58:   KSP_CG         *cgP = (KSP_CG*)ksp->data;
 60:   PetscInt        maxit = ksp->max_it;

 63:   /* 
 64:        This implementation of CG only handles left preconditioning
 65:      so generate an error otherwise.
 66:   */
 67:   if (ksp->pc_side == PC_RIGHT) {
 68:     SETERRQ(2,"No right preconditioning for KSPCG");
 69:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 70:     SETERRQ(2,"No symmetric preconditioning for KSPCG");
 71:   }

 73:   /* get work vectors needed by CG */
 74:   KSPDefaultGetWork(ksp,3);

 76:   /*
 77:      If user requested computations of eigenvalues then allocate work
 78:      work space needed
 79:   */
 80:   if (ksp->calc_sings) {
 81:     /* get space to store tridiagonal matrix for Lanczos */
 82:     PetscMalloc(2*(maxit+1)*sizeof(PetscScalar),&cgP->e);
 83:     PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
 84:     cgP->d                         = cgP->e + maxit + 1;
 85:     PetscMalloc(2*(maxit+1)*sizeof(PetscReal),&cgP->ee);
 86:     PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
 87:     cgP->dd                        = cgP->ee + maxit + 1;
 88:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 89:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 90:   }
 91:   return(0);
 92: }

 94: /*
 95:        KSPSolve_CG - This routine actually applies the conjugate gradient 
 96:     method

 98:    Input Parameter:
 99: .     ksp - the Krylov space object that was set to use conjugate gradient, by, for 
100:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);

102:    Output Parameter:
103: .     its - number of iterations used

105: */
108: PetscErrorCode  KSPSolve_CG(KSP ksp)
109: {
111:   PetscInt       i,stored_max_it,eigs;
112:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b,*e = 0,*d = 0,mone = -1.0,ma;
113:   PetscReal      dp = 0.0;
114:   Vec            X,B,Z,R,P;
115:   KSP_CG         *cg;
116:   Mat            Amat,Pmat;
117:   MatStructure   pflag;
118:   PetscTruth     diagonalscale;

121:   PCDiagonalScale(ksp->pc,&diagonalscale);
122:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);

124:   cg            = (KSP_CG*)ksp->data;
125:   eigs          = ksp->calc_sings;
126:   stored_max_it = ksp->max_it;
127:   X             = ksp->vec_sol;
128:   B             = ksp->vec_rhs;
129:   R             = ksp->work[0];
130:   Z             = ksp->work[1];
131:   P             = ksp->work[2];

133: #if !defined(PETSC_USE_COMPLEX)
134: #define VecXDot(x,y,a) VecDot(x,y,a)
135: #else
136: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
137: #endif

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

142:   ksp->its = 0;
143:   if (!ksp->guess_zero) {
144:     KSP_MatMult(ksp,Amat,X,R);             /*   r <- b - Ax       */
145:     VecAYPX(&mone,B,R);
146:   } else {
147:     VecCopy(B,R);                         /*     r <- b (x is 0) */
148:   }
149:   KSP_PCApply(ksp,R,Z);                   /*     z <- Br         */
150:   VecXDot(Z,R,&beta);
151:   if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
152:     VecNorm(Z,NORM_2,&dp);                /*    dp <- z'*z = e'*A'*B'*B*A'*e'   */
153:   } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
154:     VecNorm(R,NORM_2,&dp);                /*    dp <- r'*r = e'*A'*A*e         */
155:   } else if (ksp->normtype == KSP_NATURAL_NORM) {
156:     dp = sqrt(PetscAbsScalar(beta));                           /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
157:   } else dp = 0.0;
158:   KSPLogResidualHistory(ksp,dp);
159:   KSPMonitor(ksp,0,dp);                              /* call any registered monitor routines */
160:   ksp->rnorm = dp;

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

165:   i = 0;
166:   do {
167:      ksp->its = i+1;
168:      VecXDot(Z,R,&beta);     /*     beta <- r'z     */
169:      if (beta == 0.0) {
170:        ksp->reason = KSP_CONVERGED_ATOL;
171:        PetscLogInfo(ksp,"KSPSolve_CG:converged due to beta = 0");
172:        break;
173: #if !defined(PETSC_USE_COMPLEX)
174:      } else if (beta < 0.0) {
175:        ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
176:        PetscLogInfo(ksp,"KSPSolve_CG:diverging due to indefinite preconditioner");
177:        break;
178: #endif
179:      }
180:      if (!i) {
181:        VecCopy(Z,P);         /*     p <- z          */
182:      } else {
183:          b = beta/betaold;
184:          if (eigs) {
185:            if (ksp->max_it != stored_max_it) {
186:              SETERRQ(PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
187:            }
188:            e[i] = sqrt(PetscAbsScalar(b))/a;
189:          }
190:          VecAYPX(&b,Z,P);    /*     p <- z + b* p   */
191:      }
192:      betaold = beta;
193:      KSP_MatMult(ksp,Amat,P,Z);      /*     z <- Kp         */
194:      VecXDot(P,Z,&dpi);      /*     dpi <- z'p      */
195:      a = beta/dpi;                                 /*     a = beta/p'z    */
196:      if (eigs) {
197:        d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
198:      }
199:      VecAXPY(&a,P,X);          /*     x <- x + ap     */
200:      ma = -a; VecAXPY(&ma,Z,R);                      /*     r <- r - az     */
201:      if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
202:        KSP_PCApply(ksp,R,Z);        /*     z <- Br         */
203:        VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z       */
204:      } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
205:        VecNorm(R,NORM_2,&dp);              /*    dp <- r'*r       */
206:      } else if (ksp->normtype == KSP_NATURAL_NORM) {
207:        dp = sqrt(PetscAbsScalar(beta));
208:      } else {
209:        dp = 0.0;
210:      }
211:      ksp->rnorm = dp;
212:      KSPLogResidualHistory(ksp,dp);
213:      KSPMonitor(ksp,i+1,dp);
214:      (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
215:      if (ksp->reason) break;
216:      if (ksp->normtype != KSP_PRECONDITIONED_NORM) {
217:        KSP_PCApply(ksp,R,Z); /* z <- Br  */
218:      }
219:      i++;
220:   } while (i<ksp->max_it);
221:   if (i >= ksp->max_it) {
222:     ksp->reason = KSP_DIVERGED_ITS;
223:   }
224:   return(0);
225: }
226: /*
227:        KSPDestroy_CG - Frees all memory space used by the Krylov method

229: */
232: PetscErrorCode KSPDestroy_CG(KSP ksp)
233: {
234:   KSP_CG         *cg = (KSP_CG*)ksp->data;

238:   /* free space used for singular value calculations */
239:   if (ksp->calc_sings) {
240:     PetscFree(cg->e);
241:     PetscFree(cg->ee);
242:   }

244:   KSPDefaultFreeWork(ksp);
245: 
246:   /* free the context variable */
247:   PetscFree(cg);
248:   return(0);
249: }

251: /*
252:      KSPView_CG - Prints information about the current Krylov method being used

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

258: */
261: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
262: {
263: #if defined(PETSC_USE_COMPLEX)
264:   KSP_CG         *cg = (KSP_CG *)ksp->data;
266:   PetscTruth     iascii;

269:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
270:   if (iascii) {
271:     if (cg->type == KSP_CG_HERMITIAN) {
272:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, Hermitian system\n");
273:     } else if (cg->type == KSP_CG_SYMMETRIC) {
274:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, symmetric system\n");
275:     } else {
276:       PetscViewerASCIIPrintf(viewer,"  CG: unknown variant\n");
277:     }
278:   } else {
279:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
280:   }
281: #endif
282:   return(0);
283: }

285: /*
286:     KSPSetFromOptions_CG - Checks the options database for options related to the 
287:                            conjugate gradient method.
288: */
291: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
292: {
293: #if defined(PETSC_USE_COMPLEX)
295:   PetscTruth     flg;
296: #endif

299: #if defined(PETSC_USE_COMPLEX)
300:   PetscOptionsHead("KSP CG options");
301:     PetscOptionsLogicalGroupBegin("-ksp_cg_Hermitian","Matrix is Hermitian","KSPCGSetType",&flg);
302:     if (flg) { KSPCGSetType(ksp,KSP_CG_HERMITIAN); }
303:     PetscOptionsLogicalGroupEnd("-ksp_cg_symmetric","Matrix is complex symmetric, not Hermitian","KSPCGSetType",&flg);
304:     if (flg) { KSPCGSetType(ksp,KSP_CG_SYMMETRIC); }
305:   PetscOptionsTail();
306: #endif
307:   return(0);
308: }

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

316: */
320: PetscErrorCode KSPCGSetType_CG(KSP ksp,KSPCGType type)
321: {
322:   KSP_CG *cg;

325:   cg = (KSP_CG *)ksp->data;
326:   cg->type = type;
327:   return(0);
328: }

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

336: */
337: /*MC
338:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

340:    Options Database Keys:
341: +   -ksp_cg_Hermitian - (for complex matrices only) indicates the matrix is Hermitian
342: -   -ksp_cg_symmetric - (for complex matrices only) indicates the matrix is symmetric

344:    Level: beginner

346:    Notes: The PCG method requires both the matrix and preconditioner to 
347:           be symmetric positive (semi) definite

349: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
350:            KSPCGSetType()

352: M*/
356: PetscErrorCode KSPCreate_CG(KSP ksp)
357: {
359:   KSP_CG         *cg;

362:   PetscNew(KSP_CG,&cg);
363:   PetscMemzero(cg,sizeof(KSP_CG));
364:   PetscLogObjectMemory(ksp,sizeof(KSP_CG));
365: #if !defined(PETSC_USE_COMPLEX)
366:   cg->type                       = KSP_CG_SYMMETRIC;
367: #else
368:   cg->type                       = KSP_CG_HERMITIAN;
369: #endif
370:   ksp->data                      = (void*)cg;
371:   ksp->pc_side                   = PC_LEFT;

373:   /*
374:        Sets the functions that are associated with this data structure 
375:        (in C++ this is the same as defining virtual functions)
376:   */
377:   ksp->ops->setup                = KSPSetUp_CG;
378:   ksp->ops->solve                = KSPSolve_CG;
379:   ksp->ops->destroy              = KSPDestroy_CG;
380:   ksp->ops->view                 = KSPView_CG;
381:   ksp->ops->setfromoptions       = KSPSetFromOptions_CG;
382:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
383:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

385:   /*
386:       Attach the function KSPCGSetType_CG() to this object. The routine 
387:       KSPCGSetType() checks for this attached function and calls it if it finds
388:       it. (Sort of like a dynamic member function that can be added at run time
389:   */
390:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG",
391:                                      KSPCGSetType_CG);
392:   return(0);
393: }