xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetType_Jacobi(PC pc, PCJacobiType type)
70d71ae5a4SJacob Faibussowitsch {
71baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
724b9ad928SBarry Smith 
734b9ad928SBarry Smith   PetscFunctionBegin;
74baa89ecbSBarry Smith   j->userowmax = PETSC_FALSE;
75baa89ecbSBarry Smith   j->userowsum = PETSC_FALSE;
76baa89ecbSBarry Smith   if (type == PC_JACOBI_ROWMAX) {
774b9ad928SBarry Smith     j->userowmax = PETSC_TRUE;
78baa89ecbSBarry Smith   } else if (type == PC_JACOBI_ROWSUM) {
79baa89ecbSBarry Smith     j->userowsum = PETSC_TRUE;
80baa89ecbSBarry Smith   }
814b9ad928SBarry Smith   PetscFunctionReturn(0);
824b9ad928SBarry Smith }
834b9ad928SBarry Smith 
84d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetType_Jacobi(PC pc, PCJacobiType *type)
85d71ae5a4SJacob Faibussowitsch {
86baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
8786697f06SMatthew Knepley 
8886697f06SMatthew Knepley   PetscFunctionBegin;
89baa89ecbSBarry Smith   if (j->userowmax) {
90baa89ecbSBarry Smith     *type = PC_JACOBI_ROWMAX;
91baa89ecbSBarry Smith   } else if (j->userowsum) {
92baa89ecbSBarry Smith     *type = PC_JACOBI_ROWSUM;
93baa89ecbSBarry Smith   } else {
94baa89ecbSBarry Smith     *type = PC_JACOBI_DIAGONAL;
95baa89ecbSBarry Smith   }
9686697f06SMatthew Knepley   PetscFunctionReturn(0);
9786697f06SMatthew Knepley }
9886697f06SMatthew Knepley 
99d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc, PetscBool flg)
100d71ae5a4SJacob Faibussowitsch {
101baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
102cd47f5d9SBarry Smith 
103cd47f5d9SBarry Smith   PetscFunctionBegin;
104baa89ecbSBarry Smith   j->useabs = flg;
105baa89ecbSBarry Smith   PetscFunctionReturn(0);
106baa89ecbSBarry Smith }
107baa89ecbSBarry Smith 
108d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc, PetscBool *flg)
109d71ae5a4SJacob Faibussowitsch {
110baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi *)pc->data;
111baa89ecbSBarry Smith 
112baa89ecbSBarry Smith   PetscFunctionBegin;
113baa89ecbSBarry Smith   *flg = j->useabs;
114cd47f5d9SBarry Smith   PetscFunctionReturn(0);
115cd47f5d9SBarry Smith }
116cd47f5d9SBarry Smith 
117d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetFixDiagonal_Jacobi(PC pc, PetscBool flg)
118d71ae5a4SJacob Faibussowitsch {
11967ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
12067ed0f3bSStefano Zampini 
12167ed0f3bSStefano Zampini   PetscFunctionBegin;
12267ed0f3bSStefano Zampini   j->fixdiag = flg;
12367ed0f3bSStefano Zampini   PetscFunctionReturn(0);
12467ed0f3bSStefano Zampini }
12567ed0f3bSStefano Zampini 
126d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetFixDiagonal_Jacobi(PC pc, PetscBool *flg)
127d71ae5a4SJacob Faibussowitsch {
12867ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi *)pc->data;
12967ed0f3bSStefano Zampini 
13067ed0f3bSStefano Zampini   PetscFunctionBegin;
13167ed0f3bSStefano Zampini   *flg = j->fixdiag;
13267ed0f3bSStefano Zampini   PetscFunctionReturn(0);
13367ed0f3bSStefano Zampini }
13467ed0f3bSStefano Zampini 
1354b9ad928SBarry Smith /*
1364b9ad928SBarry Smith    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1374b9ad928SBarry Smith                     by setting data structures and options.
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith    Input Parameter:
1404b9ad928SBarry Smith .  pc - the preconditioner context
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith    Application Interface Routine: PCSetUp()
1434b9ad928SBarry Smith 
144f1580f4eSBarry Smith    Note:
1454b9ad928SBarry Smith    The interface routine PCSetUp() is not usually called directly by
1464b9ad928SBarry Smith    the user, but instead is called by PCApply() if necessary.
1474b9ad928SBarry Smith */
148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi(PC pc)
149d71ae5a4SJacob Faibussowitsch {
1504b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
1514b9ad928SBarry Smith   Vec          diag, diagsqrt;
152cd47f5d9SBarry Smith   PetscInt     n, i;
1534b9ad928SBarry Smith   PetscScalar *x;
154ace3abfcSBarry Smith   PetscBool    zeroflag = PETSC_FALSE;
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith   PetscFunctionBegin;
1574b9ad928SBarry Smith   /*
1584b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1619566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->mat,&jac->diag));
1624b9ad928SBarry Smith   }
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1654b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1664b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1674b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1684b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1694b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1704b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1714b9ad928SBarry Smith   */
1724b9ad928SBarry Smith 
1734b9ad928SBarry Smith   /*
1744b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1754b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1764b9ad928SBarry Smith   */
1774b9ad928SBarry Smith   diag     = jac->diag;
1784b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith   if (diag) {
181b94d7dedSBarry Smith     PetscBool isset, isspd;
18267ed0f3bSStefano Zampini 
1834b9ad928SBarry Smith     if (jac->userowmax) {
1849566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL));
18586697f06SMatthew Knepley     } else if (jac->userowsum) {
1869566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diag));
1874b9ad928SBarry Smith     } else {
1889566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diag));
1894b9ad928SBarry Smith     }
1909566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
1911baa6e33SBarry Smith     if (jac->useabs) PetscCall(VecAbs(diag));
192b94d7dedSBarry Smith     PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
193b94d7dedSBarry Smith     if (jac->fixdiag && (!isset || !isspd)) {
1949566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(diag, &n));
1959566063dSJacob Faibussowitsch       PetscCall(VecGetArray(diag, &x));
1964b9ad928SBarry Smith       for (i = 0; i < n; i++) {
1974b9ad928SBarry Smith         if (x[i] == 0.0) {
1984b9ad928SBarry Smith           x[i]     = 1.0;
199cd47f5d9SBarry Smith           zeroflag = PETSC_TRUE;
2004b9ad928SBarry Smith         }
2014b9ad928SBarry Smith       }
2029566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(diag, &x));
2034b9ad928SBarry Smith     }
20467ed0f3bSStefano Zampini   }
2054b9ad928SBarry Smith   if (diagsqrt) {
2064b9ad928SBarry Smith     if (jac->userowmax) {
2079566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL));
20886697f06SMatthew Knepley     } else if (jac->userowsum) {
2099566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diagsqrt));
2104b9ad928SBarry Smith     } else {
2119566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diagsqrt));
2124b9ad928SBarry Smith     }
2139566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(diagsqrt, &n));
2149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(diagsqrt, &x));
2154b9ad928SBarry Smith     for (i = 0; i < n; i++) {
2168f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i]));
2174b9ad928SBarry Smith       else {
2184b9ad928SBarry Smith         x[i]     = 1.0;
219cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2204b9ad928SBarry Smith       }
2214b9ad928SBarry Smith     }
2229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(diagsqrt, &x));
2234b9ad928SBarry Smith   }
22448a46eb9SPierre Jolivet   if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
2254b9ad928SBarry Smith   PetscFunctionReturn(0);
2264b9ad928SBarry Smith }
227f1580f4eSBarry Smith 
2284b9ad928SBarry Smith /*
2294b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2304b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2314b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith    Input Parameter:
2344b9ad928SBarry Smith .  pc - the preconditioner context
2354b9ad928SBarry Smith */
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
237d71ae5a4SJacob Faibussowitsch {
2384b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   PetscFunctionBegin;
2419566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL));
2429566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2434b9ad928SBarry Smith   PetscFunctionReturn(0);
2444b9ad928SBarry Smith }
245f1580f4eSBarry Smith 
2464b9ad928SBarry Smith /*
2474b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2484b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2494b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2504b9ad928SBarry Smith 
2514b9ad928SBarry Smith    Input Parameter:
2524b9ad928SBarry Smith .  pc - the preconditioner context
2534b9ad928SBarry Smith */
254d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
255d71ae5a4SJacob Faibussowitsch {
2564b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL));
2609566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2614b9ad928SBarry Smith   PetscFunctionReturn(0);
2624b9ad928SBarry Smith }
263f1580f4eSBarry Smith 
2644b9ad928SBarry Smith /*
2654b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2664b9ad928SBarry Smith 
2674b9ad928SBarry Smith    Input Parameters:
2684b9ad928SBarry Smith .  pc - the preconditioner context
2694b9ad928SBarry Smith .  x - input vector
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Output Parameter:
2724b9ad928SBarry Smith .  y - output vector
2734b9ad928SBarry Smith 
2744b9ad928SBarry Smith    Application Interface Routine: PCApply()
2754b9ad928SBarry Smith  */
276d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y)
277d71ae5a4SJacob Faibussowitsch {
2784b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2794b9ad928SBarry Smith 
2804b9ad928SBarry Smith   PetscFunctionBegin;
28148a46eb9SPierre Jolivet   if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
2829566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diag));
2834b9ad928SBarry Smith   PetscFunctionReturn(0);
2844b9ad928SBarry Smith }
285f1580f4eSBarry Smith 
2864b9ad928SBarry Smith /*
2874b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2884b9ad928SBarry Smith    symmetric preconditioner to a vector.
2894b9ad928SBarry Smith 
2904b9ad928SBarry Smith    Input Parameters:
2914b9ad928SBarry Smith .  pc - the preconditioner context
2924b9ad928SBarry Smith .  x - input vector
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith    Output Parameter:
2954b9ad928SBarry Smith .  y - output vector
2964b9ad928SBarry Smith 
2974b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
2984b9ad928SBarry Smith */
299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y)
300d71ae5a4SJacob Faibussowitsch {
3014b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3024b9ad928SBarry Smith 
3034b9ad928SBarry Smith   PetscFunctionBegin;
30448a46eb9SPierre Jolivet   if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc));
3059566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diagsqrt));
3064b9ad928SBarry Smith   PetscFunctionReturn(0);
3074b9ad928SBarry Smith }
30867ed0f3bSStefano Zampini 
309d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Jacobi(PC pc)
310d71ae5a4SJacob Faibussowitsch {
311a06653b4SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
312a06653b4SBarry Smith 
313a06653b4SBarry Smith   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
3159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diagsqrt));
316a06653b4SBarry Smith   PetscFunctionReturn(0);
317a06653b4SBarry Smith }
318a06653b4SBarry Smith 
3194b9ad928SBarry Smith /*
3204b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3214b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith    Input Parameter:
3244b9ad928SBarry Smith .  pc - the preconditioner context
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3274b9ad928SBarry Smith */
328d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Jacobi(PC pc)
329d71ae5a4SJacob Faibussowitsch {
3304b9ad928SBarry Smith   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(PCReset_Jacobi(pc));
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL));
3339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL));
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL));
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL));
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith   /*
3404b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3414b9ad928SBarry Smith   */
3429566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3434b9ad928SBarry Smith   PetscFunctionReturn(0);
3444b9ad928SBarry Smith }
3454b9ad928SBarry Smith 
346d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject)
347d71ae5a4SJacob Faibussowitsch {
3484b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
349baa89ecbSBarry Smith   PetscBool    flg;
350baa89ecbSBarry Smith   PCJacobiType deflt, type;
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(PCJacobiGetType(pc, &deflt));
354d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options");
3559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg));
3561baa6e33SBarry Smith   if (flg) PetscCall(PCJacobiSetType(pc, type));
3579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL));
3589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL));
359d0609cedSBarry Smith   PetscOptionsHeadEnd();
3604b9ad928SBarry Smith   PetscFunctionReturn(0);
3614b9ad928SBarry Smith }
3624b9ad928SBarry Smith 
363d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
364d71ae5a4SJacob Faibussowitsch {
365569e28a7SMatthew G. Knepley   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
366569e28a7SMatthew G. Knepley   PetscBool  iascii;
367569e28a7SMatthew G. Knepley 
368569e28a7SMatthew G. Knepley   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
370569e28a7SMatthew G. Knepley   if (iascii) {
371569e28a7SMatthew G. Knepley     PCJacobiType      type;
37267ed0f3bSStefano Zampini     PetscBool         useAbs, fixdiag;
373569e28a7SMatthew G. Knepley     PetscViewerFormat format;
374569e28a7SMatthew G. Knepley 
3759566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetType(pc, &type));
3769566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
3779566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
3789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
3799566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
38048a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(VecView(jac->diag, viewer));
381569e28a7SMatthew G. Knepley   }
382569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
383569e28a7SMatthew G. Knepley }
384569e28a7SMatthew G. Knepley 
3854b9ad928SBarry Smith /*
3864b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
3874b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
3884b9ad928SBarry Smith    context, PC, that was created within PCCreate().
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith    Input Parameter:
3914b9ad928SBarry Smith .  pc - the preconditioner context
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Application Interface Routine: PCCreate()
3944b9ad928SBarry Smith */
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith /*MC
3975a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
3984b9ad928SBarry Smith 
399f1580f4eSBarry Smith    Options Database Keys:
400967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
40167ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
402147403d9SBarry Smith -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
4034b9ad928SBarry Smith 
4044b9ad928SBarry Smith    Level: beginner
4054b9ad928SBarry Smith 
40695452b02SPatrick Sanan   Notes:
407f1580f4eSBarry Smith     By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric
4084b9ad928SBarry Smith     can scale each side of the matrix by the square root of the diagonal entries.
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith     Zero entries along the diagonal are replaced with the value 1.0
4114b9ad928SBarry Smith 
412f1580f4eSBarry Smith     See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks
413422a814eSBarry Smith 
414db781477SPatrick Sanan .seealso:  `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
415f1580f4eSBarry Smith            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`,
416db781477SPatrick Sanan            `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
417db781477SPatrick Sanan            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
4184b9ad928SBarry Smith M*/
4194b9ad928SBarry Smith 
420d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
421d71ae5a4SJacob Faibussowitsch {
4224b9ad928SBarry Smith   PC_Jacobi *jac;
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith   PetscFunctionBegin;
4254b9ad928SBarry Smith   /*
4264b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4274b9ad928SBarry Smith      attach it to the PC object.
4284b9ad928SBarry Smith   */
4294dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
4304b9ad928SBarry Smith   pc->data = (void *)jac;
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith   /*
4334b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4344b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4354b9ad928SBarry Smith   */
4360a545947SLisandro Dalcin   jac->diag      = NULL;
4370a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4384b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
43986697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
440cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
44167ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4424b9ad928SBarry Smith 
4434b9ad928SBarry Smith   /*
4444b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4454b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4464b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4474b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4484b9ad928SBarry Smith       not needed.
4494b9ad928SBarry Smith   */
4504b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4514b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4524b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
453a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4544b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4554b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
456569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4570a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4584b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4594b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4602fa5cd67SKarl Rupp 
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi));
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi));
4639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi));
4649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi));
4659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi));
4669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi));
4674b9ad928SBarry Smith   PetscFunctionReturn(0);
4684b9ad928SBarry Smith }
469cd47f5d9SBarry Smith 
470cd47f5d9SBarry Smith /*@
471f1580f4eSBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the
47267ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
473cd47f5d9SBarry Smith 
474*c3339decSBarry Smith    Logically Collective
475cd47f5d9SBarry Smith 
476cd47f5d9SBarry Smith    Input Parameters:
477baa89ecbSBarry Smith +  pc - the preconditioner context
478baa89ecbSBarry Smith -  flg - whether to use absolute values or not
479baa89ecbSBarry Smith 
480baa89ecbSBarry Smith    Options Database Key:
481147403d9SBarry Smith .  -pc_jacobi_abs <bool> - use absolute values
482baa89ecbSBarry Smith 
483f1580f4eSBarry Smith    Note:
48495452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
485baa89ecbSBarry Smith 
486baa89ecbSBarry Smith    Level: intermediate
487baa89ecbSBarry Smith 
488f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiGetUseAbs()`
489baa89ecbSBarry Smith @*/
490d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg)
491d71ae5a4SJacob Faibussowitsch {
492baa89ecbSBarry Smith   PetscFunctionBegin;
493baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
494cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg));
495baa89ecbSBarry Smith   PetscFunctionReturn(0);
496baa89ecbSBarry Smith }
497baa89ecbSBarry Smith 
498baa89ecbSBarry Smith /*@
499f1580f4eSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the
50067ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
501baa89ecbSBarry Smith 
502*c3339decSBarry Smith    Logically Collective
503baa89ecbSBarry Smith 
504baa89ecbSBarry Smith    Input Parameter:
505cd47f5d9SBarry Smith .  pc - the preconditioner context
506cd47f5d9SBarry Smith 
507baa89ecbSBarry Smith    Output Parameter:
508baa89ecbSBarry Smith .  flg - whether to use absolute values or not
509cd47f5d9SBarry Smith 
510cd47f5d9SBarry Smith    Level: intermediate
511cd47f5d9SBarry Smith 
512f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
513cd47f5d9SBarry Smith @*/
514d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg)
515d71ae5a4SJacob Faibussowitsch {
516cd47f5d9SBarry Smith   PetscFunctionBegin;
5170700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
518cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg));
519cd47f5d9SBarry Smith   PetscFunctionReturn(0);
520cd47f5d9SBarry Smith }
521cd47f5d9SBarry Smith 
5224b9ad928SBarry Smith /*@
523147403d9SBarry Smith    PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
52467ed0f3bSStefano Zampini 
525*c3339decSBarry Smith    Logically Collective
52667ed0f3bSStefano Zampini 
52767ed0f3bSStefano Zampini    Input Parameters:
52867ed0f3bSStefano Zampini +  pc - the preconditioner context
52967ed0f3bSStefano Zampini -  flg - the boolean flag
53067ed0f3bSStefano Zampini 
53167ed0f3bSStefano Zampini    Options Database Key:
532147403d9SBarry Smith .  -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
53367ed0f3bSStefano Zampini 
534f1580f4eSBarry Smith    Note:
53567ed0f3bSStefano Zampini    This takes affect at the next construction of the preconditioner
53667ed0f3bSStefano Zampini 
53767ed0f3bSStefano Zampini    Level: intermediate
53867ed0f3bSStefano Zampini 
539f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()`
54067ed0f3bSStefano Zampini @*/
541d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg)
542d71ae5a4SJacob Faibussowitsch {
54367ed0f3bSStefano Zampini   PetscFunctionBegin;
54467ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
545cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg));
54667ed0f3bSStefano Zampini   PetscFunctionReturn(0);
54767ed0f3bSStefano Zampini }
54867ed0f3bSStefano Zampini 
54967ed0f3bSStefano Zampini /*@
550f1580f4eSBarry Smith    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms
55167ed0f3bSStefano Zampini 
552*c3339decSBarry Smith    Logically Collective
55367ed0f3bSStefano Zampini 
55467ed0f3bSStefano Zampini    Input Parameter:
55567ed0f3bSStefano Zampini .  pc - the preconditioner context
55667ed0f3bSStefano Zampini 
55767ed0f3bSStefano Zampini    Output Parameter:
55867ed0f3bSStefano Zampini .  flg - the boolean flag
55967ed0f3bSStefano Zampini 
56067ed0f3bSStefano Zampini    Options Database Key:
56167b8a455SSatish Balay .  -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
56267ed0f3bSStefano Zampini 
56367ed0f3bSStefano Zampini    Level: intermediate
56467ed0f3bSStefano Zampini 
565f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
56667ed0f3bSStefano Zampini @*/
567d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg)
568d71ae5a4SJacob Faibussowitsch {
56967ed0f3bSStefano Zampini   PetscFunctionBegin;
57067ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
571cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg));
57267ed0f3bSStefano Zampini   PetscFunctionReturn(0);
57367ed0f3bSStefano Zampini }
57467ed0f3bSStefano Zampini 
57567ed0f3bSStefano Zampini /*@
576baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
577baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
5784b9ad928SBarry Smith 
579*c3339decSBarry Smith    Logically Collective
5804b9ad928SBarry Smith 
5814b9ad928SBarry Smith    Input Parameters:
582baa89ecbSBarry Smith +  pc - the preconditioner context
583f1580f4eSBarry Smith -  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
5844b9ad928SBarry Smith 
5854b9ad928SBarry Smith    Options Database Key:
586147403d9SBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
5874b9ad928SBarry Smith 
5884b9ad928SBarry Smith    Level: intermediate
5894b9ad928SBarry Smith 
590f1580f4eSBarry Smith    Developer Note:
591f1580f4eSBarry Smith    Why is there a separate function for using the absolute value?
592f1580f4eSBarry Smith 
593f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
5944b9ad928SBarry Smith @*/
595d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type)
596d71ae5a4SJacob Faibussowitsch {
5974b9ad928SBarry Smith   PetscFunctionBegin;
5980700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
599cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type));
6004b9ad928SBarry Smith   PetscFunctionReturn(0);
6014b9ad928SBarry Smith }
6024b9ad928SBarry Smith 
60386697f06SMatthew Knepley /*@
604baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
60586697f06SMatthew Knepley 
606*c3339decSBarry Smith    Not Collective
60786697f06SMatthew Knepley 
608baa89ecbSBarry Smith    Input Parameter:
60986697f06SMatthew Knepley .  pc - the preconditioner context
61086697f06SMatthew Knepley 
611baa89ecbSBarry Smith    Output Parameter:
612f1580f4eSBarry Smith .  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
61386697f06SMatthew Knepley 
61486697f06SMatthew Knepley    Level: intermediate
61586697f06SMatthew Knepley 
616f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaUseAbs()`, `PCJacobiSetType()`
61786697f06SMatthew Knepley @*/
618d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type)
619d71ae5a4SJacob Faibussowitsch {
62086697f06SMatthew Knepley   PetscFunctionBegin;
6210700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
622cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type));
62386697f06SMatthew Knepley   PetscFunctionReturn(0);
62486697f06SMatthew Knepley }
625