Actual source code: jacobi.c

  1: #define PETSCKSP_DLL

  3: /*  -------------------------------------------------------------------- 

  5:      This file implements a Jacobi preconditioner in PETSc as part of PC.
  6:      You can use this as a starting point for implementing your own 
  7:      preconditioner that is not provided with PETSc. (You might also consider
  8:      just using PCSHELL)

 10:      The following basic routines are required for each preconditioner.
 11:           PCCreate_XXX()          - Creates a preconditioner context
 12:           PCSetFromOptions_XXX()  - Sets runtime options
 13:           PCApply_XXX()           - Applies the preconditioner
 14:           PCDestroy_XXX()         - Destroys the preconditioner context
 15:      where the suffix "_XXX" denotes a particular implementation, in
 16:      this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi).
 17:      These routines are actually called via the common user interface
 18:      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 
 19:      so the application code interface remains identical for all 
 20:      preconditioners.  

 22:      Another key routine is:
 23:           PCSetUp_XXX()           - Prepares for the use of a preconditioner
 24:      by setting data structures and options.   The interface routine PCSetUp()
 25:      is not usually called directly by the user, but instead is called by
 26:      PCApply() if necessary.

 28:      Additional basic routines are:
 29:           PCView_XXX()            - Prints details of runtime options that
 30:                                     have actually been used.
 31:      These are called by application codes via the interface routines
 32:      PCView().

 34:      The various types of solvers (preconditioners, Krylov subspace methods,
 35:      nonlinear solvers, timesteppers) are all organized similarly, so the
 36:      above description applies to these categories also.  One exception is
 37:      that the analogues of PCApply() for these components are KSPSolve(), 
 38:      SNESSolve(), and TSSolve().

 40:      Additional optional functionality unique to preconditioners is left and
 41:      right symmetric preconditioner application via PCApplySymmetricLeft() 
 42:      and PCApplySymmetricRight().  The Jacobi implementation is 
 43:      PCApplySymmetricLeftOrRight_Jacobi().

 45:     -------------------------------------------------------------------- */

 47: /* 
 48:    Include files needed for the Jacobi preconditioner:
 49:      pcimpl.h - private include file intended for use by all preconditioners 
 50: */

 52:  #include private/pcimpl.h

 54: /* 
 55:    Private context (data structure) for the Jacobi preconditioner.  
 56: */
 57: typedef struct {
 58:   Vec        diag;               /* vector containing the reciprocals of the diagonal elements
 59:                                     of the preconditioner matrix */
 60:   Vec        diagsqrt;           /* vector containing the reciprocals of the square roots of
 61:                                     the diagonal elements of the preconditioner matrix (used 
 62:                                     only for symmetric preconditioner application) */
 63:   PetscTruth userowmax;
 64:   PetscTruth userowsum;
 65:   PetscTruth useabs;             /* use the absolute values of the diagonal entries */
 66: } PC_Jacobi;

 71: PetscErrorCode  PCJacobiSetUseRowMax_Jacobi(PC pc)
 72: {
 73:   PC_Jacobi *j;

 76:   j            = (PC_Jacobi*)pc->data;
 77:   j->userowmax = PETSC_TRUE;
 78:   return(0);
 79: }

 85: PetscErrorCode  PCJacobiSetUseRowSum_Jacobi(PC pc)
 86: {
 87:   PC_Jacobi *j;

 90:   j            = (PC_Jacobi*)pc->data;
 91:   j->userowsum = PETSC_TRUE;
 92:   return(0);
 93: }

 99: PetscErrorCode  PCJacobiSetUseAbs_Jacobi(PC pc)
100: {
101:   PC_Jacobi *j;

104:   j         = (PC_Jacobi*)pc->data;
105:   j->useabs = PETSC_TRUE;
106:   return(0);
107: }

