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