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