110: /* -------------------------------------------------------------------------- */
111: /*
112:    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
113:                     by setting data structures and options.   

115:    Input Parameter:
116: .  pc - the preconditioner context

118:    Application Interface Routine: PCSetUp()

120:    Notes:
121:    The interface routine PCSetUp() is not usually called directly by
122:    the user, but instead is called by PCApply() if necessary.
123: */
126: static PetscErrorCode PCSetUp_Jacobi(PC pc)
127: {
128:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
129:   Vec            diag,diagsqrt;
131:   PetscInt       n,i;
132:   PetscScalar    *x;
133:   PetscTruth     zeroflag = PETSC_FALSE;

136:   /*
137:        For most preconditioners the code would begin here something like

139:   if (pc->setupcalled == 0) { allocate space the first time this is ever called
140:     MatGetVecs(pc->mat,&jac->diag);
141:     PetscLogObjectParent(pc,jac->diag);
142:   }

144:     But for this preconditioner we want to support use of both the matrix' diagonal
145:     elements (for left or right preconditioning) and square root of diagonal elements
146:     (for symmetric preconditioning).  Hence we do not allocate space here, since we
147:     don't know at this point which will be needed (diag and/or diagsqrt) until the user
148:     applies the preconditioner, and we don't want to allocate BOTH unless we need
149:     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
150:     and PCSetUp_Jacobi_Symmetric(), respectively.
151:   */

153:   /*
154:     Here we set up the preconditioner; that is, we copy the diagonal values from
155:     the matrix and put them into a format to make them quick to apply as a preconditioner.
156:   */
157:   diag     = jac->diag;
158:   diagsqrt = jac->diagsqrt;

160:   if (diag) {
161:     if (jac->userowmax) {
162:       MatGetRowMaxAbs(pc->pmat,diag,PETSC_NULL);
163:     } else if (jac->userowsum) {
164:       MatGetRowSum(pc->pmat,diag);
165:     } else {
166:       MatGetDiagonal(pc->pmat,diag);
167:     }
168:     VecReciprocal(diag);
169:     VecGetLocalSize(diag,&n);
170:     VecGetArray(diag,&x);
171:     if (jac->useabs) {
172:       for (i=0; i<n; i++) {
173:         x[i]     = PetscAbsScalar(x[i]);
174:       }
175:     }
176:     for (i=0; i<n; i++) {
177:       if (x[i] == 0.0) {
178:         x[i]     = 1.0;
179:         zeroflag = PETSC_TRUE;
180:       }
181:     }
182:     VecRestoreArray(diag,&x);
183:   }
184:   if (diagsqrt) {
185:     if (jac->userowmax) {
186:       MatGetRowMaxAbs(pc->pmat,diagsqrt,PETSC_NULL);
187:     } else if (jac->userowsum) {
188:       MatGetRowSum(pc->pmat,diagsqrt);
189:     } else {
190:       MatGetDiagonal(pc->pmat,diagsqrt);
191:     }
192:     VecGetLocalSize(diagsqrt,&n);
193:     VecGetArray(diagsqrt,&x);
194:     for (i=0; i<n; i++) {
195:       if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i]));
196:       else {
197:         x[i]     = 1.0;
198:         zeroflag = PETSC_TRUE;
199:       }
200:     }
201:     VecRestoreArray(diagsqrt,&x);
202:   }
203:   if (zeroflag) {
204:     PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");
205:   }
206:   return(0);
207: }
208: /* -------------------------------------------------------------------------- */
209: /*
210:    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
211:    inverse of the square root of the diagonal entries of the matrix.  This
212:    is used for symmetric application of the Jacobi preconditioner.

214:    Input Parameter:
215: .  pc - the preconditioner context
216: */
219: static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
220: {
222:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

225:   MatGetVecs(pc->pmat,&jac->diagsqrt,0);
226:   PetscLogObjectParent(pc,jac->diagsqrt);
227:   PCSetUp_Jacobi(pc);
228:   return(0);
229: }
230: /* -------------------------------------------------------------------------- */
231: /*
232:    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
233:    inverse of the diagonal entries of the matrix.  This is used for left of
234:    right application of the Jacobi preconditioner.

236:    Input Parameter:
237: .  pc - the preconditioner context
238: */
241: static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
242: {
244:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

247:   MatGetVecs(pc->pmat,&jac->diag,0);
248:   PetscLogObjectParent(pc,jac->diag);
249:   PCSetUp_Jacobi(pc);
250:   return(0);
251: }
252: /* -------------------------------------------------------------------------- */
253: /*
254:    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.

256:    Input Parameters:
257: .  pc - the preconditioner context
258: .  x - input vector

260:    Output Parameter:
261: .  y - output vector

263:    Application Interface Routine: PCApply()
264:  */
267: static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y)
268: {
269:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

273:   if (!jac->diag) {
274:     PCSetUp_Jacobi_NonSymmetric(pc);
275:   }
276:   VecPointwiseMult(y,x,jac->diag);
277:   return(0);
278: }
279: /* -------------------------------------------------------------------------- */
280: /*
281:    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
282:    symmetric preconditioner to a vector.

284:    Input Parameters:
285: .  pc - the preconditioner context
286: .  x - input vector

288:    Output Parameter:
289: .  y - output vector

291:    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
292: */
295: static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)
296: {
298:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

301:   if (!jac->diagsqrt) {
302:     PCSetUp_Jacobi_Symmetric(pc);
303:   }
304:   VecPointwiseMult(y,x,jac->diagsqrt);
305:   return(0);
306: }
307: /* -------------------------------------------------------------------------- */
308: /*
309:    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
310:    that was created with PCCreate_Jacobi().

312:    Input Parameter:
313: .  pc - the preconditioner context

315:    Application Interface Routine: PCDestroy()
316: */
319: static PetscErrorCode PCDestroy_Jacobi(PC pc)
320: {
321:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

325:   if (jac->diag)     {VecDestroy(jac->diag);}
326:   if (jac->diagsqrt) {VecDestroy(jac->diagsqrt);}

328:   /*
329:       Free the private data structure that was hanging off the PC
330:   */
331:   PetscFree(jac);
332:   return(0);
333: }

