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