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