Actual source code: cg.c

  1: /*$Id: cg.c,v 1.117 2001/08/07 03:03:50 balay Exp $*/

  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_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/sles/ksp/impls/cg/cgctx.h
 46: EXTERN int KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
 47: EXTERN int KSPComputeEigenvalues_CG(KSP,int,PetscReal *,PetscReal *,int *);

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

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

 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: int  KSPSolve_CG(KSP ksp,int *its)
109: {
110:   int          ierr,i,maxit,eigs;
111:   PetscScalar  dpi,a = 1.0,beta,betaold = 1.0,b,*e = 0,*d = 0,mone = -1.0,ma;
112:   PetscReal    dp = 0.0;
113:   Vec          X,B,Z,R,P;
114:   KSP_CG       *cg;
115:   Mat          Amat,Pmat;
116:   MatStructure pflag;
117:   PetscTruth   diagonalscale;

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

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

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

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

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

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

223: */
226: int KSPDestroy_CG(KSP ksp)
227: {
228:   KSP_CG *cg = (KSP_CG*)ksp->data;
229:   int    ierr;

232:   /* free space used for singular value calculations */
233:   if (ksp->calc_sings) {
234:     PetscFree(cg->e);
235:     PetscFree(cg->ee);
236:   }

238:   KSPDefaultFreeWork(ksp);
239: 
240:   /* free the context variable */
241:   PetscFree(cg);
242:   return(0);
243: }

245: /*
246:      KSPView_CG - Prints information about the current Krylov method being used

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

252: */
255: int KSPView_CG(KSP ksp,PetscViewer viewer)
256: {
257: #if defined(PETSC_USE_COMPLEX)
258:   KSP_CG     *cg = (KSP_CG *)ksp->data;
259:   int        ierr;
260:   PetscTruth isascii;

263:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
264:   if (isascii) {
265:     if (cg->type == KSP_CG_HERMITIAN) {
266:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, Hermitian system\n");
267:     } else if (cg->type == KSP_CG_SYMMETRIC) {
268:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, symmetric system\n");
269:     } else {
270:       PetscViewerASCIIPrintf(viewer,"  CG: unknown variant\n");
271:     }
272:   } else {
273:     SETERRQ1(1,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
274:   }
275: #endif
276:   return(0);
277: }

279: /*
280:     KSPSetFromOptions_CG - Checks the options database for options related to the 
281:                            conjugate gradient method.
282: */
285: int KSPSetFromOptions_CG(KSP ksp)
286: {
287: #if defined(PETSC_USE_COMPLEX)
288:   int        ierr;
289:   PetscTruth flg;
290: #endif

293: #if defined(PETSC_USE_COMPLEX)
294:   PetscOptionsHead("KSP CG options");
295:     PetscOptionsLogicalGroupBegin("-ksp_cg_Hermitian","Matrix is Hermitian","KSPCGSetType",&flg);
296:     if (flg) { KSPCGSetType(ksp,KSP_CG_HERMITIAN); }
297:     PetscOptionsLogicalGroupEnd("-ksp_cg_symmetric","Matrix is complex symmetric, not Hermitian","KSPCGSetType",&flg);
298:     if (flg) { KSPCGSetType(ksp,KSP_CG_SYMMETRIC); }
299:   PetscOptionsTail();
300: #endif
301:   return(0);
302: }

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

309:         This must be wrapped in an EXTERN_C_BEGIN to be dynamically linkable in C++
310: */
311: EXTERN_C_BEGIN
314: int KSPCGSetType_CG(KSP ksp,KSPCGType type)
315: {
316:   KSP_CG *cg;

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

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

329:     It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
330: */
331: EXTERN_C_BEGIN
334: int KSPCreate_CG(KSP ksp)
335: {
336:   int    ierr;
337:   KSP_CG *cg;

340:   PetscNew(KSP_CG,&cg);
341:   PetscMemzero(cg,sizeof(KSP_CG));
342:   PetscLogObjectMemory(ksp,sizeof(KSP_CG));
343: #if !defined(PETSC_USE_COMPLEX)
344:   cg->type                       = KSP_CG_SYMMETRIC;
345: #else
346:   cg->type                       = KSP_CG_HERMITIAN;
347: #endif
348:   ksp->data                      = (void*)cg;
349:   ksp->pc_side                   = PC_LEFT;

351:   /*
352:        Sets the functions that are associated with this data structure 
353:        (in C++ this is the same as defining virtual functions)
354:   */
355:   ksp->ops->setup                = KSPSetUp_CG;
356:   ksp->ops->solve                = KSPSolve_CG;
357:   ksp->ops->destroy              = KSPDestroy_CG;
358:   ksp->ops->view                 = KSPView_CG;
359:   ksp->ops->setfromoptions       = KSPSetFromOptions_CG;
360:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
361:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

363:   /*
364:       Attach the function KSPCGSetType_CG() to this object. The routine 
365:       KSPCGSetType() checks for this attached function and calls it if it finds
366:       it. (Sort of like a dynamic member function that can be added at run time
367:   */
368:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG",
369:                                      KSPCGSetType_CG);
370:   return(0);
371: }
372: EXTERN_C_END