xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision 04c3f3b8f563eff258c1de90d4e02d41e6027585)
14b9ad928SBarry Smith 
2*04c3f3b8SBarry Smith /*
35a46d3faSBarry Smith      This file implements a Jacobi preconditioner in PETSc as part of PC.
45a46d3faSBarry Smith      You can use this as a starting point for implementing your own
55a46d3faSBarry Smith      preconditioner that is not provided with PETSc. (You might also consider
65a46d3faSBarry Smith      just using PCSHELL)
74b9ad928SBarry Smith 
84b9ad928SBarry Smith      The following basic routines are required for each preconditioner.
94b9ad928SBarry Smith           PCCreate_XXX()          - Creates a preconditioner context
104b9ad928SBarry Smith           PCSetFromOptions_XXX()  - Sets runtime options
114b9ad928SBarry Smith           PCApply_XXX()           - Applies the preconditioner
124b9ad928SBarry Smith           PCDestroy_XXX()         - Destroys the preconditioner context
134b9ad928SBarry Smith      where the suffix "_XXX" denotes a particular implementation, in
144b9ad928SBarry Smith      this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi).
154b9ad928SBarry Smith      These routines are actually called via the common user interface
164b9ad928SBarry Smith      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
174b9ad928SBarry Smith      so the application code interface remains identical for all
184b9ad928SBarry Smith      preconditioners.
194b9ad928SBarry Smith 
204b9ad928SBarry Smith      Another key routine is:
214b9ad928SBarry Smith           PCSetUp_XXX()           - Prepares for the use of a preconditioner
224b9ad928SBarry Smith      by setting data structures and options.   The interface routine PCSetUp()
234b9ad928SBarry Smith      is not usually called directly by the user, but instead is called by
244b9ad928SBarry Smith      PCApply() if necessary.
254b9ad928SBarry Smith 
264b9ad928SBarry Smith      Additional basic routines are:
274b9ad928SBarry Smith           PCView_XXX()            - Prints details of runtime options that
284b9ad928SBarry Smith                                     have actually been used.
294b9ad928SBarry Smith      These are called by application codes via the interface routines
304b9ad928SBarry Smith      PCView().
314b9ad928SBarry Smith 
324b9ad928SBarry Smith      The various types of solvers (preconditioners, Krylov subspace methods,
334b9ad928SBarry Smith      nonlinear solvers, timesteppers) are all organized similarly, so the
344b9ad928SBarry Smith      above description applies to these categories also.  One exception is
354b9ad928SBarry Smith      that the analogues of PCApply() for these components are KSPSolve(),
364b9ad928SBarry Smith      SNESSolve(), and TSSolve().
374b9ad928SBarry Smith 
384b9ad928SBarry Smith      Additional optional functionality unique to preconditioners is left and
394b9ad928SBarry Smith      right symmetric preconditioner application via PCApplySymmetricLeft()
404b9ad928SBarry Smith      and PCApplySymmetricRight().  The Jacobi implementation is
414b9ad928SBarry Smith      PCApplySymmetricLeftOrRight_Jacobi().
42*04c3f3b8SBarry Smith */
434b9ad928SBarry Smith 
444b9ad928SBarry Smith /*
454b9ad928SBarry Smith    Include files needed for the Jacobi preconditioner:
464b9ad928SBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
474b9ad928SBarry Smith */
484b9ad928SBarry Smith 
49af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
504b9ad928SBarry Smith 
510a545947SLisandro Dalcin const char *const PCJacobiTypes[] = {"DIAGONAL", "ROWMAX", "ROWSUM", "PCJacobiType", "PC_JACOBI_", NULL};
52baa89ecbSBarry Smith 
534b9ad928SBarry Smith /*
544b9ad928SBarry Smith    Private context (data structure) for the Jacobi preconditioner.
554b9ad928SBarry Smith */
564b9ad928SBarry Smith typedef struct {
572fa5cd67SKarl Rupp   Vec       diag;      /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */
584b9ad928SBarry Smith   Vec       diagsqrt;  /* vector containing the reciprocals of the square roots of
594b9ad928SBarry Smith                                     the diagonal elements of the preconditioner matrix (used
604b9ad928SBarry Smith                                     only for symmetric preconditioner application) */
61baa89ecbSBarry Smith   PetscBool userowmax; /* set with PCJacobiSetType() */
62ace3abfcSBarry Smith   PetscBool userowsum;
63ace3abfcSBarry Smith   PetscBool useabs;  /* use the absolute values of the diagonal entries */
6467ed0f3bSStefano Zampini   PetscBool fixdiag; /* fix zero diagonal terms */
654b9ad928SBarry Smith } PC_Jacobi;
664b9ad928SBarry Smith 
677a660b1eSJed Brown static PetscErrorCode PCReset_Jacobi(PC);
687a660b1eSJed Brown 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetType_Jacobi(PC pc, PCJacobiType type)
70d71ae5a4SJacob Faibussowitsch {
71baa89ecbSBarry Smith   PC_Jacobi   *j = (PC_Jacobi *)pc->data;
727a660b1eSJed Brown   PCJacobiType old_type;
734b9ad928SBarry Smith 
744b9ad928SBarry Smith   PetscFunctionBegin;
757a660b1eSJed Brown   PetscCall(PCJacobiGetType(pc, &old_type));
763ba16761SJacob Faibussowitsch   if (old_type == type) PetscFunctionReturn(PETSC_SUCCESS);
777a660b1eSJed Brown   PetscCall(PCReset_Jacobi(pc));
78baa89ecbSBarry Smith   j->userowmax = PETSC_FALSE;
79baa89ecbSBarry Smith   j->userowsum = PETSC_FALSE;
80baa89ecbSBarry Smith   if (type == PC_JACOBI_ROWMAX) {
814b9ad928SBarry Smith     j->userowmax = PETSC_TRUE;
82baa89ecbSBarry Smith   } else if (type == PC_JACOBI_ROWSUM) {
83baa89ecbSBarry Smith     j->userowsum = PETSC_TRUE;
84baa89ecbSBarry Smith   }
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
864b9ad928SBarry Smith }
874b9ad928SBarry Smith 
88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetType_Jacobi(PC pc, PCJacobiType *type)
89d71ae5a4SJacob Faibussowitsch {
90baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
9186697f06SMatthew Knepley 
9286697f06SMatthew Knepley   PetscFunctionBegin;
93baa89ecbSBarry Smith   if (j->userowmax) {
94baa89ecbSBarry Smith     *type = PC_JACOBI_ROWMAX;
95baa89ecbSBarry Smith   } else if (j->userowsum) {
96baa89ecbSBarry Smith     *type = PC_JACOBI_ROWSUM;
97baa89ecbSBarry Smith   } else {
98baa89ecbSBarry Smith     *type = PC_JACOBI_DIAGONAL;
99baa89ecbSBarry Smith   }
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10186697f06SMatthew Knepley }
10286697f06SMatthew Knepley 
103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc, PetscBool flg)
104d71ae5a4SJacob Faibussowitsch {
105baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
106cd47f5d9SBarry Smith 
107cd47f5d9SBarry Smith   PetscFunctionBegin;
108baa89ecbSBarry Smith   j->useabs = flg;
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110baa89ecbSBarry Smith }
111baa89ecbSBarry Smith 
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc, PetscBool *flg)
113d71ae5a4SJacob Faibussowitsch {
114baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
115baa89ecbSBarry Smith 
116baa89ecbSBarry Smith   PetscFunctionBegin;
117baa89ecbSBarry Smith   *flg = j->useabs;
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119cd47f5d9SBarry Smith }
120cd47f5d9SBarry Smith 
121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetFixDiagonal_Jacobi(PC pc, PetscBool flg)
122d71ae5a4SJacob Faibussowitsch {
12367ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
12467ed0f3bSStefano Zampini 
12567ed0f3bSStefano Zampini   PetscFunctionBegin;
12667ed0f3bSStefano Zampini   j->fixdiag = flg;
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12867ed0f3bSStefano Zampini }
12967ed0f3bSStefano Zampini 
130d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetFixDiagonal_Jacobi(PC pc, PetscBool *flg)
131d71ae5a4SJacob Faibussowitsch {
13267ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
13367ed0f3bSStefano Zampini 
13467ed0f3bSStefano Zampini   PetscFunctionBegin;
13567ed0f3bSStefano Zampini   *flg = j->fixdiag;
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13767ed0f3bSStefano Zampini }
13867ed0f3bSStefano Zampini 
1394b9ad928SBarry Smith /*
1404b9ad928SBarry Smith    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1414b9ad928SBarry Smith                     by setting data structures and options.
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith    Input Parameter:
1444b9ad928SBarry Smith .  pc - the preconditioner context
1454b9ad928SBarry Smith 
1464b9ad928SBarry Smith    Application Interface Routine: PCSetUp()
1474b9ad928SBarry Smith 
148f1580f4eSBarry Smith    Note:
1494b9ad928SBarry Smith    The interface routine PCSetUp() is not usually called directly by
1504b9ad928SBarry Smith    the user, but instead is called by PCApply() if necessary.
1514b9ad928SBarry Smith */
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi(PC pc)
153d71ae5a4SJacob Faibussowitsch {
1544b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
1554b9ad928SBarry Smith   Vec          diag, diagsqrt;
156cd47f5d9SBarry Smith   PetscInt     n, i;
1574b9ad928SBarry Smith   PetscScalar *x;
158ace3abfcSBarry Smith   PetscBool    zeroflag = PETSC_FALSE;
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith   PetscFunctionBegin;
1614b9ad928SBarry Smith   /*
1624b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1659566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->mat,&jac->diag));
1664b9ad928SBarry Smith   }
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1694b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1704b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1714b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1724b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1734b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1744b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1754b9ad928SBarry Smith   */
1764b9ad928SBarry Smith 
1774b9ad928SBarry Smith   /*
1784b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1794b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1804b9ad928SBarry Smith   */
1814b9ad928SBarry Smith   diag     = jac->diag;
1824b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   if (diag) {
185b94d7dedSBarry Smith     PetscBool isset, isspd;
18667ed0f3bSStefano Zampini 
1874b9ad928SBarry Smith     if (jac->userowmax) {
1889566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL));
18986697f06SMatthew Knepley     } else if (jac->userowsum) {
1909566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diag));
1914b9ad928SBarry Smith     } else {
1929566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diag));
1934b9ad928SBarry Smith     }
1949566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
1951baa6e33SBarry Smith     if (jac->useabs) PetscCall(VecAbs(diag));
196b94d7dedSBarry Smith     PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
197b94d7dedSBarry Smith     if (jac->fixdiag && (!isset || !isspd)) {
1989566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(diag, &n));
1999566063dSJacob Faibussowitsch       PetscCall(VecGetArray(diag, &x));
2004b9ad928SBarry Smith       for (i = 0; i < n; i++) {
2014b9ad928SBarry Smith         if (x[i] == 0.0) {
2024b9ad928SBarry Smith           x[i]     = 1.0;
203cd47f5d9SBarry Smith           zeroflag = PETSC_TRUE;
2044b9ad928SBarry Smith         }
2054b9ad928SBarry Smith       }
2069566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(diag, &x));
2074b9ad928SBarry Smith     }
20867ed0f3bSStefano Zampini   }
2094b9ad928SBarry Smith   if (diagsqrt) {
2104b9ad928SBarry Smith     if (jac->userowmax) {
2119566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL));
21286697f06SMatthew Knepley     } else if (jac->userowsum) {
2139566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diagsqrt));
2144b9ad928SBarry Smith     } else {
2159566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diagsqrt));
2164b9ad928SBarry Smith     }
2179566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(diagsqrt, &n));
2189566063dSJacob Faibussowitsch     PetscCall(VecGetArray(diagsqrt, &x));
2194b9ad928SBarry Smith     for (i = 0; i < n; i++) {
2208f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i]));
2214b9ad928SBarry Smith       else {
2224b9ad928SBarry Smith         x[i]     = 1.0;
223cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2244b9ad928SBarry Smith       }
2254b9ad928SBarry Smith     }
2269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(diagsqrt, &x));
2274b9ad928SBarry Smith   }
22848a46eb9SPierre Jolivet   if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2304b9ad928SBarry Smith }
231f1580f4eSBarry Smith 
2324b9ad928SBarry Smith /*
2334b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2344b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2354b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith    Input Parameter:
2384b9ad928SBarry Smith .  pc - the preconditioner context
2394b9ad928SBarry Smith */
240d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
241d71ae5a4SJacob Faibussowitsch {
2424b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2434b9ad928SBarry Smith 
2444b9ad928SBarry Smith   PetscFunctionBegin;
2459566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL));
2469566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2484b9ad928SBarry Smith }
249f1580f4eSBarry Smith 
2504b9ad928SBarry Smith /*
2514b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2524b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2534b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith    Input Parameter:
2564b9ad928SBarry Smith .  pc - the preconditioner context
2574b9ad928SBarry Smith */
258d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
259d71ae5a4SJacob Faibussowitsch {
2604b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL));
2649566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2664b9ad928SBarry Smith }
267f1580f4eSBarry Smith 
2684b9ad928SBarry Smith /*
2694b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Input Parameters:
2724b9ad928SBarry Smith .  pc - the preconditioner context
2734b9ad928SBarry Smith .  x - input vector
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith    Output Parameter:
2764b9ad928SBarry Smith .  y - output vector
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith    Application Interface Routine: PCApply()
2794b9ad928SBarry Smith  */
280d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y)
281d71ae5a4SJacob Faibussowitsch {
2824b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2834b9ad928SBarry Smith 
2844b9ad928SBarry Smith   PetscFunctionBegin;
28548a46eb9SPierre Jolivet   if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
2869566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diag));
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2884b9ad928SBarry Smith }
289f1580f4eSBarry Smith 
2904b9ad928SBarry Smith /*
2914b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2924b9ad928SBarry Smith    symmetric preconditioner to a vector.
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith    Input Parameters:
2954b9ad928SBarry Smith .  pc - the preconditioner context
2964b9ad928SBarry Smith .  x - input vector
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith    Output Parameter:
2994b9ad928SBarry Smith .  y - output vector
3004b9ad928SBarry Smith 
3014b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
3024b9ad928SBarry Smith */
303d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y)
304d71ae5a4SJacob Faibussowitsch {
3054b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith   PetscFunctionBegin;
30848a46eb9SPierre Jolivet   if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc));
3099566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diagsqrt));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3114b9ad928SBarry Smith }
31267ed0f3bSStefano Zampini 
313d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Jacobi(PC pc)
314d71ae5a4SJacob Faibussowitsch {
315a06653b4SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
316a06653b4SBarry Smith 
317a06653b4SBarry Smith   PetscFunctionBegin;
3189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diagsqrt));
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
321a06653b4SBarry Smith }
322a06653b4SBarry Smith 
3234b9ad928SBarry Smith /*
3244b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3254b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3264b9ad928SBarry Smith 
3274b9ad928SBarry Smith    Input Parameter:
3284b9ad928SBarry Smith .  pc - the preconditioner context
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3314b9ad928SBarry Smith */
332d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Jacobi(PC pc)
333d71ae5a4SJacob Faibussowitsch {
3344b9ad928SBarry Smith   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(PCReset_Jacobi(pc));
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL));
3399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL));
3409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL));
3419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL));
3424b9ad928SBarry Smith 
3434b9ad928SBarry Smith   /*
3444b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3454b9ad928SBarry Smith   */
3469566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3484b9ad928SBarry Smith }
3494b9ad928SBarry Smith 
350d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject)
351d71ae5a4SJacob Faibussowitsch {
3524b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
353baa89ecbSBarry Smith   PetscBool    flg;
354baa89ecbSBarry Smith   PCJacobiType deflt, type;
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(PCJacobiGetType(pc, &deflt));
358d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options");
3599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg));
3601baa6e33SBarry Smith   if (flg) PetscCall(PCJacobiSetType(pc, type));
3619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL));
3629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL));
363d0609cedSBarry Smith   PetscOptionsHeadEnd();
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3654b9ad928SBarry Smith }
3664b9ad928SBarry Smith 
367d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
368d71ae5a4SJacob Faibussowitsch {
369569e28a7SMatthew G. Knepley   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
370569e28a7SMatthew G. Knepley   PetscBool  iascii;
371569e28a7SMatthew G. Knepley 
372569e28a7SMatthew G. Knepley   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
374569e28a7SMatthew G. Knepley   if (iascii) {
375569e28a7SMatthew G. Knepley     PCJacobiType      type;
37667ed0f3bSStefano Zampini     PetscBool         useAbs, fixdiag;
377569e28a7SMatthew G. Knepley     PetscViewerFormat format;
378569e28a7SMatthew G. Knepley 
3799566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetType(pc, &type));
3809566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
3819566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
3829566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
3839566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
384b24842d1SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL && jac->diag) PetscCall(VecView(jac->diag, viewer));
385569e28a7SMatthew G. Knepley   }
3863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
387569e28a7SMatthew G. Knepley }
388569e28a7SMatthew G. Knepley 
3894b9ad928SBarry Smith /*
3904b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
3914b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
3924b9ad928SBarry Smith    context, PC, that was created within PCCreate().
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith    Input Parameter:
3954b9ad928SBarry Smith .  pc - the preconditioner context
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Application Interface Routine: PCCreate()
3984b9ad928SBarry Smith */
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith /*MC
4015a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
4024b9ad928SBarry Smith 
403f1580f4eSBarry Smith    Options Database Keys:
404967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
40567ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
406147403d9SBarry Smith -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith    Level: beginner
4094b9ad928SBarry Smith 
41095452b02SPatrick Sanan   Notes:
411f1580f4eSBarry Smith     By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric
4124b9ad928SBarry Smith     can scale each side of the matrix by the square root of the diagonal entries.
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith     Zero entries along the diagonal are replaced with the value 1.0
4154b9ad928SBarry Smith 
416f1580f4eSBarry Smith     See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks
417422a814eSBarry Smith 
418db781477SPatrick Sanan .seealso:  `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
419f1580f4eSBarry Smith            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`,
420db781477SPatrick Sanan            `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
421db781477SPatrick Sanan            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
4224b9ad928SBarry Smith M*/
4234b9ad928SBarry Smith 
424d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
425d71ae5a4SJacob Faibussowitsch {
4264b9ad928SBarry Smith   PC_Jacobi *jac;
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith   PetscFunctionBegin;
4294b9ad928SBarry Smith   /*
4304b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4314b9ad928SBarry Smith      attach it to the PC object.
4324b9ad928SBarry Smith   */
4334dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
4344b9ad928SBarry Smith   pc->data = (void *)jac;
4354b9ad928SBarry Smith 
4364b9ad928SBarry Smith   /*
4374b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4384b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4394b9ad928SBarry Smith   */
4400a545947SLisandro Dalcin   jac->diag      = NULL;
4410a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4424b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
44386697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
444cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
44567ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith   /*
4484b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4494b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4504b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4514b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4524b9ad928SBarry Smith       not needed.
4534b9ad928SBarry Smith   */
4544b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4554b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4564b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
457a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4584b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4594b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
460569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4610a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4624b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4634b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4642fa5cd67SKarl Rupp 
4659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi));
4669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi));
4679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi));
4689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi));
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi));
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi));
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4724b9ad928SBarry Smith }
473cd47f5d9SBarry Smith 
474cd47f5d9SBarry Smith /*@
475f1580f4eSBarry Smith   PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the
47667ed0f3bSStefano Zampini   absolute values of the diagonal divisors in the preconditioner
477cd47f5d9SBarry Smith 
478c3339decSBarry Smith   Logically Collective
479cd47f5d9SBarry Smith 
480cd47f5d9SBarry Smith   Input Parameters:
481baa89ecbSBarry Smith + pc  - the preconditioner context
482baa89ecbSBarry Smith - flg - whether to use absolute values or not
483baa89ecbSBarry Smith 
484baa89ecbSBarry Smith   Options Database Key:
485147403d9SBarry Smith . -pc_jacobi_abs <bool> - use absolute values
486baa89ecbSBarry Smith 
487f1580f4eSBarry Smith   Note:
48895452b02SPatrick Sanan   This takes affect at the next construction of the preconditioner
489baa89ecbSBarry Smith 
490baa89ecbSBarry Smith   Level: intermediate
491baa89ecbSBarry Smith 
492f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiGetUseAbs()`
493baa89ecbSBarry Smith @*/
494d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg)
495d71ae5a4SJacob Faibussowitsch {
496baa89ecbSBarry Smith   PetscFunctionBegin;
497baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
498cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500baa89ecbSBarry Smith }
501baa89ecbSBarry Smith 
502baa89ecbSBarry Smith /*@
503f1580f4eSBarry Smith   PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the
50467ed0f3bSStefano Zampini   absolute values of the diagonal divisors in the preconditioner
505baa89ecbSBarry Smith 
506c3339decSBarry Smith   Logically Collective
507baa89ecbSBarry Smith 
508baa89ecbSBarry Smith   Input Parameter:
509cd47f5d9SBarry Smith . pc - the preconditioner context
510cd47f5d9SBarry Smith 
511baa89ecbSBarry Smith   Output Parameter:
512baa89ecbSBarry Smith . flg - whether to use absolute values or not
513cd47f5d9SBarry Smith 
514cd47f5d9SBarry Smith   Level: intermediate
515cd47f5d9SBarry Smith 
516f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
517cd47f5d9SBarry Smith @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg)
519d71ae5a4SJacob Faibussowitsch {
520cd47f5d9SBarry Smith   PetscFunctionBegin;
5210700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
522cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg));
5233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
524cd47f5d9SBarry Smith }
525cd47f5d9SBarry Smith 
5264b9ad928SBarry Smith /*@
527147403d9SBarry Smith   PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
52867ed0f3bSStefano Zampini 
529c3339decSBarry Smith   Logically Collective
53067ed0f3bSStefano Zampini 
53167ed0f3bSStefano Zampini   Input Parameters:
53267ed0f3bSStefano Zampini + pc  - the preconditioner context
53367ed0f3bSStefano Zampini - flg - the boolean flag
53467ed0f3bSStefano Zampini 
53567ed0f3bSStefano Zampini   Options Database Key:
536147403d9SBarry Smith . -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
53767ed0f3bSStefano Zampini 
538f1580f4eSBarry Smith   Note:
53967ed0f3bSStefano Zampini   This takes affect at the next construction of the preconditioner
54067ed0f3bSStefano Zampini 
54167ed0f3bSStefano Zampini   Level: intermediate
54267ed0f3bSStefano Zampini 
543f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()`
54467ed0f3bSStefano Zampini @*/
545d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg)
546d71ae5a4SJacob Faibussowitsch {
54767ed0f3bSStefano Zampini   PetscFunctionBegin;
54867ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
549cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg));
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55167ed0f3bSStefano Zampini }
55267ed0f3bSStefano Zampini 
55367ed0f3bSStefano Zampini /*@
554f1580f4eSBarry Smith   PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms
55567ed0f3bSStefano Zampini 
556c3339decSBarry Smith   Logically Collective
55767ed0f3bSStefano Zampini 
55867ed0f3bSStefano Zampini   Input Parameter:
55967ed0f3bSStefano Zampini . pc - the preconditioner context
56067ed0f3bSStefano Zampini 
56167ed0f3bSStefano Zampini   Output Parameter:
56267ed0f3bSStefano Zampini . flg - the boolean flag
56367ed0f3bSStefano Zampini 
56467ed0f3bSStefano Zampini   Options Database Key:
56567b8a455SSatish Balay . -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
56667ed0f3bSStefano Zampini 
56767ed0f3bSStefano Zampini   Level: intermediate
56867ed0f3bSStefano Zampini 
569f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
57067ed0f3bSStefano Zampini @*/
571d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg)
572d71ae5a4SJacob Faibussowitsch {
57367ed0f3bSStefano Zampini   PetscFunctionBegin;
57467ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
575cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg));
5763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57767ed0f3bSStefano Zampini }
57867ed0f3bSStefano Zampini 
57967ed0f3bSStefano Zampini /*@
580baa89ecbSBarry Smith   PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
581baa89ecbSBarry Smith   of the sum of rows entries for the diagonal preconditioner
5824b9ad928SBarry Smith 
583c3339decSBarry Smith   Logically Collective
5844b9ad928SBarry Smith 
5854b9ad928SBarry Smith   Input Parameters:
586baa89ecbSBarry Smith + pc   - the preconditioner context
587f1580f4eSBarry Smith - type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
5884b9ad928SBarry Smith 
5894b9ad928SBarry Smith   Options Database Key:
590147403d9SBarry Smith . -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
5914b9ad928SBarry Smith 
5924b9ad928SBarry Smith   Level: intermediate
5934b9ad928SBarry Smith 
594feefa0e1SJacob Faibussowitsch   Developer Notes:
595f1580f4eSBarry Smith   Why is there a separate function for using the absolute value?
596f1580f4eSBarry Smith 
597f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
5984b9ad928SBarry Smith @*/
599d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type)
600d71ae5a4SJacob Faibussowitsch {
6014b9ad928SBarry Smith   PetscFunctionBegin;
6020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
603cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type));
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6054b9ad928SBarry Smith }
6064b9ad928SBarry Smith 
60786697f06SMatthew Knepley /*@
608baa89ecbSBarry Smith   PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
60986697f06SMatthew Knepley 
610c3339decSBarry Smith   Not Collective
61186697f06SMatthew Knepley 
612baa89ecbSBarry Smith   Input Parameter:
61386697f06SMatthew Knepley . pc - the preconditioner context
61486697f06SMatthew Knepley 
615baa89ecbSBarry Smith   Output Parameter:
616f1580f4eSBarry Smith . type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
61786697f06SMatthew Knepley 
61886697f06SMatthew Knepley   Level: intermediate
61986697f06SMatthew Knepley 
620f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaUseAbs()`, `PCJacobiSetType()`
62186697f06SMatthew Knepley @*/
622d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type)
623d71ae5a4SJacob Faibussowitsch {
62486697f06SMatthew Knepley   PetscFunctionBegin;
6250700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
626cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62886697f06SMatthew Knepley }
629