xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
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 
138*f1580f4eSBarry 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));
1553bb1ff40SBarry Smith     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
1564b9ad928SBarry Smith   }
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1594b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1604b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1614b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1624b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1634b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1644b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1654b9ad928SBarry Smith   */
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith   /*
1684b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1694b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1704b9ad928SBarry Smith   */
1714b9ad928SBarry Smith   diag     = jac->diag;
1724b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith   if (diag) {
175b94d7dedSBarry Smith     PetscBool isset, isspd;
17667ed0f3bSStefano Zampini 
1774b9ad928SBarry Smith     if (jac->userowmax) {
1789566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL));
17986697f06SMatthew Knepley     } else if (jac->userowsum) {
1809566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diag));
1814b9ad928SBarry Smith     } else {
1829566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diag));
1834b9ad928SBarry Smith     }
1849566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
1851baa6e33SBarry Smith     if (jac->useabs) PetscCall(VecAbs(diag));
186b94d7dedSBarry Smith     PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
187b94d7dedSBarry Smith     if (jac->fixdiag && (!isset || !isspd)) {
1889566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(diag, &n));
1899566063dSJacob Faibussowitsch       PetscCall(VecGetArray(diag, &x));
1904b9ad928SBarry Smith       for (i = 0; i < n; i++) {
1914b9ad928SBarry Smith         if (x[i] == 0.0) {
1924b9ad928SBarry Smith           x[i]     = 1.0;
193cd47f5d9SBarry Smith           zeroflag = PETSC_TRUE;
1944b9ad928SBarry Smith         }
1954b9ad928SBarry Smith       }
1969566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(diag, &x));
1974b9ad928SBarry Smith     }
19867ed0f3bSStefano Zampini   }
1994b9ad928SBarry Smith   if (diagsqrt) {
2004b9ad928SBarry Smith     if (jac->userowmax) {
2019566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL));
20286697f06SMatthew Knepley     } else if (jac->userowsum) {
2039566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat, diagsqrt));
2044b9ad928SBarry Smith     } else {
2059566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat, diagsqrt));
2064b9ad928SBarry Smith     }
2079566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(diagsqrt, &n));
2089566063dSJacob Faibussowitsch     PetscCall(VecGetArray(diagsqrt, &x));
2094b9ad928SBarry Smith     for (i = 0; i < n; i++) {
2108f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i]));
2114b9ad928SBarry Smith       else {
2124b9ad928SBarry Smith         x[i]     = 1.0;
213cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2144b9ad928SBarry Smith       }
2154b9ad928SBarry Smith     }
2169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(diagsqrt, &x));
2174b9ad928SBarry Smith   }
21848a46eb9SPierre Jolivet   if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
2194b9ad928SBarry Smith   PetscFunctionReturn(0);
2204b9ad928SBarry Smith }
221*f1580f4eSBarry Smith 
2224b9ad928SBarry Smith /*
2234b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2244b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2254b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith    Input Parameter:
2284b9ad928SBarry Smith .  pc - the preconditioner context
2294b9ad928SBarry Smith */
2309371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) {
2314b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith   PetscFunctionBegin;
2349566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->diagsqrt));
2369566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2374b9ad928SBarry Smith   PetscFunctionReturn(0);
2384b9ad928SBarry Smith }
239*f1580f4eSBarry Smith 
2404b9ad928SBarry Smith /*
2414b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2424b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2434b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith    Input Parameter:
2464b9ad928SBarry Smith .  pc - the preconditioner context
2474b9ad928SBarry Smith */
2489371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) {
2494b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2504b9ad928SBarry Smith 
2514b9ad928SBarry Smith   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL));
2539566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->diag));
2549566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2554b9ad928SBarry Smith   PetscFunctionReturn(0);
2564b9ad928SBarry Smith }
257*f1580f4eSBarry Smith 
2584b9ad928SBarry Smith /*
2594b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2604b9ad928SBarry Smith 
2614b9ad928SBarry Smith    Input Parameters:
2624b9ad928SBarry Smith .  pc - the preconditioner context
2634b9ad928SBarry Smith .  x - input vector
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Output Parameter:
2664b9ad928SBarry Smith .  y - output vector
2674b9ad928SBarry Smith 
2684b9ad928SBarry Smith    Application Interface Routine: PCApply()
2694b9ad928SBarry Smith  */
2709371c9d4SSatish Balay static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y) {
2714b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith   PetscFunctionBegin;
27448a46eb9SPierre Jolivet   if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
2759566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diag));
2764b9ad928SBarry Smith   PetscFunctionReturn(0);
2774b9ad928SBarry Smith }
278*f1580f4eSBarry Smith 
2794b9ad928SBarry Smith /*
2804b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2814b9ad928SBarry Smith    symmetric preconditioner to a vector.
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith    Input Parameters:
2844b9ad928SBarry Smith .  pc - the preconditioner context
2854b9ad928SBarry Smith .  x - input vector
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith    Output Parameter:
2884b9ad928SBarry Smith .  y - output vector
2894b9ad928SBarry Smith 
2904b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
2914b9ad928SBarry Smith */
2929371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y) {
2934b9ad928SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
2944b9ad928SBarry Smith 
2954b9ad928SBarry Smith   PetscFunctionBegin;
29648a46eb9SPierre Jolivet   if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc));
2979566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y, x, jac->diagsqrt));
2984b9ad928SBarry Smith   PetscFunctionReturn(0);
2994b9ad928SBarry Smith }
30067ed0f3bSStefano Zampini 
3019371c9d4SSatish Balay static PetscErrorCode PCReset_Jacobi(PC pc) {
302a06653b4SBarry Smith   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
303a06653b4SBarry Smith 
304a06653b4SBarry Smith   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
3069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diagsqrt));
307a06653b4SBarry Smith   PetscFunctionReturn(0);
308a06653b4SBarry Smith }
309a06653b4SBarry Smith 
3104b9ad928SBarry Smith /*
3114b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3124b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith    Input Parameter:
3154b9ad928SBarry Smith .  pc - the preconditioner context
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3184b9ad928SBarry Smith */
3199371c9d4SSatish Balay static PetscErrorCode PCDestroy_Jacobi(PC pc) {
3204b9ad928SBarry Smith   PetscFunctionBegin;
3219566063dSJacob Faibussowitsch   PetscCall(PCReset_Jacobi(pc));
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL));
3239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL));
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL));
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL));
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL));
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL));
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith   /*
3304b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3314b9ad928SBarry Smith   */
3329566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3334b9ad928SBarry Smith   PetscFunctionReturn(0);
3344b9ad928SBarry Smith }
3354b9ad928SBarry Smith 
3369371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject) {
3374b9ad928SBarry Smith   PC_Jacobi   *jac = (PC_Jacobi *)pc->data;
338baa89ecbSBarry Smith   PetscBool    flg;
339baa89ecbSBarry Smith   PCJacobiType deflt, type;
3404b9ad928SBarry Smith 
3414b9ad928SBarry Smith   PetscFunctionBegin;
3429566063dSJacob Faibussowitsch   PetscCall(PCJacobiGetType(pc, &deflt));
343d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options");
3449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg));
3451baa6e33SBarry Smith   if (flg) PetscCall(PCJacobiSetType(pc, type));
3469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL));
3479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL));
348d0609cedSBarry Smith   PetscOptionsHeadEnd();
3494b9ad928SBarry Smith   PetscFunctionReturn(0);
3504b9ad928SBarry Smith }
3514b9ad928SBarry Smith 
3529371c9d4SSatish Balay static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer) {
353569e28a7SMatthew G. Knepley   PC_Jacobi *jac = (PC_Jacobi *)pc->data;
354569e28a7SMatthew G. Knepley   PetscBool  iascii;
355569e28a7SMatthew G. Knepley 
356569e28a7SMatthew G. Knepley   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
358569e28a7SMatthew G. Knepley   if (iascii) {
359569e28a7SMatthew G. Knepley     PCJacobiType      type;
36067ed0f3bSStefano Zampini     PetscBool         useAbs, fixdiag;
361569e28a7SMatthew G. Knepley     PetscViewerFormat format;
362569e28a7SMatthew G. Knepley 
3639566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetType(pc, &type));
3649566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
3659566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
3669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
3679566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
36848a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(VecView(jac->diag, viewer));
369569e28a7SMatthew G. Knepley   }
370569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
371569e28a7SMatthew G. Knepley }
372569e28a7SMatthew G. Knepley 
3734b9ad928SBarry Smith /*
3744b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
3754b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
3764b9ad928SBarry Smith    context, PC, that was created within PCCreate().
3774b9ad928SBarry Smith 
3784b9ad928SBarry Smith    Input Parameter:
3794b9ad928SBarry Smith .  pc - the preconditioner context
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith    Application Interface Routine: PCCreate()
3824b9ad928SBarry Smith */
3834b9ad928SBarry Smith 
3844b9ad928SBarry Smith /*MC
3855a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
3864b9ad928SBarry Smith 
387*f1580f4eSBarry Smith    Options Database Keys:
388967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
38967ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
390147403d9SBarry Smith -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Level: beginner
3934b9ad928SBarry Smith 
39495452b02SPatrick Sanan   Notes:
395*f1580f4eSBarry Smith     By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric
3964b9ad928SBarry Smith     can scale each side of the matrix by the square root of the diagonal entries.
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith     Zero entries along the diagonal are replaced with the value 1.0
3994b9ad928SBarry Smith 
400*f1580f4eSBarry Smith     See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks
401422a814eSBarry Smith 
402db781477SPatrick Sanan .seealso:  `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
403*f1580f4eSBarry Smith            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`,
404db781477SPatrick Sanan            `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
405db781477SPatrick Sanan            `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
4064b9ad928SBarry Smith M*/
4074b9ad928SBarry Smith 
4089371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) {
4094b9ad928SBarry Smith   PC_Jacobi *jac;
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith   PetscFunctionBegin;
4124b9ad928SBarry Smith   /*
4134b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4144b9ad928SBarry Smith      attach it to the PC object.
4154b9ad928SBarry Smith   */
4169566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &jac));
4174b9ad928SBarry Smith   pc->data = (void *)jac;
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith   /*
4204b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4214b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4224b9ad928SBarry Smith   */
4230a545947SLisandro Dalcin   jac->diag      = NULL;
4240a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4254b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
42686697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
427cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
42867ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith   /*
4314b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4324b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4334b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4344b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4354b9ad928SBarry Smith       not needed.
4364b9ad928SBarry Smith   */
4374b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4384b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4394b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
440a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4414b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4424b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
443569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4440a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4454b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4464b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4472fa5cd67SKarl Rupp 
4489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi));
4499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi));
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi));
4519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi));
4529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi));
4539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi));
4544b9ad928SBarry Smith   PetscFunctionReturn(0);
4554b9ad928SBarry Smith }
456cd47f5d9SBarry Smith 
457cd47f5d9SBarry Smith /*@
458*f1580f4eSBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the
45967ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
460cd47f5d9SBarry Smith 
461*f1580f4eSBarry Smith    Logically Collective on pc
462cd47f5d9SBarry Smith 
463cd47f5d9SBarry Smith    Input Parameters:
464baa89ecbSBarry Smith +  pc - the preconditioner context
465baa89ecbSBarry Smith -  flg - whether to use absolute values or not
466baa89ecbSBarry Smith 
467baa89ecbSBarry Smith    Options Database Key:
468147403d9SBarry Smith .  -pc_jacobi_abs <bool> - use absolute values
469baa89ecbSBarry Smith 
470*f1580f4eSBarry Smith    Note:
47195452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
472baa89ecbSBarry Smith 
473baa89ecbSBarry Smith    Level: intermediate
474baa89ecbSBarry Smith 
475*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiGetUseAbs()`
476baa89ecbSBarry Smith @*/
4779371c9d4SSatish Balay PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg) {
478baa89ecbSBarry Smith   PetscFunctionBegin;
479baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
480cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg));
481baa89ecbSBarry Smith   PetscFunctionReturn(0);
482baa89ecbSBarry Smith }
483baa89ecbSBarry Smith 
484baa89ecbSBarry Smith /*@
485*f1580f4eSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the
48667ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
487baa89ecbSBarry Smith 
488*f1580f4eSBarry Smith    Logically Collective on pc
489baa89ecbSBarry Smith 
490baa89ecbSBarry Smith    Input Parameter:
491cd47f5d9SBarry Smith .  pc - the preconditioner context
492cd47f5d9SBarry Smith 
493baa89ecbSBarry Smith    Output Parameter:
494baa89ecbSBarry Smith .  flg - whether to use absolute values or not
495cd47f5d9SBarry Smith 
496cd47f5d9SBarry Smith    Level: intermediate
497cd47f5d9SBarry Smith 
498*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
499cd47f5d9SBarry Smith @*/
5009371c9d4SSatish Balay PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg) {
501cd47f5d9SBarry Smith   PetscFunctionBegin;
5020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
503cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg));
504cd47f5d9SBarry Smith   PetscFunctionReturn(0);
505cd47f5d9SBarry Smith }
506cd47f5d9SBarry Smith 
5074b9ad928SBarry Smith /*@
508147403d9SBarry Smith    PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
50967ed0f3bSStefano Zampini 
510*f1580f4eSBarry Smith    Logically Collective on pc
51167ed0f3bSStefano Zampini 
51267ed0f3bSStefano Zampini    Input Parameters:
51367ed0f3bSStefano Zampini +  pc - the preconditioner context
51467ed0f3bSStefano Zampini -  flg - the boolean flag
51567ed0f3bSStefano Zampini 
51667ed0f3bSStefano Zampini    Options Database Key:
517147403d9SBarry Smith .  -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
51867ed0f3bSStefano Zampini 
519*f1580f4eSBarry Smith    Note:
52067ed0f3bSStefano Zampini    This takes affect at the next construction of the preconditioner
52167ed0f3bSStefano Zampini 
52267ed0f3bSStefano Zampini    Level: intermediate
52367ed0f3bSStefano Zampini 
524*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()`
52567ed0f3bSStefano Zampini @*/
5269371c9d4SSatish Balay PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg) {
52767ed0f3bSStefano Zampini   PetscFunctionBegin;
52867ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
529cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg));
53067ed0f3bSStefano Zampini   PetscFunctionReturn(0);
53167ed0f3bSStefano Zampini }
53267ed0f3bSStefano Zampini 
53367ed0f3bSStefano Zampini /*@
534*f1580f4eSBarry Smith    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms
53567ed0f3bSStefano Zampini 
536*f1580f4eSBarry Smith    Logically Collective on pc
53767ed0f3bSStefano Zampini 
53867ed0f3bSStefano Zampini    Input Parameter:
53967ed0f3bSStefano Zampini .  pc - the preconditioner context
54067ed0f3bSStefano Zampini 
54167ed0f3bSStefano Zampini    Output Parameter:
54267ed0f3bSStefano Zampini .  flg - the boolean flag
54367ed0f3bSStefano Zampini 
54467ed0f3bSStefano Zampini    Options Database Key:
54567b8a455SSatish Balay .  -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
54667ed0f3bSStefano Zampini 
54767ed0f3bSStefano Zampini    Level: intermediate
54867ed0f3bSStefano Zampini 
549*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
55067ed0f3bSStefano Zampini @*/
5519371c9d4SSatish Balay PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg) {
55267ed0f3bSStefano Zampini   PetscFunctionBegin;
55367ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
554cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg));
55567ed0f3bSStefano Zampini   PetscFunctionReturn(0);
55667ed0f3bSStefano Zampini }
55767ed0f3bSStefano Zampini 
55867ed0f3bSStefano Zampini /*@
559baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
560baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
5614b9ad928SBarry Smith 
562*f1580f4eSBarry Smith    Logically Collective on pc
5634b9ad928SBarry Smith 
5644b9ad928SBarry Smith    Input Parameters:
565baa89ecbSBarry Smith +  pc - the preconditioner context
566*f1580f4eSBarry Smith -  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
5674b9ad928SBarry Smith 
5684b9ad928SBarry Smith    Options Database Key:
569147403d9SBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
5704b9ad928SBarry Smith 
5714b9ad928SBarry Smith    Level: intermediate
5724b9ad928SBarry Smith 
573*f1580f4eSBarry Smith    Developer Note:
574*f1580f4eSBarry Smith    Why is there a separate function for using the absolute value?
575*f1580f4eSBarry Smith 
576*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
5774b9ad928SBarry Smith @*/
5789371c9d4SSatish Balay PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type) {
5794b9ad928SBarry Smith   PetscFunctionBegin;
5800700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
581cac4c232SBarry Smith   PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type));
5824b9ad928SBarry Smith   PetscFunctionReturn(0);
5834b9ad928SBarry Smith }
5844b9ad928SBarry Smith 
58586697f06SMatthew Knepley /*@
586baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
58786697f06SMatthew Knepley 
588*f1580f4eSBarry Smith    Not Collective on pc
58986697f06SMatthew Knepley 
590baa89ecbSBarry Smith    Input Parameter:
59186697f06SMatthew Knepley .  pc - the preconditioner context
59286697f06SMatthew Knepley 
593baa89ecbSBarry Smith    Output Parameter:
594*f1580f4eSBarry Smith .  type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
59586697f06SMatthew Knepley 
59686697f06SMatthew Knepley    Level: intermediate
59786697f06SMatthew Knepley 
598*f1580f4eSBarry Smith .seealso: `PCJACOBI`, `PCJacobiaUseAbs()`, `PCJacobiSetType()`
59986697f06SMatthew Knepley @*/
6009371c9d4SSatish Balay PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type) {
60186697f06SMatthew Knepley   PetscFunctionBegin;
6020700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
603cac4c232SBarry Smith   PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type));
60486697f06SMatthew Knepley   PetscFunctionReturn(0);
60586697f06SMatthew Knepley }
606