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