xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
14b9ad928SBarry Smith 
24b9ad928SBarry Smith /*  --------------------------------------------------------------------
34b9ad928SBarry Smith 
45a46d3faSBarry Smith      This file implements a Jacobi preconditioner in PETSc as part of PC.
55a46d3faSBarry Smith      You can use this as a starting point for implementing your own
65a46d3faSBarry Smith      preconditioner that is not provided with PETSc. (You might also consider
75a46d3faSBarry Smith      just using PCSHELL)
84b9ad928SBarry Smith 
94b9ad928SBarry Smith      The following basic routines are required for each preconditioner.
104b9ad928SBarry Smith           PCCreate_XXX()          - Creates a preconditioner context
114b9ad928SBarry Smith           PCSetFromOptions_XXX()  - Sets runtime options
124b9ad928SBarry Smith           PCApply_XXX()           - Applies the preconditioner
134b9ad928SBarry Smith           PCDestroy_XXX()         - Destroys the preconditioner context
144b9ad928SBarry Smith      where the suffix "_XXX" denotes a particular implementation, in
154b9ad928SBarry Smith      this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi).
164b9ad928SBarry Smith      These routines are actually called via the common user interface
174b9ad928SBarry Smith      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
184b9ad928SBarry Smith      so the application code interface remains identical for all
194b9ad928SBarry Smith      preconditioners.
204b9ad928SBarry Smith 
214b9ad928SBarry Smith      Another key routine is:
224b9ad928SBarry Smith           PCSetUp_XXX()           - Prepares for the use of a preconditioner
234b9ad928SBarry Smith      by setting data structures and options.   The interface routine PCSetUp()
244b9ad928SBarry Smith      is not usually called directly by the user, but instead is called by
254b9ad928SBarry Smith      PCApply() if necessary.
264b9ad928SBarry Smith 
274b9ad928SBarry Smith      Additional basic routines are:
284b9ad928SBarry Smith           PCView_XXX()            - Prints details of runtime options that
294b9ad928SBarry Smith                                     have actually been used.
304b9ad928SBarry Smith      These are called by application codes via the interface routines
314b9ad928SBarry Smith      PCView().
324b9ad928SBarry Smith 
334b9ad928SBarry Smith      The various types of solvers (preconditioners, Krylov subspace methods,
344b9ad928SBarry Smith      nonlinear solvers, timesteppers) are all organized similarly, so the
354b9ad928SBarry Smith      above description applies to these categories also.  One exception is
364b9ad928SBarry Smith      that the analogues of PCApply() for these components are KSPSolve(),
374b9ad928SBarry Smith      SNESSolve(), and TSSolve().
384b9ad928SBarry Smith 
394b9ad928SBarry Smith      Additional optional functionality unique to preconditioners is left and
404b9ad928SBarry Smith      right symmetric preconditioner application via PCApplySymmetricLeft()
414b9ad928SBarry Smith      and PCApplySymmetricRight().  The Jacobi implementation is
424b9ad928SBarry Smith      PCApplySymmetricLeftOrRight_Jacobi().
434b9ad928SBarry Smith 
444b9ad928SBarry Smith     -------------------------------------------------------------------- */
454b9ad928SBarry Smith 
464b9ad928SBarry Smith /*
474b9ad928SBarry Smith    Include files needed for the Jacobi preconditioner:
484b9ad928SBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
494b9ad928SBarry Smith */
504b9ad928SBarry Smith 
51af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
524b9ad928SBarry Smith 
530a545947SLisandro Dalcin const char *const PCJacobiTypes[] = {"DIAGONAL", "ROWMAX", "ROWSUM", "PCJacobiType", "PC_JACOBI_", NULL};
54baa89ecbSBarry Smith 
554b9ad928SBarry Smith /*
564b9ad928SBarry Smith    Private context (data structure) for the Jacobi preconditioner.
574b9ad928SBarry Smith */
584b9ad928SBarry Smith typedef struct {
592fa5cd67SKarl Rupp   Vec       diag;      /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */
604b9ad928SBarry Smith   Vec       diagsqrt;  /* vector containing the reciprocals of the square roots of
614b9ad928SBarry Smith                                     the diagonal elements of the preconditioner matrix (used
624b9ad928SBarry Smith                                     only for symmetric preconditioner application) */
63baa89ecbSBarry Smith   PetscBool userowmax; /* set with PCJacobiSetType() */
64ace3abfcSBarry Smith   PetscBool userowsum;
65ace3abfcSBarry Smith   PetscBool useabs;  /* use the absolute values of the diagonal entries */
6667ed0f3bSStefano Zampini   PetscBool fixdiag; /* fix zero diagonal terms */
674b9ad928SBarry Smith } PC_Jacobi;
684b9ad928SBarry Smith 
699371c9d4SSatish Balay static PetscErrorCode PCJacobiSetType_Jacobi(PC pc, PCJacobiType type) {
70baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
714b9ad928SBarry Smith 
724b9ad928SBarry Smith   PetscFunctionBegin;
73baa89ecbSBarry Smith   j->userowmax = PETSC_FALSE;
74baa89ecbSBarry Smith   j->userowsum = PETSC_FALSE;
75baa89ecbSBarry Smith   if (type == PC_JACOBI_ROWMAX) {
764b9ad928SBarry Smith     j->userowmax = PETSC_TRUE;
77baa89ecbSBarry Smith   } else if (type == PC_JACOBI_ROWSUM) {
78baa89ecbSBarry Smith     j->userowsum = PETSC_TRUE;
79baa89ecbSBarry Smith   }
804b9ad928SBarry Smith   PetscFunctionReturn(0);
814b9ad928SBarry Smith }
824b9ad928SBarry Smith 
839371c9d4SSatish Balay static PetscErrorCode PCJacobiGetType_Jacobi(PC pc, PCJacobiType *type) {
84baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
8586697f06SMatthew Knepley 
8686697f06SMatthew Knepley   PetscFunctionBegin;
87baa89ecbSBarry Smith   if (j->userowmax) {
88baa89ecbSBarry Smith     *type = PC_JACOBI_ROWMAX;
89baa89ecbSBarry Smith   } else if (j->userowsum) {
90baa89ecbSBarry Smith     *type = PC_JACOBI_ROWSUM;
91baa89ecbSBarry Smith   } else {
92baa89ecbSBarry Smith     *type = PC_JACOBI_DIAGONAL;
93baa89ecbSBarry Smith   }
9486697f06SMatthew Knepley   PetscFunctionReturn(0);
9586697f06SMatthew Knepley }
9686697f06SMatthew Knepley 
979371c9d4SSatish Balay static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc, PetscBool flg) {
98baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
99cd47f5d9SBarry Smith 
100cd47f5d9SBarry Smith   PetscFunctionBegin;
101baa89ecbSBarry Smith   j->useabs = flg;
102baa89ecbSBarry Smith   PetscFunctionReturn(0);
103baa89ecbSBarry Smith }
104baa89ecbSBarry Smith 
1059371c9d4SSatish Balay static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc, PetscBool *flg) {
106baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
107baa89ecbSBarry Smith 
108baa89ecbSBarry Smith   PetscFunctionBegin;
109baa89ecbSBarry Smith   *flg = j->useabs;
110cd47f5d9SBarry Smith   PetscFunctionReturn(0);
111cd47f5d9SBarry Smith }
112cd47f5d9SBarry Smith 
1139371c9d4SSatish Balay static PetscErrorCode PCJacobiSetFixDiagonal_Jacobi(PC pc, PetscBool flg) {
11467ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
11567ed0f3bSStefano Zampini 
11667ed0f3bSStefano Zampini   PetscFunctionBegin;
11767ed0f3bSStefano Zampini   j->fixdiag = flg;
11867ed0f3bSStefano Zampini   PetscFunctionReturn(0);
11967ed0f3bSStefano Zampini }
12067ed0f3bSStefano Zampini 
1219371c9d4SSatish Balay static PetscErrorCode PCJacobiGetFixDiagonal_Jacobi(PC pc, PetscBool *flg) {
12267ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
12367ed0f3bSStefano Zampini 
12467ed0f3bSStefano Zampini   PetscFunctionBegin;
12567ed0f3bSStefano Zampini   *flg = j->fixdiag;
12667ed0f3bSStefano Zampini   PetscFunctionReturn(0);
12767ed0f3bSStefano Zampini }
12867ed0f3bSStefano Zampini 
1294b9ad928SBarry Smith /*
1304b9ad928SBarry Smith    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1314b9ad928SBarry Smith                     by setting data structures and options.
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith    Input Parameter:
1344b9ad928SBarry Smith .  pc - the preconditioner context
1354b9ad928SBarry Smith 
1364b9ad928SBarry Smith    Application Interface Routine: PCSetUp()
1374b9ad928SBarry Smith 
138f1580f4eSBarry Smith    Note:
1394b9ad928SBarry Smith    The interface routine PCSetUp() is not usually called directly by
1404b9ad928SBarry Smith    the user, but instead is called by PCApply() if necessary.
1414b9ad928SBarry Smith */
1429371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi(PC pc) {
1434b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
1444b9ad928SBarry Smith   Vec          diag, diagsqrt;
145cd47f5d9SBarry Smith   PetscInt     n, i;
1464b9ad928SBarry Smith   PetscScalar *x;
147ace3abfcSBarry Smith   PetscBool    zeroflag = PETSC_FALSE;
1484b9ad928SBarry Smith 
1494b9ad928SBarry Smith   PetscFunctionBegin;
1504b9ad928SBarry Smith   /*
1514b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1524b9ad928SBarry Smith 
1534b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1549566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->mat,&jac->diag));
1554b9ad928SBarry Smith   }
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1584b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1594b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1604b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1614b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1624b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1634b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1644b9ad928SBarry Smith   */
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith   /*
1674b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1684b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1694b9ad928SBarry Smith   */
1704b9ad928SBarry Smith   diag     = jac->diag;
1714b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1724b9ad928SBarry Smith 
1734b9ad928SBarry Smith   if (diag) {
174b94d7dedSBarry Smith     PetscBool isset, isspd;
17567ed0f3bSStefano Zampini 
1764b9ad928SBarry Smith     if (jac->userowmax) {
1779566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL));
17886697f06SMatthew Knepley     } else if (jac->userowsum) {
1799566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diag));
1804b9ad928SBarry Smith     } else {
1819566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diag));
1824b9ad928SBarry Smith     }
1839566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
1841baa6e33SBarry Smith     if (jac->useabs) PetscCall(VecAbs(diag));
185b94d7dedSBarry Smith     PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
186b94d7dedSBarry Smith     if (jac->fixdiag && (!isset || !isspd)) {
1879566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(diag, &n));
1889566063dSJacob Faibussowitsch       PetscCall(VecGetArray(diag, &x));
1894b9ad928SBarry Smith       for (i = 0; i < n; i++) {
1904b9ad928SBarry Smith         if (x[i] == 0.0) {
1914b9ad928SBarry Smith           x[i]     = 1.0;
192cd47f5d9SBarry Smith           zeroflag = PETSC_TRUE;
1934b9ad928SBarry Smith         }
1944b9ad928SBarry Smith       }
1959566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(diag, &x));
1964b9ad928SBarry Smith     }
19767ed0f3bSStefano Zampini   }
1984b9ad928SBarry Smith   if (diagsqrt) {
1994b9ad928SBarry Smith     if (jac->userowmax) {
2009566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL));
20186697f06SMatthew Knepley     } else if (jac->userowsum) {
2029566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diagsqrt));
2034b9ad928SBarry Smith     } else {
2049566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diagsqrt));
2054b9ad928SBarry Smith     }
2069566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(diagsqrt, &n));
2079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(diagsqrt, &x));
2084b9ad928SBarry Smith     for (i = 0; i < n; i++) {
2098f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i]));
2104b9ad928SBarry Smith       else {
2114b9ad928SBarry Smith         x[i]     = 1.0;
212cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2134b9ad928SBarry Smith       }
2144b9ad928SBarry Smith     }
2159566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(diagsqrt, &x));
2164b9ad928SBarry Smith   }
21748a46eb9SPierre Jolivet   if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
2184b9ad928SBarry Smith   PetscFunctionReturn(0);
2194b9ad928SBarry Smith }
220f1580f4eSBarry Smith 
2214b9ad928SBarry Smith /*
2224b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2234b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2244b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2254b9ad928SBarry Smith 
2264b9ad928SBarry Smith    Input Parameter:
2274b9ad928SBarry Smith .  pc - the preconditioner context
2284b9ad928SBarry Smith */
2299371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) {
2304b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2354b9ad928SBarry Smith   PetscFunctionReturn(0);
2364b9ad928SBarry Smith }
237f1580f4eSBarry Smith 
2384b9ad928SBarry Smith /*
2394b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2404b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2414b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2424b9ad928SBarry Smith 
2434b9ad928SBarry Smith    Input Parameter:
2444b9ad928SBarry Smith .  pc - the preconditioner context
2454b9ad928SBarry Smith */
2469371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) {
2474b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL));
2519566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2524b9ad928SBarry Smith   PetscFunctionReturn(0);
2534b9ad928SBarry Smith }
254f1580f4eSBarry Smith 
2554b9ad928SBarry Smith /*
2564b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith    Input Parameters:
2594b9ad928SBarry Smith .  pc - the preconditioner context
2604b9ad928SBarry Smith .  x - input vector
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith    Output Parameter:
2634b9ad928SBarry Smith .  y - output vector
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Application Interface Routine: PCApply()
2664b9ad928SBarry Smith  */
2679371c9d4SSatish Balay static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y) {
2684b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith   PetscFunctionBegin;
27148a46eb9SPierre Jolivet   if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
2729566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diag));
2734b9ad928SBarry Smith   PetscFunctionReturn(0);
2744b9ad928SBarry Smith }
275f1580f4eSBarry Smith 
2764b9ad928SBarry Smith /*
2774b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2784b9ad928SBarry Smith    symmetric preconditioner to a vector.
2794b9ad928SBarry Smith 
2804b9ad928SBarry Smith    Input Parameters:
2814b9ad928SBarry Smith .  pc - the preconditioner context
2824b9ad928SBarry Smith .  x - input vector
2834b9ad928SBarry Smith 
2844b9ad928SBarry Smith    Output Parameter:
2854b9ad928SBarry Smith .  y - output vector
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
2884b9ad928SBarry Smith */
2899371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y) {
2904b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith   PetscFunctionBegin;
29348a46eb9SPierre Jolivet   if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc));
2949566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diagsqrt));
2954b9ad928SBarry Smith   PetscFunctionReturn(0);
2964b9ad928SBarry Smith }
29767ed0f3bSStefano Zampini 
2989371c9d4SSatish Balay static PetscErrorCode PCReset_Jacobi(PC pc) {
299a06653b4SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
300a06653b4SBarry Smith 
301a06653b4SBarry Smith   PetscFunctionBegin;
3029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
3039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diagsqrt));
304a06653b4SBarry Smith   PetscFunctionReturn(0);
305a06653b4SBarry Smith }
306a06653b4SBarry Smith 
3074b9ad928SBarry Smith /*
3084b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3094b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3104b9ad928SBarry Smith 
3114b9ad928SBarry Smith    Input Parameter:
3124b9ad928SBarry Smith .  pc - the preconditioner context
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3154b9ad928SBarry Smith */
3169371c9d4SSatish Balay static PetscErrorCode PCDestroy_Jacobi(PC pc) {
3174b9ad928SBarry Smith   PetscFunctionBegin;
3189566063dSJacob Faibussowitsch   PetscCall(PCReset_Jacobi(pc));
3199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL));
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL));
3219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL));
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL));
3239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL));
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL));
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith   /*
3274b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3284b9ad928SBarry Smith   */
3299566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3304b9ad928SBarry Smith   PetscFunctionReturn(0);
3314b9ad928SBarry Smith }
3324b9ad928SBarry Smith 
3339371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject) {
3344b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
335baa89ecbSBarry Smith   PetscBool    flg;
336baa89ecbSBarry Smith   PCJacobiType deflt, type;
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(PCJacobiGetType(pc, &deflt));
340d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options");
3419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg));
3421baa6e33SBarry Smith   if (flg) PetscCall(PCJacobiSetType(pc, type));
3439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL));
3449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL));
345d0609cedSBarry Smith   PetscOptionsHeadEnd();
3464b9ad928SBarry Smith   PetscFunctionReturn(0);
3474b9ad928SBarry Smith }
3484b9ad928SBarry Smith 
3499371c9d4SSatish Balay static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer) {
350569e28a7SMatthew G. Knepley   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
351569e28a7SMatthew G. Knepley   PetscBool  iascii;
352569e28a7SMatthew G. Knepley 
353569e28a7SMatthew G. Knepley   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
355569e28a7SMatthew G. Knepley   if (iascii) {
356569e28a7SMatthew G. Knepley     PCJacobiType      type;
35767ed0f3bSStefano Zampini     PetscBool         useAbs, fixdiag;
358569e28a7SMatthew G. Knepley     PetscViewerFormat format;
359569e28a7SMatthew G. Knepley 
3609566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetType(pc, &type));
3619566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
3629566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
3639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
3649566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
36548a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(VecView(jac->diag, viewer));
366569e28a7SMatthew G. Knepley   }
367569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
368569e28a7SMatthew G. Knepley }
369569e28a7SMatthew G. Knepley 
3704b9ad928SBarry Smith /*
3714b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
3724b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
3734b9ad928SBarry Smith    context, PC, that was created within PCCreate().
3744b9ad928SBarry Smith 
3754b9ad928SBarry Smith    Input Parameter:
3764b9ad928SBarry Smith .  pc - the preconditioner context
3774b9ad928SBarry Smith 
3784b9ad928SBarry Smith    Application Interface Routine: PCCreate()
3794b9ad928SBarry Smith */
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith /*MC
3825a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
3834b9ad928SBarry Smith 
384f1580f4eSBarry Smith    Options Database Keys:
385967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
38667ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
387147403d9SBarry Smith -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
3884b9ad928SBarry Smith 
3894b9ad928SBarry Smith    Level: beginner
3904b9ad928SBarry Smith 
39195452b02SPatrick Sanan   Notes:
392f1580f4eSBarry Smith     By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric
3934b9ad928SBarry Smith     can scale each side of the matrix by the square root of the diagonal entries.
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith     Zero entries along the diagonal are replaced with the value 1.0
3964b9ad928SBarry Smith 
397f1580f4eSBarry Smith     See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks
398422a814eSBarry Smith 
399db781477SPatrick Sanan .seealso:  `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
400f1580f4eSBarry Smith            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`,
401db781477SPatrick Sanan            `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
402db781477SPatrick Sanan            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
4034b9ad928SBarry Smith M*/
4044b9ad928SBarry Smith 
4059371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) {
4064b9ad928SBarry Smith   PC_Jacobi *jac;
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith   PetscFunctionBegin;
4094b9ad928SBarry Smith   /*
4104b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4114b9ad928SBarry Smith      attach it to the PC object.
4124b9ad928SBarry Smith   */
413*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
4144b9ad928SBarry Smith   pc->data = (void *)jac;
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith   /*
4174b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4184b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4194b9ad928SBarry Smith   */
4200a545947SLisandro Dalcin   jac->diag      = NULL;
4210a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4224b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
42386697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
424cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
42567ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith   /*
4284b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4294b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4304b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4314b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4324b9ad928SBarry Smith       not needed.
4334b9ad928SBarry Smith   */
4344b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4354b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4364b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
437a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4384b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4394b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
440569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4410a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4424b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4434b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4442fa5cd67SKarl Rupp 
4459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi));
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi));
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi));
4489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi));
4499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi));
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi));
4514b9ad928SBarry Smith   PetscFunctionReturn(0);
4524b9ad928SBarry Smith }
453cd47f5d9SBarry Smith 
454cd47f5d9SBarry Smith /*@
455f1580f4eSBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the
45667ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
457cd47f5d9SBarry Smith 
458f1580f4eSBarry Smith    Logically Collective on pc
459cd47f5d9SBarry Smith 
460cd47f5d9SBarry Smith    Input Parameters:
461baa89ecbSBarry Smith +  pc - the preconditioner context
462baa89ecbSBarry Smith -  flg - whether to use absolute values or not
463baa89ecbSBarry Smith 
464baa89ecbSBarry Smith    Options Database Key:
465147403d9SBarry Smith .  -pc_jacobi_abs <bool> - use absolute values
466baa89ecbSBarry Smith 
467f1580f4eSBarry Smith    Note:
46895452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
469baa89ecbSBarry Smith 
470baa89ecbSBarry Smith    Level: intermediate
471baa89ecbSBarry Smith 
472f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiGetUseAbs()`
473baa89ecbSBarry Smith @*/
4749371c9d4SSatish Balay PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg) {
475baa89ecbSBarry Smith   PetscFunctionBegin;
476baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
477cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg));
478baa89ecbSBarry Smith   PetscFunctionReturn(0);
479baa89ecbSBarry Smith }
480baa89ecbSBarry Smith 
481baa89ecbSBarry Smith /*@
482f1580f4eSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the
48367ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
484baa89ecbSBarry Smith 
485f1580f4eSBarry Smith    Logically Collective on pc
486baa89ecbSBarry Smith 
487baa89ecbSBarry Smith    Input Parameter:
488cd47f5d9SBarry Smith .  pc - the preconditioner context
489cd47f5d9SBarry Smith 
490baa89ecbSBarry Smith    Output Parameter:
491baa89ecbSBarry Smith .  flg - whether to use absolute values or not
492cd47f5d9SBarry Smith 
493cd47f5d9SBarry Smith    Level: intermediate
494cd47f5d9SBarry Smith 
495f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
496cd47f5d9SBarry Smith @*/
4979371c9d4SSatish Balay PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg) {
498cd47f5d9SBarry Smith   PetscFunctionBegin;
4990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
500cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg));
501cd47f5d9SBarry Smith   PetscFunctionReturn(0);
502cd47f5d9SBarry Smith }
503cd47f5d9SBarry Smith 
5044b9ad928SBarry Smith /*@
505147403d9SBarry Smith    PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
50667ed0f3bSStefano Zampini 
507f1580f4eSBarry Smith    Logically Collective on pc
50867ed0f3bSStefano Zampini 
50967ed0f3bSStefano Zampini    Input Parameters:
51067ed0f3bSStefano Zampini +  pc - the preconditioner context
51167ed0f3bSStefano Zampini -  flg - the boolean flag
51267ed0f3bSStefano Zampini 
51367ed0f3bSStefano Zampini    Options Database Key:
514147403d9SBarry Smith .  -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
51567ed0f3bSStefano Zampini 
516f1580f4eSBarry Smith    Note:
51767ed0f3bSStefano Zampini    This takes affect at the next construction of the preconditioner
51867ed0f3bSStefano Zampini 
51967ed0f3bSStefano Zampini    Level: intermediate
52067ed0f3bSStefano Zampini 
521f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()`
52267ed0f3bSStefano Zampini @*/
5239371c9d4SSatish Balay PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg) {
52467ed0f3bSStefano Zampini   PetscFunctionBegin;
52567ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
526cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg));
52767ed0f3bSStefano Zampini   PetscFunctionReturn(0);
52867ed0f3bSStefano Zampini }
52967ed0f3bSStefano Zampini 
53067ed0f3bSStefano Zampini /*@
531f1580f4eSBarry Smith    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms
53267ed0f3bSStefano Zampini 
533f1580f4eSBarry Smith    Logically Collective on pc
53467ed0f3bSStefano Zampini 
53567ed0f3bSStefano Zampini    Input Parameter:
53667ed0f3bSStefano Zampini .  pc - the preconditioner context
53767ed0f3bSStefano Zampini 
53867ed0f3bSStefano Zampini    Output Parameter:
53967ed0f3bSStefano Zampini .  flg - the boolean flag
54067ed0f3bSStefano Zampini 
54167ed0f3bSStefano Zampini    Options Database Key:
54267b8a455SSatish Balay .  -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
54367ed0f3bSStefano Zampini 
54467ed0f3bSStefano Zampini    Level: intermediate
54567ed0f3bSStefano Zampini 
546f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
54767ed0f3bSStefano Zampini @*/
5489371c9d4SSatish Balay PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg) {
54967ed0f3bSStefano Zampini   PetscFunctionBegin;
55067ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
551cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg));
55267ed0f3bSStefano Zampini   PetscFunctionReturn(0);
55367ed0f3bSStefano Zampini }
55467ed0f3bSStefano Zampini 
55567ed0f3bSStefano Zampini /*@
556baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
557baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
5584b9ad928SBarry Smith 
559f1580f4eSBarry Smith    Logically Collective on pc
5604b9ad928SBarry Smith 
5614b9ad928SBarry Smith    Input Parameters:
562baa89ecbSBarry Smith +  pc - the preconditioner context
563f1580f4eSBarry Smith -  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith    Options Database Key:
566147403d9SBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
5674b9ad928SBarry Smith 
5684b9ad928SBarry Smith    Level: intermediate
5694b9ad928SBarry Smith 
570f1580f4eSBarry Smith    Developer Note:
571f1580f4eSBarry Smith    Why is there a separate function for using the absolute value?
572f1580f4eSBarry Smith 
573f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
5744b9ad928SBarry Smith @*/
5759371c9d4SSatish Balay PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type) {
5764b9ad928SBarry Smith   PetscFunctionBegin;
5770700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
578cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type));
5794b9ad928SBarry Smith   PetscFunctionReturn(0);
5804b9ad928SBarry Smith }
5814b9ad928SBarry Smith 
58286697f06SMatthew Knepley /*@
583baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
58486697f06SMatthew Knepley 
585f1580f4eSBarry Smith    Not Collective on pc
58686697f06SMatthew Knepley 
587baa89ecbSBarry Smith    Input Parameter:
58886697f06SMatthew Knepley .  pc - the preconditioner context
58986697f06SMatthew Knepley 
590baa89ecbSBarry Smith    Output Parameter:
591f1580f4eSBarry Smith .  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
59286697f06SMatthew Knepley 
59386697f06SMatthew Knepley    Level: intermediate
59486697f06SMatthew Knepley 
595f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaUseAbs()`, `PCJacobiSetType()`
59686697f06SMatthew Knepley @*/
5979371c9d4SSatish Balay PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type) {
59886697f06SMatthew Knepley   PetscFunctionBegin;
5990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
600cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type));
60186697f06SMatthew Knepley   PetscFunctionReturn(0);
60286697f06SMatthew Knepley }
603