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: */
 55: int KSPSetUp_CG(KSP ksp)
 56: {
 57:   KSP_CG *cgP = (KSP_CG*)ksp->data;
 58:   int    maxit = ksp->max_it,ierr;

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

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

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

 92: /*
 93:        KSPSolve_CG - This routine actually applies the conjugate gradient 
 94:     method

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

100:    Output Parameter:
101: .     its - number of iterations used

103: */
104: int  KSPSolve_CG(KSP ksp,int *its)
105: {
106:   int          ierr,i,maxit,eigs;
107:   PetscScalar  dpi,a = 1.0,beta,betaold = 1.0,b,*e = 0,*d = 0,mone = -1.0,ma;
108:   PetscReal    dp = 0.0;
109:   Vec          X,B,Z,R,P;
110:   KSP_CG       *cg;
111:   Mat          Amat,Pmat;
112:   MatStructure pflag;
113:   PetscTruth   diagonalscale;

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

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

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

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

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

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

215: */
216: int KSPDestroy_CG(KSP ksp)
217: {
218:   KSP_CG *cg = (KSP_CG*)ksp->data;
219:   int    ierr;

222:   /* free space used for singular value calculations */
223:   if (ksp->calc_sings) {
224:     PetscFree(cg->e);
225:     PetscFree(cg->ee);
226:   }

228:   KSPDefaultFreeWork(ksp);
229: 
230:   /* free the context variable */
231:   PetscFree(cg);
232:   return(0);
233: }

235: /*
236:      KSPView_CG - Prints information about the current Krylov method being used

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

242: */
243: int KSPView_CG(KSP ksp,PetscViewer viewer)
244: {
245: #if defined(PETSC_USE_COMPLEX)
246:   KSP_CG     *cg = (KSP_CG *)ksp->data;
247:   int        ierr;
248:   PetscTruth isascii;

251:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
252:   if (isascii) {
253:     if (cg->type == KSP_CG_HERMITIAN) {
254:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, Hermitian systemn");
255:     } else if (cg->type == KSP_CG_SYMMETRIC) {
256:       PetscViewerASCIIPrintf(viewer,"  CG: variant for complex, symmetric systemn");
257:     } else {
258:       PetscViewerASCIIPrintf(viewer,"  CG: unknown variantn");
259:     }
260:   } else {
261:     SETERRQ1(1,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
262:   }
263: #endif
264:   return(0);
265: }

267: /*
268:     KSPSetFromOptions_CG - Checks the options database for options related to the 
269:                            conjugate gradient method.
270: */
271: int KSPSetFromOptions_CG(KSP ksp)
272: {
273: #if defined(PETSC_USE_COMPLEX)
274:   int        ierr;
275:   PetscTruth flg;
276: #endif

279: #if defined(PETSC_USE_COMPLEX)
280:   PetscOptionsHead("KSP CG options");
281:     PetscOptionsLogicalGroupBegin("-ksp_cg_Hermitian","Matrix is Hermitian","KSPCGSetType",&flg);
282:     if (flg) { KSPCGSetType(ksp,KSP_CG_HERMITIAN); }
283:     PetscOptionsLogicalGroupEnd("-ksp_cg_symmetric","Matrix is complex symmetric, not Hermitian","KSPCGSetType",&flg);
284:     if (flg) { KSPCGSetType(ksp,KSP_CG_SYMMETRIC); }
285:   PetscOptionsTail();
286: #endif
287:   return(0);
288: }

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

295:         This must be wrapped in an EXTERN_C_BEGIN to be dynamically linkable in C++
296: */
297: EXTERN_C_BEGIN
298: int KSPCGSetType_CG(KSP ksp,KSPCGType type)
299: {
300:   KSP_CG *cg;

303:   cg = (KSP_CG *)ksp->data;
304:   cg->type = type;
305:   return(0);
306: }
307: EXTERN_C_END

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

313:     It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
314: */
315: EXTERN_C_BEGIN
316: int KSPCreate_CG(KSP ksp)
317: {
318:   int    ierr;
319:   KSP_CG *cg;

322:   PetscNew(KSP_CG,&cg);
323:   PetscMemzero(cg,sizeof(KSP_CG));
324:   PetscLogObjectMemory(ksp,sizeof(KSP_CG));
325: #if !defined(PETSC_USE_COMPLEX)
326:   cg->type                       = KSP_CG_SYMMETRIC;
327: #else
328:   cg->type                       = KSP_CG_HERMITIAN;
329: #endif
330:   ksp->data                      = (void*)cg;
331:   ksp->pc_side                   = PC_LEFT;

333:   /*
334:        Sets the functions that are associated with this data structure 
335:        (in C++ this is the same as defining virtual functions)
336:   */
337:   ksp->ops->setup                = KSPSetUp_CG;
338:   ksp->ops->solve                = KSPSolve_CG;
339:   ksp->ops->destroy              = KSPDestroy_CG;
340:   ksp->ops->view                 = KSPView_CG;
341:   ksp->ops->setfromoptions       = KSPSetFromOptions_CG;
342:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
343:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

345:   /*
346:       Attach the function KSPCGSetType_CG() to this object. The routine 
347:       KSPCGSetType() checks for this attached function and calls it if it finds
348:       it. (Sort of like a dynamic member function that can be added at run time
349:   */
350:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG",
351:                                      KSPCGSetType_CG);
352:   return(0);
353: }
354: EXTERN_C_END