Actual source code: cg.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     This file implements the conjugate gradient method in PETSc as part of
  5:     KSP. You can use this as a starting point for implementing your own 
  6:     Krylov method that is not provided with PETSc.

  8:     The following basic routines are required for each Krylov method.
  9:         KSPCreate_XXX()          - Creates the Krylov context
 10:         KSPSetFromOptions_XXX()  - Sets runtime options
 11:         KSPSolve_XXX()           - Runs the Krylov method
 12:         KSPDestroy_XXX()         - Destroys the Krylov context, freeing all 
 13:                                    memory it needed
 14:     Here the "_XXX" denotes a particular implementation, in this case 
 15:     we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines are 
 16:     are actually called vai the common user interface routines
 17:     KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the
 18:     application code interface remains identical for all preconditioners.

 20:     Other basic routines for the KSP objects include
 21:         KSPSetUp_XXX()
 22:         KSPView_XXX()             - Prints details of solver being used.

 24:     Detailed notes:                         
 25:     By default, this code implements the CG (Conjugate Gradient) method,
 26:     which is valid for real symmetric (and complex Hermitian) positive
 27:     definite matrices. Note that for the complex Hermitian case, the
 28:     VecDot() arguments within the code MUST remain in the order given
 29:     for correct computation of inner products.

 31:     Reference: Hestenes and Steifel, 1952.

 33:     By switching to the indefinite vector inner product, VecTDot(), the
 34:     same code is used for the complex symmetric case as well.  The user
 35:     must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option 
 36:     -ksp_cg_type symmetric to invoke this variant for the complex case.
 37:     Note, however, that the complex symmetric code is NOT valid for
 38:     all such matrices ... and thus we don't recommend using this method.
 39: */
 40: /*
 41:        cgctx.h defines the simple data structured used to store information
 42:     related to the type of matrix (e.g. complex symmetric) being solved and
 43:     data used during the optional Lanczo process used to compute eigenvalues
 44: */
 45:  #include src/ksp/ksp/impls/cg/cgctx.h
 46: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
 47: EXTERN PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);

 49: /*
 50:      KSPSetUp_CG - Sets up the workspace needed by the CG method. 

 52:       This is called once, usually automatically by KSPSolve() or KSPSetUp()
 53:      but can be called directly by KSPSetUp()
 54: */
 57: PetscErrorCode KSPSetUp_CG(KSP ksp)
 58: {
 59:   KSP_CG         *cgP = (KSP_CG*)ksp->data;
 61:   PetscInt        maxit = ksp->max_it;

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

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

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

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

 99:    Input Parameter:
100: .     ksp - the Krylov space object that was set to use conjugate gradient, by, for 
101:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
102: */
105: PetscErrorCode  KSPSolve_CG(KSP ksp)
106: {
108:   PetscInt       i,stored_max_it,eigs;
109:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0;
110:   PetscReal      dp = 0.0;
111:   Vec            X,B,Z,R,P;
112:   KSP_CG         *cg;
113:   Mat            Amat,Pmat;
114:   MatStructure   pflag;
115:   PetscTruth     diagonalscale;

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

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

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

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

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

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

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

232: */
235: PetscErrorCode KSPDestroy_CG(KSP ksp)
236: {
237:   KSP_CG         *cg = (KSP_CG*)ksp->data;

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

247:   KSPDefaultFreeWork(ksp);
248: 
249:   /* free the context variable */
250:   PetscFree(cg);
251:   return(0);
252: }

254: /*
255:      KSPView_CG - Prints information about the current Krylov method being used

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

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

272:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
273:   if (iascii) {
274:     PetscViewerASCIIPrintf(viewer,"  CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
275:   } else {
276:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
277:   }
278: #endif
279:   return(0);
280: }

282: /*
283:     KSPSetFromOptions_CG - Checks the options database for options related to the 
284:                            conjugate gradient method.
285: */
288: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
289: {
290: #if defined(PETSC_USE_COMPLEX)
292:   KSP_CG         *cg = (KSP_CG *)ksp->data;
293: #endif

296: #if defined(PETSC_USE_COMPLEX)
297:   PetscOptionsHead("KSP CG and CGNE options");
298:   PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
299:                           (PetscEnum*)&cg->type,PETSC_NULL);
300:   PetscOptionsTail();
301: #endif
302:   return(0);
303: }

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

311: */
315: PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
316: {
317:   KSP_CG *cg;

320:   cg = (KSP_CG *)ksp->data;
321:   cg->type = type;
322:   return(0);
323: }

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

331: */
332: /*MC
333:      KSPCG - The preconditioned conjugate gradient (PCG) iterative method

335:    Options Database Keys:
336: +   -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian
337: -   -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric

339:    Level: beginner

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

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

347: M*/
351: PetscErrorCode  KSPCreate_CG(KSP ksp)
352: {
354:   KSP_CG         *cg;

357:   PetscNew(KSP_CG,&cg);
358:   PetscLogObjectMemory(ksp,sizeof(KSP_CG));
359: #if !defined(PETSC_USE_COMPLEX)
360:   cg->type                       = KSP_CG_SYMMETRIC;
361: #else
362:   cg->type                       = KSP_CG_HERMITIAN;
363: #endif
364:   ksp->data                      = (void*)cg;
365:   ksp->pc_side                   = PC_LEFT;

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

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