337: static PetscErrorCode PCSetFromOptions_Jacobi(PC pc)
338: {
339:   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;

343:   PetscOptionsHead("Jacobi options");
344:     PetscOptionsTruth("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax,
345:                           &jac->userowmax,PETSC_NULL);
346:     PetscOptionsTruth("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum,
347:                           &jac->userowsum,PETSC_NULL);
348:     PetscOptionsTruth("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs,
349:                           &jac->useabs,PETSC_NULL);
350:   PetscOptionsTail();
351:   return(0);
352: }

354: /* -------------------------------------------------------------------------- */
355: /*
356:    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 
357:    and sets this as the private data within the generic preconditioning 
358:    context, PC, that was created within PCCreate().

360:    Input Parameter:
361: .  pc - the preconditioner context

363:    Application Interface Routine: PCCreate()
364: */

366: /*MC
367:      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)

369:    Options Database Key:
370: +    -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor,
371:                         rather than the diagonal
372: .    -pc_jacobi_rowsum - use the maximum absolute value in each row as the scaling factor,
373:                         rather than the diagonal
374: -    -pc_jacobi_abs - use the absolute value of the diagaonl entry

376:    Level: beginner

378:   Concepts: Jacobi, diagonal scaling, preconditioners

380:   Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 
381:          can scale each side of the matrix by the squareroot of the diagonal entries.

383:          Zero entries along the diagonal are replaced with the value 1.0

385: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
386:            PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs()
387: M*/

392: PetscErrorCode  PCCreate_Jacobi(PC pc)
393: {
394:   PC_Jacobi      *jac;

398:   /*
399:      Creates the private data structure for this preconditioner and
400:      attach it to the PC object.
401:   */
402:   PetscNewLog(pc,PC_Jacobi,&jac);
403:   pc->data  = (void*)jac;

405:   /*
406:      Initialize the pointers to vectors to ZERO; these will be used to store
407:      diagonal entries of the matrix for fast preconditioner application.
408:   */
409:   jac->diag          = 0;
410:   jac->diagsqrt      = 0;
411:   jac->userowmax     = PETSC_FALSE;
412:   jac->userowsum     = PETSC_FALSE;
413:   jac->useabs        = PETSC_FALSE;

415:   /*
416:       Set the pointers for the functions that are provided above.
417:       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
418:       are called, they will automatically call these functions.  Note we
419:       choose not to provide a couple of these functions since they are
420:       not needed.
421:   */
422:   pc->ops->apply               = PCApply_Jacobi;
423:   pc->ops->applytranspose      = PCApply_Jacobi;
424:   pc->ops->setup               = PCSetUp_Jacobi;
425:   pc->ops->destroy             = PCDestroy_Jacobi;
426:   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
427:   pc->ops->view                = 0;
428:   pc->ops->applyrichardson     = 0;
429:   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
430:   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
431:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);
432:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);
433:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);
434:   return(0);
435: }


441: /*@
442:    PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 
443:       absolute value of the diagonal to for the preconditioner

445:    Collective on PC

447:    Input Parameters:
448: .  pc - the preconditioner context


451:    Options Database Key:
452: .  -pc_jacobi_abs

454:    Level: intermediate

456:    Concepts: Jacobi preconditioner

458: .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum()

460: @*/
461: PetscErrorCode  PCJacobiSetUseAbs(PC pc)
462: {
463:   PetscErrorCode ierr,(*f)(PC);

467:   PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",(void (**)(void))&f);
468:   if (f) {
469:     (*f)(pc);
470:   }
471:   return(0);
472: }

476: /*@
477:    PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 
478:       maximum entry in each row as the diagonal preconditioner, instead of
479:       the diagonal entry

481:    Collective on PC

483:    Input Parameters:
484: .  pc - the preconditioner context


487:    Options Database Key:
488: .  -pc_jacobi_rowmax 

490:    Level: intermediate

492:    Concepts: Jacobi preconditioner

494: .seealso: PCJacobiaUseAbs()
495: @*/
496: PetscErrorCode  PCJacobiSetUseRowMax(PC pc)
497: {
498:   PetscErrorCode ierr,(*f)(PC);

502:   PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseRowMax_C",(void (**)(void))&f);
503:   if (f) {
504:     (*f)(pc);
505:   }
506:   return(0);
507: }

511: /*@
512:    PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 
513:       sum of each row as the diagonal preconditioner, instead of
514:       the diagonal entry

516:    Collective on PC

518:    Input Parameters:
519: .  pc - the preconditioner context


522:    Options Database Key:
523: .  -pc_jacobi_rowsum

525:    Level: intermediate

527:    Concepts: Jacobi preconditioner

529: .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum()
530: @*/
531: PetscErrorCode  PCJacobiSetUseRowSum(PC pc)
532: {
533:   PetscErrorCode ierr,(*f)(PC);

537:   PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseRowSum_C",(void (**)(void))&f);
538:   if (f) {
539:     (*f)(pc);
540:   }
541:   return(0);
542: }