xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision fdf6c4e30aafdbc795e4f76379caa977fd5cdf5a)
1 
2 /*  --------------------------------------------------------------------
3 
4      This file implements a Jacobi preconditioner in PETSc as part of PC.
5      You can use this as a starting point for implementing your own
6      preconditioner that is not provided with PETSc. (You might also consider
7      just using PCSHELL)
8 
9      The following basic routines are required for each preconditioner.
10           PCCreate_XXX()          - Creates a preconditioner context
11           PCSetFromOptions_XXX()  - Sets runtime options
12           PCApply_XXX()           - Applies the preconditioner
13           PCDestroy_XXX()         - Destroys the preconditioner context
14      where the suffix "_XXX" denotes a particular implementation, in
15      this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi).
16      These routines are actually called via the common user interface
17      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
18      so the application code interface remains identical for all
19      preconditioners.
20 
21      Another key routine is:
22           PCSetUp_XXX()           - Prepares for the use of a preconditioner
23      by setting data structures and options.   The interface routine PCSetUp()
24      is not usually called directly by the user, but instead is called by
25      PCApply() if necessary.
26 
27      Additional basic routines are:
28           PCView_XXX()            - Prints details of runtime options that
29                                     have actually been used.
30      These are called by application codes via the interface routines
31      PCView().
32 
33      The various types of solvers (preconditioners, Krylov subspace methods,
34      nonlinear solvers, timesteppers) are all organized similarly, so the
35      above description applies to these categories also.  One exception is
36      that the analogues of PCApply() for these components are KSPSolve(),
37      SNESSolve(), and TSSolve().
38 
39      Additional optional functionality unique to preconditioners is left and
40      right symmetric preconditioner application via PCApplySymmetricLeft()
41      and PCApplySymmetricRight().  The Jacobi implementation is
42      PCApplySymmetricLeftOrRight_Jacobi().
43 
44     -------------------------------------------------------------------- */
45 
46 /*
47    Include files needed for the Jacobi preconditioner:
48      pcimpl.h - private include file intended for use by all preconditioners
49 */
50 
51 #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
52 
53 const char *const PCJacobiTypes[]    = {"DIAGONAL","ROWMAX","ROWSUM","PCJacobiType","PC_JACOBI_",NULL};
54 
55 /*
56    Private context (data structure) for the Jacobi preconditioner.
57 */
58 typedef struct {
59   Vec diag;                      /* vector containing the reciprocals of the diagonal elements 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   PetscBool userowmax;           /* set with PCJacobiSetType() */
64   PetscBool userowsum;
65   PetscBool useabs;              /* use the absolute values of the diagonal entries */
66   PetscBool fixdiag;             /* fix zero diagonal terms */
67 } PC_Jacobi;
68 
69 static PetscErrorCode  PCJacobiSetType_Jacobi(PC pc,PCJacobiType type)
70 {
71   PC_Jacobi *j = (PC_Jacobi*)pc->data;
72 
73   PetscFunctionBegin;
74   j->userowmax = PETSC_FALSE;
75   j->userowsum = PETSC_FALSE;
76   if (type == PC_JACOBI_ROWMAX) {
77     j->userowmax = PETSC_TRUE;
78   } else if (type == PC_JACOBI_ROWSUM) {
79     j->userowsum = PETSC_TRUE;
80   }
81   PetscFunctionReturn(0);
82 }
83 
84 static PetscErrorCode  PCJacobiGetType_Jacobi(PC pc,PCJacobiType *type)
85 {
86   PC_Jacobi *j = (PC_Jacobi*)pc->data;
87 
88   PetscFunctionBegin;
89   if (j->userowmax) {
90     *type = PC_JACOBI_ROWMAX;
91   } else if (j->userowsum) {
92     *type = PC_JACOBI_ROWSUM;
93   } else {
94     *type = PC_JACOBI_DIAGONAL;
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode  PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg)
100 {
101   PC_Jacobi *j = (PC_Jacobi*)pc->data;
102 
103   PetscFunctionBegin;
104   j->useabs = flg;
105   PetscFunctionReturn(0);
106 }
107 
108 static PetscErrorCode  PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool *flg)
109 {
110   PC_Jacobi *j = (PC_Jacobi*)pc->data;
111 
112   PetscFunctionBegin;
113   *flg = j->useabs;
114   PetscFunctionReturn(0);
115 }
116 
117 static PetscErrorCode  PCJacobiSetFixDiagonal_Jacobi(PC pc,PetscBool flg)
118 {
119   PC_Jacobi *j = (PC_Jacobi*)pc->data;
120 
121   PetscFunctionBegin;
122   j->fixdiag = flg;
123   PetscFunctionReturn(0);
124 }
125 
126 static PetscErrorCode  PCJacobiGetFixDiagonal_Jacobi(PC pc,PetscBool *flg)
127 {
128   PC_Jacobi *j = (PC_Jacobi*)pc->data;
129 
130   PetscFunctionBegin;
131   *flg = j->fixdiag;
132   PetscFunctionReturn(0);
133 }
134 
135 /* -------------------------------------------------------------------------- */
136 /*
137    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
138                     by setting data structures and options.
139 
140    Input Parameter:
141 .  pc - the preconditioner context
142 
143    Application Interface Routine: PCSetUp()
144 
145    Notes:
146    The interface routine PCSetUp() is not usually called directly by
147    the user, but instead is called by PCApply() if necessary.
148 */
149 static PetscErrorCode PCSetUp_Jacobi(PC pc)
150 {
151   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
152   Vec            diag,diagsqrt;
153   PetscInt       n,i;
154   PetscScalar    *x;
155   PetscBool      zeroflag = PETSC_FALSE;
156 
157   PetscFunctionBegin;
158   /*
159        For most preconditioners the code would begin here something like
160 
161   if (pc->setupcalled == 0) { allocate space the first time this is ever called
162     PetscCall(MatCreateVecs(pc->mat,&jac->diag));
163     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
164   }
165 
166     But for this preconditioner we want to support use of both the matrix' diagonal
167     elements (for left or right preconditioning) and square root of diagonal elements
168     (for symmetric preconditioning).  Hence we do not allocate space here, since we
169     don't know at this point which will be needed (diag and/or diagsqrt) until the user
170     applies the preconditioner, and we don't want to allocate BOTH unless we need
171     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
172     and PCSetUp_Jacobi_Symmetric(), respectively.
173   */
174 
175   /*
176     Here we set up the preconditioner; that is, we copy the diagonal values from
177     the matrix and put them into a format to make them quick to apply as a preconditioner.
178   */
179   diag     = jac->diag;
180   diagsqrt = jac->diagsqrt;
181 
182   if (diag) {
183     PetscBool isset,isspd;
184 
185     if (jac->userowmax) {
186       PetscCall(MatGetRowMaxAbs(pc->pmat,diag,NULL));
187     } else if (jac->userowsum) {
188       PetscCall(MatGetRowSum(pc->pmat,diag));
189     } else {
190       PetscCall(MatGetDiagonal(pc->pmat,diag));
191     }
192     PetscCall(VecReciprocal(diag));
193     if (jac->useabs) PetscCall(VecAbs(diag));
194     PetscCall(MatIsSPDKnown(pc->pmat,&isset,&isspd));
195     if (jac->fixdiag && (!isset || !isspd)) {
196       PetscCall(VecGetLocalSize(diag,&n));
197       PetscCall(VecGetArray(diag,&x));
198       for (i=0; i<n; i++) {
199         if (x[i] == 0.0) {
200           x[i]     = 1.0;
201           zeroflag = PETSC_TRUE;
202         }
203       }
204       PetscCall(VecRestoreArray(diag,&x));
205     }
206   }
207   if (diagsqrt) {
208     if (jac->userowmax) {
209       PetscCall(MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL));
210     } else if (jac->userowsum) {
211       PetscCall(MatGetRowSum(pc->pmat,diagsqrt));
212     } else {
213       PetscCall(MatGetDiagonal(pc->pmat,diagsqrt));
214     }
215     PetscCall(VecGetLocalSize(diagsqrt,&n));
216     PetscCall(VecGetArray(diagsqrt,&x));
217     for (i=0; i<n; i++) {
218       if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i]));
219       else {
220         x[i]     = 1.0;
221         zeroflag = PETSC_TRUE;
222       }
223     }
224     PetscCall(VecRestoreArray(diagsqrt,&x));
225   }
226   if (zeroflag) {
227     PetscCall(PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n"));
228   }
229   PetscFunctionReturn(0);
230 }
231 /* -------------------------------------------------------------------------- */
232 /*
233    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
234    inverse of the square root of the diagonal entries of the matrix.  This
235    is used for symmetric application of the Jacobi preconditioner.
236 
237    Input Parameter:
238 .  pc - the preconditioner context
239 */
240 static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
241 {
242   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
243 
244   PetscFunctionBegin;
245   PetscCall(MatCreateVecs(pc->pmat,&jac->diagsqrt,NULL));
246   PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt));
247   PetscCall(PCSetUp_Jacobi(pc));
248   PetscFunctionReturn(0);
249 }
250 /* -------------------------------------------------------------------------- */
251 /*
252    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
253    inverse of the diagonal entries of the matrix.  This is used for left of
254    right application of the Jacobi preconditioner.
255 
256    Input Parameter:
257 .  pc - the preconditioner context
258 */
259 static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
260 {
261   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
262 
263   PetscFunctionBegin;
264   PetscCall(MatCreateVecs(pc->pmat,&jac->diag,NULL));
265   PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag));
266   PetscCall(PCSetUp_Jacobi(pc));
267   PetscFunctionReturn(0);
268 }
269 /* -------------------------------------------------------------------------- */
270 /*
271    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
272 
273    Input Parameters:
274 .  pc - the preconditioner context
275 .  x - input vector
276 
277    Output Parameter:
278 .  y - output vector
279 
280    Application Interface Routine: PCApply()
281  */
282 static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y)
283 {
284   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
285 
286   PetscFunctionBegin;
287   if (!jac->diag) {
288     PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
289   }
290   PetscCall(VecPointwiseMult(y,x,jac->diag));
291   PetscFunctionReturn(0);
292 }
293 /* -------------------------------------------------------------------------- */
294 /*
295    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
296    symmetric preconditioner to a vector.
297 
298    Input Parameters:
299 .  pc - the preconditioner context
300 .  x - input vector
301 
302    Output Parameter:
303 .  y - output vector
304 
305    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
306 */
307 static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)
308 {
309   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
310 
311   PetscFunctionBegin;
312   if (!jac->diagsqrt) {
313     PetscCall(PCSetUp_Jacobi_Symmetric(pc));
314   }
315   PetscCall(VecPointwiseMult(y,x,jac->diagsqrt));
316   PetscFunctionReturn(0);
317 }
318 
319 /* -------------------------------------------------------------------------- */
320 static PetscErrorCode PCReset_Jacobi(PC pc)
321 {
322   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
323 
324   PetscFunctionBegin;
325   PetscCall(VecDestroy(&jac->diag));
326   PetscCall(VecDestroy(&jac->diagsqrt));
327   PetscFunctionReturn(0);
328 }
329 
330 /*
331    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
332    that was created with PCCreate_Jacobi().
333 
334    Input Parameter:
335 .  pc - the preconditioner context
336 
337    Application Interface Routine: PCDestroy()
338 */
339 static PetscErrorCode PCDestroy_Jacobi(PC pc)
340 {
341   PetscFunctionBegin;
342   PetscCall(PCReset_Jacobi(pc));
343   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",NULL));
344   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",NULL));
345   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",NULL));
346   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",NULL));
347   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",NULL));
348   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",NULL));
349 
350   /*
351       Free the private data structure that was hanging off the PC
352   */
353   PetscCall(PetscFree(pc->data));
354   PetscFunctionReturn(0);
355 }
356 
357 static PetscErrorCode PCSetFromOptions_Jacobi(PetscOptionItems *PetscOptionsObject,PC pc)
358 {
359   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
360   PetscBool      flg;
361   PCJacobiType   deflt,type;
362 
363   PetscFunctionBegin;
364   PetscCall(PCJacobiGetType(pc,&deflt));
365   PetscOptionsHeadBegin(PetscOptionsObject,"Jacobi options");
366   PetscCall(PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg));
367   if (flg) PetscCall(PCJacobiSetType(pc,type));
368   PetscCall(PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL));
369   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal","Fix null terms on diagonal","PCJacobiSetFixDiagonal",jac->fixdiag,&jac->fixdiag,NULL));
370   PetscOptionsHeadEnd();
371   PetscFunctionReturn(0);
372 }
373 
374 static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
375 {
376   PC_Jacobi     *jac = (PC_Jacobi *) pc->data;
377   PetscBool      iascii;
378 
379   PetscFunctionBegin;
380   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
381   if (iascii) {
382     PCJacobiType      type;
383     PetscBool         useAbs,fixdiag;
384     PetscViewerFormat format;
385 
386     PetscCall(PCJacobiGetType(pc, &type));
387     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
388     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
389     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
390     PetscCall(PetscViewerGetFormat(viewer, &format));
391     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
392       PetscCall(VecView(jac->diag, viewer));
393     }
394   }
395   PetscFunctionReturn(0);
396 }
397 
398 /* -------------------------------------------------------------------------- */
399 /*
400    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
401    and sets this as the private data within the generic preconditioning
402    context, PC, that was created within PCCreate().
403 
404    Input Parameter:
405 .  pc - the preconditioner context
406 
407    Application Interface Routine: PCCreate()
408 */
409 
410 /*MC
411      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
412 
413    Options Database Key:
414 +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
415 .    -pc_jacobi_abs - use the absolute value of the diagonal entry
416 -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
417 
418    Level: beginner
419 
420   Notes:
421     By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric
422     can scale each side of the matrix by the square root of the diagonal entries.
423 
424     Zero entries along the diagonal are replaced with the value 1.0
425 
426     See PCPBJACOBI for fixed-size point block, PCVPBJACOBI for variable-sized point block, and PCBJACOBI for large size blocks
427 
428 .seealso:  `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
429            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`,
430            `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
431            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
432 M*/
433 
434 PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
435 {
436   PC_Jacobi      *jac;
437 
438   PetscFunctionBegin;
439   /*
440      Creates the private data structure for this preconditioner and
441      attach it to the PC object.
442   */
443   PetscCall(PetscNewLog(pc,&jac));
444   pc->data = (void*)jac;
445 
446   /*
447      Initialize the pointers to vectors to ZERO; these will be used to store
448      diagonal entries of the matrix for fast preconditioner application.
449   */
450   jac->diag      = NULL;
451   jac->diagsqrt  = NULL;
452   jac->userowmax = PETSC_FALSE;
453   jac->userowsum = PETSC_FALSE;
454   jac->useabs    = PETSC_FALSE;
455   jac->fixdiag   = PETSC_TRUE;
456 
457   /*
458       Set the pointers for the functions that are provided above.
459       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
460       are called, they will automatically call these functions.  Note we
461       choose not to provide a couple of these functions since they are
462       not needed.
463   */
464   pc->ops->apply               = PCApply_Jacobi;
465   pc->ops->applytranspose      = PCApply_Jacobi;
466   pc->ops->setup               = PCSetUp_Jacobi;
467   pc->ops->reset               = PCReset_Jacobi;
468   pc->ops->destroy             = PCDestroy_Jacobi;
469   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
470   pc->ops->view                = PCView_Jacobi;
471   pc->ops->applyrichardson     = NULL;
472   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
473   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
474 
475   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi));
476   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi));
477   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi));
478   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi));
479   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",PCJacobiSetFixDiagonal_Jacobi));
480   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",PCJacobiGetFixDiagonal_Jacobi));
481   PetscFunctionReturn(0);
482 }
483 
484 /*@
485    PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the
486       absolute values of the diagonal divisors in the preconditioner
487 
488    Logically Collective on PC
489 
490    Input Parameters:
491 +  pc - the preconditioner context
492 -  flg - whether to use absolute values or not
493 
494    Options Database Key:
495 .  -pc_jacobi_abs <bool> - use absolute values
496 
497    Notes:
498     This takes affect at the next construction of the preconditioner
499 
500    Level: intermediate
501 
502 .seealso: `PCJacobiaSetType()`, `PCJacobiGetUseAbs()`
503 
504 @*/
505 PetscErrorCode  PCJacobiSetUseAbs(PC pc,PetscBool flg)
506 {
507   PetscFunctionBegin;
508   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
509   PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));
510   PetscFunctionReturn(0);
511 }
512 
513 /*@
514    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the
515       absolute values of the diagonal divisors in the preconditioner
516 
517    Logically Collective on PC
518 
519    Input Parameter:
520 .  pc - the preconditioner context
521 
522    Output Parameter:
523 .  flg - whether to use absolute values or not
524 
525    Level: intermediate
526 
527 .seealso: `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
528 
529 @*/
530 PetscErrorCode  PCJacobiGetUseAbs(PC pc,PetscBool *flg)
531 {
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
534   PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));
535   PetscFunctionReturn(0);
536 }
537 
538 /*@
539    PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
540 
541    Logically Collective on PC
542 
543    Input Parameters:
544 +  pc - the preconditioner context
545 -  flg - the boolean flag
546 
547    Options Database Key:
548 .  -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
549 
550    Notes:
551     This takes affect at the next construction of the preconditioner
552 
553    Level: intermediate
554 
555 .seealso: `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`
556 
557 @*/
558 PetscErrorCode  PCJacobiSetFixDiagonal(PC pc,PetscBool flg)
559 {
560   PetscFunctionBegin;
561   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
562   PetscTryMethod(pc,"PCJacobiSetFixDiagonal_C",(PC,PetscBool),(pc,flg));
563   PetscFunctionReturn(0);
564 }
565 
566 /*@
567    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner checks for zero diagonal terms
568 
569    Logically Collective on PC
570 
571    Input Parameter:
572 .  pc - the preconditioner context
573 
574    Output Parameter:
575 .  flg - the boolean flag
576 
577    Options Database Key:
578 .  -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
579 
580    Level: intermediate
581 
582 .seealso: `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
583 
584 @*/
585 PetscErrorCode  PCJacobiGetFixDiagonal(PC pc,PetscBool *flg)
586 {
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
589   PetscUseMethod(pc,"PCJacobiGetFixDiagonal_C",(PC,PetscBool*),(pc,flg));
590   PetscFunctionReturn(0);
591 }
592 
593 /*@
594    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
595       of the sum of rows entries for the diagonal preconditioner
596 
597    Logically Collective on PC
598 
599    Input Parameters:
600 +  pc - the preconditioner context
601 -  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
602 
603    Options Database Key:
604 .  -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
605 
606    Level: intermediate
607 
608 .seealso: `PCJacobiaUseAbs()`, `PCJacobiGetType()`
609 @*/
610 PetscErrorCode  PCJacobiSetType(PC pc,PCJacobiType type)
611 {
612   PetscFunctionBegin;
613   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
614   PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));
615   PetscFunctionReturn(0);
616 }
617 
618 /*@
619    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
620 
621    Not Collective on PC
622 
623    Input Parameter:
624 .  pc - the preconditioner context
625 
626    Output Parameter:
627 .  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
628 
629    Level: intermediate
630 
631 .seealso: `PCJacobiaUseAbs()`, `PCJacobiSetType()`
632 @*/
633 PetscErrorCode  PCJacobiGetType(PC pc,PCJacobiType *type)
634 {
635   PetscFunctionBegin;
636   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
637   PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));
638   PetscFunctionReturn(0);
639 }
640