xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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 
69baa89ecbSBarry Smith static PetscErrorCode  PCJacobiSetType_Jacobi(PC pc,PCJacobiType type)
704b9ad928SBarry Smith {
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 
84baa89ecbSBarry Smith static PetscErrorCode  PCJacobiGetType_Jacobi(PC pc,PCJacobiType *type)
8586697f06SMatthew Knepley {
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 
99baa89ecbSBarry Smith static PetscErrorCode  PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg)
100cd47f5d9SBarry Smith {
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 
108baa89ecbSBarry Smith static PetscErrorCode  PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool *flg)
109baa89ecbSBarry Smith {
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 
11767ed0f3bSStefano Zampini static PetscErrorCode  PCJacobiSetFixDiagonal_Jacobi(PC pc,PetscBool flg)
11867ed0f3bSStefano Zampini {
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 
12667ed0f3bSStefano Zampini static PetscErrorCode  PCJacobiGetFixDiagonal_Jacobi(PC pc,PetscBool *flg)
12767ed0f3bSStefano Zampini {
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 /*
1374b9ad928SBarry Smith    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1384b9ad928SBarry Smith                     by setting data structures and options.
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith    Input Parameter:
1414b9ad928SBarry Smith .  pc - the preconditioner context
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith    Application Interface Routine: PCSetUp()
1444b9ad928SBarry Smith 
1454b9ad928SBarry Smith    Notes:
1464b9ad928SBarry Smith    The interface routine PCSetUp() is not usually called directly by
1474b9ad928SBarry Smith    the user, but instead is called by PCApply() if necessary.
1484b9ad928SBarry Smith */
1496849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc)
1504b9ad928SBarry Smith {
1514b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
1524b9ad928SBarry Smith   Vec            diag,diagsqrt;
153cd47f5d9SBarry Smith   PetscInt       n,i;
1544b9ad928SBarry Smith   PetscScalar    *x;
155ace3abfcSBarry Smith   PetscBool      zeroflag = PETSC_FALSE;
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith   PetscFunctionBegin;
1584b9ad928SBarry Smith   /*
1594b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1629566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->mat,&jac->diag));
1633bb1ff40SBarry Smith     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
1644b9ad928SBarry Smith   }
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1674b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1684b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1694b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1704b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1714b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1724b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1734b9ad928SBarry Smith   */
1744b9ad928SBarry Smith 
1754b9ad928SBarry Smith   /*
1764b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1774b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1784b9ad928SBarry Smith   */
1794b9ad928SBarry Smith   diag     = jac->diag;
1804b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1814b9ad928SBarry Smith 
1824b9ad928SBarry Smith   if (diag) {
18367ed0f3bSStefano Zampini     PetscBool isspd;
18467ed0f3bSStefano Zampini 
1854b9ad928SBarry Smith     if (jac->userowmax) {
1869566063dSJacob Faibussowitsch       PetscCall(MatGetRowMaxAbs(pc->pmat,diag,NULL));
18786697f06SMatthew Knepley     } else if (jac->userowsum) {
1889566063dSJacob Faibussowitsch       PetscCall(MatGetRowSum(pc->pmat,diag));
1894b9ad928SBarry Smith     } else {
1909566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(pc->pmat,diag));
1914b9ad928SBarry Smith     }
1929566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
193cd47f5d9SBarry Smith     if (jac->useabs) {
1949566063dSJacob Faibussowitsch       PetscCall(VecAbs(diag));
195cd47f5d9SBarry Smith     }
1969566063dSJacob Faibussowitsch     PetscCall(MatGetOption(pc->pmat,MAT_SPD,&isspd));
19767ed0f3bSStefano Zampini     if (jac->fixdiag && !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   }
2284b9ad928SBarry Smith   if (zeroflag) {
2299566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n"));
2304b9ad928SBarry Smith   }
2314b9ad928SBarry Smith   PetscFunctionReturn(0);
2324b9ad928SBarry Smith }
2334b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2344b9ad928SBarry Smith /*
2354b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2364b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2374b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2384b9ad928SBarry Smith 
2394b9ad928SBarry Smith    Input Parameter:
2404b9ad928SBarry Smith .  pc - the preconditioner context
2414b9ad928SBarry Smith */
2426849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
2434b9ad928SBarry Smith {
2444b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2454b9ad928SBarry Smith 
2464b9ad928SBarry Smith   PetscFunctionBegin;
2479566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat,&jac->diagsqrt,NULL));
2489566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt));
2499566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2504b9ad928SBarry Smith   PetscFunctionReturn(0);
2514b9ad928SBarry Smith }
2524b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2534b9ad928SBarry Smith /*
2544b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2554b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2564b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith    Input Parameter:
2594b9ad928SBarry Smith .  pc - the preconditioner context
2604b9ad928SBarry Smith */
2616849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
2624b9ad928SBarry Smith {
2634b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(pc->pmat,&jac->diag,NULL));
2679566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag));
2689566063dSJacob Faibussowitsch   PetscCall(PCSetUp_Jacobi(pc));
2694b9ad928SBarry Smith   PetscFunctionReturn(0);
2704b9ad928SBarry Smith }
2714b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2724b9ad928SBarry Smith /*
2734b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith    Input Parameters:
2764b9ad928SBarry Smith .  pc - the preconditioner context
2774b9ad928SBarry Smith .  x - input vector
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith    Output Parameter:
2804b9ad928SBarry Smith .  y - output vector
2814b9ad928SBarry Smith 
2824b9ad928SBarry Smith    Application Interface Routine: PCApply()
2834b9ad928SBarry Smith  */
2846849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y)
2854b9ad928SBarry Smith {
2864b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith   PetscFunctionBegin;
2894b9ad928SBarry Smith   if (!jac->diag) {
2909566063dSJacob Faibussowitsch     PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
2914b9ad928SBarry Smith   }
2929566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y,x,jac->diag));
2934b9ad928SBarry Smith   PetscFunctionReturn(0);
2944b9ad928SBarry Smith }
2954b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2964b9ad928SBarry Smith /*
2974b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2984b9ad928SBarry Smith    symmetric preconditioner to a vector.
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith    Input Parameters:
3014b9ad928SBarry Smith .  pc - the preconditioner context
3024b9ad928SBarry Smith .  x - input vector
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith    Output Parameter:
3054b9ad928SBarry Smith .  y - output vector
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
3084b9ad928SBarry Smith */
3096849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)
3104b9ad928SBarry Smith {
3114b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
3124b9ad928SBarry Smith 
3134b9ad928SBarry Smith   PetscFunctionBegin;
3144b9ad928SBarry Smith   if (!jac->diagsqrt) {
3159566063dSJacob Faibussowitsch     PetscCall(PCSetUp_Jacobi_Symmetric(pc));
3164b9ad928SBarry Smith   }
3179566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(y,x,jac->diagsqrt));
3184b9ad928SBarry Smith   PetscFunctionReturn(0);
3194b9ad928SBarry Smith }
32067ed0f3bSStefano Zampini 
3214b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
322a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc)
323a06653b4SBarry Smith {
324a06653b4SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
325a06653b4SBarry Smith 
326a06653b4SBarry Smith   PetscFunctionBegin;
3279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diag));
3289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->diagsqrt));
329a06653b4SBarry Smith   PetscFunctionReturn(0);
330a06653b4SBarry Smith }
331a06653b4SBarry Smith 
3324b9ad928SBarry Smith /*
3334b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3344b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3354b9ad928SBarry Smith 
3364b9ad928SBarry Smith    Input Parameter:
3374b9ad928SBarry Smith .  pc - the preconditioner context
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3404b9ad928SBarry Smith */
3416849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc)
3424b9ad928SBarry Smith {
3434b9ad928SBarry Smith   PetscFunctionBegin;
3449566063dSJacob Faibussowitsch   PetscCall(PCReset_Jacobi(pc));
3459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",NULL));
3469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",NULL));
3479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",NULL));
3489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",NULL));
3499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",NULL));
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",NULL));
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith   /*
3534b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3544b9ad928SBarry Smith   */
3559566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
3564b9ad928SBarry Smith   PetscFunctionReturn(0);
3574b9ad928SBarry Smith }
3584b9ad928SBarry Smith 
3594416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PetscOptionItems *PetscOptionsObject,PC pc)
3604b9ad928SBarry Smith {
3614b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
362baa89ecbSBarry Smith   PetscBool      flg;
363baa89ecbSBarry Smith   PCJacobiType   deflt,type;
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith   PetscFunctionBegin;
3669566063dSJacob Faibussowitsch   PetscCall(PCJacobiGetType(pc,&deflt));
3679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Jacobi options"));
3689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg));
369baa89ecbSBarry Smith   if (flg) {
3709566063dSJacob Faibussowitsch     PetscCall(PCJacobiSetType(pc,type));
371baa89ecbSBarry Smith   }
3729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL));
3739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal","Fix null terms on diagonal","PCJacobiSetFixDiagonal",jac->fixdiag,&jac->fixdiag,NULL));
3749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
3754b9ad928SBarry Smith   PetscFunctionReturn(0);
3764b9ad928SBarry Smith }
3774b9ad928SBarry Smith 
378569e28a7SMatthew G. Knepley static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
379569e28a7SMatthew G. Knepley {
380569e28a7SMatthew G. Knepley   PC_Jacobi     *jac = (PC_Jacobi *) pc->data;
381569e28a7SMatthew G. Knepley   PetscBool      iascii;
382569e28a7SMatthew G. Knepley 
383569e28a7SMatthew G. Knepley   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
385569e28a7SMatthew G. Knepley   if (iascii) {
386569e28a7SMatthew G. Knepley     PCJacobiType      type;
38767ed0f3bSStefano Zampini     PetscBool         useAbs,fixdiag;
388569e28a7SMatthew G. Knepley     PetscViewerFormat format;
389569e28a7SMatthew G. Knepley 
3909566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetType(pc, &type));
3919566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
3929566063dSJacob Faibussowitsch     PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
3939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
3949566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
395569e28a7SMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
3969566063dSJacob Faibussowitsch       PetscCall(VecView(jac->diag, viewer));
397569e28a7SMatthew G. Knepley     }
398569e28a7SMatthew G. Knepley   }
399569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
400569e28a7SMatthew G. Knepley }
401569e28a7SMatthew G. Knepley 
4024b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
4034b9ad928SBarry Smith /*
4044b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
4054b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
4064b9ad928SBarry Smith    context, PC, that was created within PCCreate().
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith    Input Parameter:
4094b9ad928SBarry Smith .  pc - the preconditioner context
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith    Application Interface Routine: PCCreate()
4124b9ad928SBarry Smith */
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith /*MC
4155a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
4164b9ad928SBarry Smith 
4174b9ad928SBarry Smith    Options Database Key:
418967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
41967ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
420147403d9SBarry Smith -    -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith    Level: beginner
4234b9ad928SBarry Smith 
42495452b02SPatrick Sanan   Notes:
42595452b02SPatrick Sanan     By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric
4264b9ad928SBarry Smith     can scale each side of the matrix by the square root of the diagonal entries.
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith     Zero entries along the diagonal are replaced with the value 1.0
4294b9ad928SBarry Smith 
430422a814eSBarry Smith     See PCPBJACOBI for a point-block Jacobi preconditioner
431422a814eSBarry Smith 
4324b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
43367ed0f3bSStefano Zampini            PCJacobiSetType(), PCJacobiSetUseAbs(), PCJacobiGetUseAbs(),
43467ed0f3bSStefano Zampini            PCJacobiSetFixDiagonal(), PCJacobiGetFixDiagonal(), PCPBJACOBI
4354b9ad928SBarry Smith M*/
4364b9ad928SBarry Smith 
4378cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
4384b9ad928SBarry Smith {
4394b9ad928SBarry Smith   PC_Jacobi      *jac;
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith   PetscFunctionBegin;
4424b9ad928SBarry Smith   /*
4434b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4444b9ad928SBarry Smith      attach it to the PC object.
4454b9ad928SBarry Smith   */
4469566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&jac));
4474b9ad928SBarry Smith   pc->data = (void*)jac;
4484b9ad928SBarry Smith 
4494b9ad928SBarry Smith   /*
4504b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4514b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4524b9ad928SBarry Smith   */
4530a545947SLisandro Dalcin   jac->diag      = NULL;
4540a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4554b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
45686697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
457cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
45867ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith   /*
4614b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4624b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4634b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4644b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4654b9ad928SBarry Smith       not needed.
4664b9ad928SBarry Smith   */
4674b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4684b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4694b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
470a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4714b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4724b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
473569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4740a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4754b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4764b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4772fa5cd67SKarl Rupp 
4789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi));
4799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi));
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi));
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi));
4829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",PCJacobiSetFixDiagonal_Jacobi));
4839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",PCJacobiGetFixDiagonal_Jacobi));
4844b9ad928SBarry Smith   PetscFunctionReturn(0);
4854b9ad928SBarry Smith }
486cd47f5d9SBarry Smith 
487cd47f5d9SBarry Smith /*@
488cd47f5d9SBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the
48967ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
490cd47f5d9SBarry Smith 
491ad4df100SBarry Smith    Logically Collective on PC
492cd47f5d9SBarry Smith 
493cd47f5d9SBarry Smith    Input Parameters:
494baa89ecbSBarry Smith +  pc - the preconditioner context
495baa89ecbSBarry Smith -  flg - whether to use absolute values or not
496baa89ecbSBarry Smith 
497baa89ecbSBarry Smith    Options Database Key:
498147403d9SBarry Smith .  -pc_jacobi_abs <bool> - use absolute values
499baa89ecbSBarry Smith 
50095452b02SPatrick Sanan    Notes:
50195452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
502baa89ecbSBarry Smith 
503baa89ecbSBarry Smith    Level: intermediate
504baa89ecbSBarry Smith 
505baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiGetUseAbs()
506baa89ecbSBarry Smith 
507baa89ecbSBarry Smith @*/
508baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetUseAbs(PC pc,PetscBool flg)
509baa89ecbSBarry Smith {
510baa89ecbSBarry Smith   PetscFunctionBegin;
511baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
512*cac4c232SBarry Smith   PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));
513baa89ecbSBarry Smith   PetscFunctionReturn(0);
514baa89ecbSBarry Smith }
515baa89ecbSBarry Smith 
516baa89ecbSBarry Smith /*@
517baa89ecbSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the
51867ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
519baa89ecbSBarry Smith 
520baa89ecbSBarry Smith    Logically Collective on PC
521baa89ecbSBarry Smith 
522baa89ecbSBarry Smith    Input Parameter:
523cd47f5d9SBarry Smith .  pc - the preconditioner context
524cd47f5d9SBarry Smith 
525baa89ecbSBarry Smith    Output Parameter:
526baa89ecbSBarry Smith .  flg - whether to use absolute values or not
527cd47f5d9SBarry Smith 
528cd47f5d9SBarry Smith    Level: intermediate
529cd47f5d9SBarry Smith 
530baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiSetUseAbs(), PCJacobiGetType()
531cd47f5d9SBarry Smith 
532cd47f5d9SBarry Smith @*/
533baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetUseAbs(PC pc,PetscBool *flg)
534cd47f5d9SBarry Smith {
535cd47f5d9SBarry Smith   PetscFunctionBegin;
5360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
537*cac4c232SBarry Smith   PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));
538cd47f5d9SBarry Smith   PetscFunctionReturn(0);
539cd47f5d9SBarry Smith }
540cd47f5d9SBarry Smith 
5414b9ad928SBarry Smith /*@
542147403d9SBarry Smith    PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
54367ed0f3bSStefano Zampini 
54467ed0f3bSStefano Zampini    Logically Collective on PC
54567ed0f3bSStefano Zampini 
54667ed0f3bSStefano Zampini    Input Parameters:
54767ed0f3bSStefano Zampini +  pc - the preconditioner context
54867ed0f3bSStefano Zampini -  flg - the boolean flag
54967ed0f3bSStefano Zampini 
55067ed0f3bSStefano Zampini    Options Database Key:
551147403d9SBarry Smith .  -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
55267ed0f3bSStefano Zampini 
55367ed0f3bSStefano Zampini    Notes:
55467ed0f3bSStefano Zampini     This takes affect at the next construction of the preconditioner
55567ed0f3bSStefano Zampini 
55667ed0f3bSStefano Zampini    Level: intermediate
55767ed0f3bSStefano Zampini 
55867ed0f3bSStefano Zampini .seealso: PCJacobiSetType(), PCJacobiGetFixDiagonal()
55967ed0f3bSStefano Zampini 
56067ed0f3bSStefano Zampini @*/
56167ed0f3bSStefano Zampini PetscErrorCode  PCJacobiSetFixDiagonal(PC pc,PetscBool flg)
56267ed0f3bSStefano Zampini {
56367ed0f3bSStefano Zampini   PetscFunctionBegin;
56467ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
565*cac4c232SBarry Smith   PetscTryMethod(pc,"PCJacobiSetFixDiagonal_C",(PC,PetscBool),(pc,flg));
56667ed0f3bSStefano Zampini   PetscFunctionReturn(0);
56767ed0f3bSStefano Zampini }
56867ed0f3bSStefano Zampini 
56967ed0f3bSStefano Zampini /*@
57067ed0f3bSStefano Zampini    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner checks for zero diagonal terms
57167ed0f3bSStefano Zampini 
57267ed0f3bSStefano Zampini    Logically Collective on PC
57367ed0f3bSStefano Zampini 
57467ed0f3bSStefano Zampini    Input Parameter:
57567ed0f3bSStefano Zampini .  pc - the preconditioner context
57667ed0f3bSStefano Zampini 
57767ed0f3bSStefano Zampini    Output Parameter:
57867ed0f3bSStefano Zampini .  flg - the boolean flag
57967ed0f3bSStefano Zampini 
58067ed0f3bSStefano Zampini    Options Database Key:
58167b8a455SSatish Balay .  -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
58267ed0f3bSStefano Zampini 
58367ed0f3bSStefano Zampini    Level: intermediate
58467ed0f3bSStefano Zampini 
58567ed0f3bSStefano Zampini .seealso: PCJacobiSetType(), PCJacobiSetFixDiagonal()
58667ed0f3bSStefano Zampini 
58767ed0f3bSStefano Zampini @*/
58867ed0f3bSStefano Zampini PetscErrorCode  PCJacobiGetFixDiagonal(PC pc,PetscBool *flg)
58967ed0f3bSStefano Zampini {
59067ed0f3bSStefano Zampini   PetscFunctionBegin;
59167ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
592*cac4c232SBarry Smith   PetscUseMethod(pc,"PCJacobiGetFixDiagonal_C",(PC,PetscBool*),(pc,flg));
59367ed0f3bSStefano Zampini   PetscFunctionReturn(0);
59467ed0f3bSStefano Zampini }
59567ed0f3bSStefano Zampini 
59667ed0f3bSStefano Zampini /*@
597baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
598baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
5994b9ad928SBarry Smith 
600ad4df100SBarry Smith    Logically Collective on PC
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith    Input Parameters:
603baa89ecbSBarry Smith +  pc - the preconditioner context
604baa89ecbSBarry Smith -  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
6054b9ad928SBarry Smith 
6064b9ad928SBarry Smith    Options Database Key:
607147403d9SBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith    Level: intermediate
6104b9ad928SBarry Smith 
611baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiGetType()
6124b9ad928SBarry Smith @*/
613baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetType(PC pc,PCJacobiType type)
6144b9ad928SBarry Smith {
6154b9ad928SBarry Smith   PetscFunctionBegin;
6160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
617*cac4c232SBarry Smith   PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));
6184b9ad928SBarry Smith   PetscFunctionReturn(0);
6194b9ad928SBarry Smith }
6204b9ad928SBarry Smith 
62186697f06SMatthew Knepley /*@
622baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
62386697f06SMatthew Knepley 
624baa89ecbSBarry Smith    Not Collective on PC
62586697f06SMatthew Knepley 
626baa89ecbSBarry Smith    Input Parameter:
62786697f06SMatthew Knepley .  pc - the preconditioner context
62886697f06SMatthew Knepley 
629baa89ecbSBarry Smith    Output Parameter:
630baa89ecbSBarry Smith .  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
63186697f06SMatthew Knepley 
63286697f06SMatthew Knepley    Level: intermediate
63386697f06SMatthew Knepley 
634baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiSetType()
63586697f06SMatthew Knepley @*/
636baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetType(PC pc,PCJacobiType *type)
63786697f06SMatthew Knepley {
63886697f06SMatthew Knepley   PetscFunctionBegin;
6390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
640*cac4c232SBarry Smith   PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));
64186697f06SMatthew Knepley   PetscFunctionReturn(0);
64286697f06SMatthew Knepley }
643