xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision 569e28a7502c934530efa6903dfdf087ebd3ebcd)
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 */
664b9ad928SBarry Smith } PC_Jacobi;
674b9ad928SBarry Smith 
68baa89ecbSBarry Smith static PetscErrorCode  PCJacobiSetType_Jacobi(PC pc,PCJacobiType type)
694b9ad928SBarry Smith {
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 
83baa89ecbSBarry Smith static PetscErrorCode  PCJacobiGetType_Jacobi(PC pc,PCJacobiType *type)
8486697f06SMatthew Knepley {
85baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi*)pc->data;
8686697f06SMatthew Knepley 
8786697f06SMatthew Knepley   PetscFunctionBegin;
88baa89ecbSBarry Smith   if (j->userowmax) {
89baa89ecbSBarry Smith     *type = PC_JACOBI_ROWMAX;
90baa89ecbSBarry Smith   } else if (j->userowsum) {
91baa89ecbSBarry Smith     *type = PC_JACOBI_ROWSUM;
92baa89ecbSBarry Smith   } else {
93baa89ecbSBarry Smith     *type = PC_JACOBI_DIAGONAL;
94baa89ecbSBarry Smith   }
9586697f06SMatthew Knepley   PetscFunctionReturn(0);
9686697f06SMatthew Knepley }
9786697f06SMatthew Knepley 
98baa89ecbSBarry Smith static PetscErrorCode  PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg)
99cd47f5d9SBarry Smith {
100baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi*)pc->data;
101cd47f5d9SBarry Smith 
102cd47f5d9SBarry Smith   PetscFunctionBegin;
103baa89ecbSBarry Smith   j->useabs = flg;
104baa89ecbSBarry Smith   PetscFunctionReturn(0);
105baa89ecbSBarry Smith }
106baa89ecbSBarry Smith 
107baa89ecbSBarry Smith static PetscErrorCode  PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool *flg)
108baa89ecbSBarry Smith {
109baa89ecbSBarry Smith   PC_Jacobi *j = (PC_Jacobi*)pc->data;
110baa89ecbSBarry Smith 
111baa89ecbSBarry Smith   PetscFunctionBegin;
112baa89ecbSBarry Smith   *flg = j->useabs;
113cd47f5d9SBarry Smith   PetscFunctionReturn(0);
114cd47f5d9SBarry Smith }
115cd47f5d9SBarry Smith 
1164b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1174b9ad928SBarry Smith /*
1184b9ad928SBarry Smith    PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1194b9ad928SBarry Smith                     by setting data structures and options.
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith    Input Parameter:
1224b9ad928SBarry Smith .  pc - the preconditioner context
1234b9ad928SBarry Smith 
1244b9ad928SBarry Smith    Application Interface Routine: PCSetUp()
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith    Notes:
1274b9ad928SBarry Smith    The interface routine PCSetUp() is not usually called directly by
1284b9ad928SBarry Smith    the user, but instead is called by PCApply() if necessary.
1294b9ad928SBarry Smith */
1306849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc)
1314b9ad928SBarry Smith {
1324b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
1334b9ad928SBarry Smith   Vec            diag,diagsqrt;
134dfbe8321SBarry Smith   PetscErrorCode ierr;
135cd47f5d9SBarry Smith   PetscInt       n,i;
1364b9ad928SBarry Smith   PetscScalar    *x;
137ace3abfcSBarry Smith   PetscBool      zeroflag = PETSC_FALSE;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
1404b9ad928SBarry Smith   /*
1414b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1442a7a6963SBarry Smith     ierr = MatCreateVecs(pc->mat,&jac->diag);CHKERRQ(ierr);
1453bb1ff40SBarry Smith     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
1464b9ad928SBarry Smith   }
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1494b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1504b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1514b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1524b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1534b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1544b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1554b9ad928SBarry Smith   */
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith   /*
1584b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1594b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1604b9ad928SBarry Smith   */
1614b9ad928SBarry Smith   diag     = jac->diag;
1624b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith   if (diag) {
1654b9ad928SBarry Smith     if (jac->userowmax) {
1660298fd71SBarry Smith       ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr);
16786697f06SMatthew Knepley     } else if (jac->userowsum) {
16886697f06SMatthew Knepley       ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr);
1694b9ad928SBarry Smith     } else {
1704b9ad928SBarry Smith       ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr);
1714b9ad928SBarry Smith     }
1724b9ad928SBarry Smith     ierr = VecReciprocal(diag);CHKERRQ(ierr);
1734b9ad928SBarry Smith     ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr);
1744b9ad928SBarry Smith     ierr = VecGetArray(diag,&x);CHKERRQ(ierr);
175cd47f5d9SBarry Smith     if (jac->useabs) {
1762fa5cd67SKarl Rupp       for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]);
177cd47f5d9SBarry Smith     }
1784b9ad928SBarry Smith     for (i=0; i<n; i++) {
1794b9ad928SBarry Smith       if (x[i] == 0.0) {
1804b9ad928SBarry Smith         x[i]     = 1.0;
181cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
1824b9ad928SBarry Smith       }
1834b9ad928SBarry Smith     }
1844b9ad928SBarry Smith     ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr);
1854b9ad928SBarry Smith   }
1864b9ad928SBarry Smith   if (diagsqrt) {
1874b9ad928SBarry Smith     if (jac->userowmax) {
1880298fd71SBarry Smith       ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr);
18986697f06SMatthew Knepley     } else if (jac->userowsum) {
19086697f06SMatthew Knepley       ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr);
1914b9ad928SBarry Smith     } else {
1924b9ad928SBarry Smith       ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr);
1934b9ad928SBarry Smith     }
1944b9ad928SBarry Smith     ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr);
1954b9ad928SBarry Smith     ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr);
1964b9ad928SBarry Smith     for (i=0; i<n; i++) {
1978f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i]));
1984b9ad928SBarry Smith       else {
1994b9ad928SBarry Smith         x[i]     = 1.0;
200cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2014b9ad928SBarry Smith       }
2024b9ad928SBarry Smith     }
2034b9ad928SBarry Smith     ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr);
2044b9ad928SBarry Smith   }
2054b9ad928SBarry Smith   if (zeroflag) {
206ae15b995SBarry Smith     ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr);
2074b9ad928SBarry Smith   }
2084b9ad928SBarry Smith   PetscFunctionReturn(0);
2094b9ad928SBarry Smith }
2104b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2114b9ad928SBarry Smith /*
2124b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2134b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2144b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith    Input Parameter:
2174b9ad928SBarry Smith .  pc - the preconditioner context
2184b9ad928SBarry Smith */
2196849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
2204b9ad928SBarry Smith {
221dfbe8321SBarry Smith   PetscErrorCode ierr;
2224b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith   PetscFunctionBegin;
2250a545947SLisandro Dalcin   ierr = MatCreateVecs(pc->pmat,&jac->diagsqrt,NULL);CHKERRQ(ierr);
2263bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt);CHKERRQ(ierr);
2274b9ad928SBarry Smith   ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr);
2284b9ad928SBarry Smith   PetscFunctionReturn(0);
2294b9ad928SBarry Smith }
2304b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2314b9ad928SBarry Smith /*
2324b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2334b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2344b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2354b9ad928SBarry Smith 
2364b9ad928SBarry Smith    Input Parameter:
2374b9ad928SBarry Smith .  pc - the preconditioner context
2384b9ad928SBarry Smith */
2396849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
2404b9ad928SBarry Smith {
241dfbe8321SBarry Smith   PetscErrorCode ierr;
2424b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2434b9ad928SBarry Smith 
2444b9ad928SBarry Smith   PetscFunctionBegin;
2450a545947SLisandro Dalcin   ierr = MatCreateVecs(pc->pmat,&jac->diag,NULL);CHKERRQ(ierr);
2463bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);CHKERRQ(ierr);
2474b9ad928SBarry Smith   ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr);
2484b9ad928SBarry Smith   PetscFunctionReturn(0);
2494b9ad928SBarry Smith }
2504b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2514b9ad928SBarry Smith /*
2524b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2534b9ad928SBarry Smith 
2544b9ad928SBarry Smith    Input Parameters:
2554b9ad928SBarry Smith .  pc - the preconditioner context
2564b9ad928SBarry Smith .  x - input vector
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith    Output Parameter:
2594b9ad928SBarry Smith .  y - output vector
2604b9ad928SBarry Smith 
2614b9ad928SBarry Smith    Application Interface Routine: PCApply()
2624b9ad928SBarry Smith  */
2636849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y)
2644b9ad928SBarry Smith {
2654b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
266dfbe8321SBarry Smith   PetscErrorCode ierr;
2674b9ad928SBarry Smith 
2684b9ad928SBarry Smith   PetscFunctionBegin;
2694b9ad928SBarry Smith   if (!jac->diag) {
2704b9ad928SBarry Smith     ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr);
2714b9ad928SBarry Smith   }
2722dcb1b2aSMatthew Knepley   ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr);
2734b9ad928SBarry Smith   PetscFunctionReturn(0);
2744b9ad928SBarry Smith }
2754b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2764b9ad928SBarry Smith /*
2774b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
2784b9ad928SBarry Smith    symmetric preconditioner to a vector.
2794b9ad928SBarry Smith 
2804b9ad928SBarry Smith    Input Parameters:
2814b9ad928SBarry Smith .  pc - the preconditioner context
2824b9ad928SBarry Smith .  x - input vector
2834b9ad928SBarry Smith 
2844b9ad928SBarry Smith    Output Parameter:
2854b9ad928SBarry Smith .  y - output vector
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
2884b9ad928SBarry Smith */
2896849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)
2904b9ad928SBarry Smith {
291dfbe8321SBarry Smith   PetscErrorCode ierr;
2924b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith   PetscFunctionBegin;
2954b9ad928SBarry Smith   if (!jac->diagsqrt) {
2964b9ad928SBarry Smith     ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr);
2974b9ad928SBarry Smith   }
2982dcb1b2aSMatthew Knepley   VecPointwiseMult(y,x,jac->diagsqrt);
2994b9ad928SBarry Smith   PetscFunctionReturn(0);
3004b9ad928SBarry Smith }
3014b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
302a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc)
303a06653b4SBarry Smith {
304a06653b4SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
305a06653b4SBarry Smith   PetscErrorCode ierr;
306a06653b4SBarry Smith 
307a06653b4SBarry Smith   PetscFunctionBegin;
3086bf464f9SBarry Smith   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
3096bf464f9SBarry Smith   ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr);
310a06653b4SBarry Smith   PetscFunctionReturn(0);
311a06653b4SBarry Smith }
312a06653b4SBarry Smith 
3134b9ad928SBarry Smith /*
3144b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3154b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith    Input Parameter:
3184b9ad928SBarry Smith .  pc - the preconditioner context
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3214b9ad928SBarry Smith */
3226849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc)
3234b9ad928SBarry Smith {
324dfbe8321SBarry Smith   PetscErrorCode ierr;
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith   PetscFunctionBegin;
327a06653b4SBarry Smith   ierr = PCReset_Jacobi(pc);CHKERRQ(ierr);
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith   /*
3304b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3314b9ad928SBarry Smith   */
332c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
3334b9ad928SBarry Smith   PetscFunctionReturn(0);
3344b9ad928SBarry Smith }
3354b9ad928SBarry Smith 
3364416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PetscOptionItems *PetscOptionsObject,PC pc)
3374b9ad928SBarry Smith {
3384b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
339dfbe8321SBarry Smith   PetscErrorCode ierr;
340baa89ecbSBarry Smith   PetscBool      flg;
341baa89ecbSBarry Smith   PCJacobiType   deflt,type;
3424b9ad928SBarry Smith 
3434b9ad928SBarry Smith   PetscFunctionBegin;
344baa89ecbSBarry Smith   ierr = PCJacobiGetType(pc,&deflt);CHKERRQ(ierr);
345e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Jacobi options");CHKERRQ(ierr);
346baa89ecbSBarry Smith   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
347baa89ecbSBarry Smith   if (flg) {
348baa89ecbSBarry Smith     ierr = PCJacobiSetType(pc,type);CHKERRQ(ierr);
349baa89ecbSBarry Smith   }
3508afaa268SBarry Smith   ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL);CHKERRQ(ierr);
3514b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3524b9ad928SBarry Smith   PetscFunctionReturn(0);
3534b9ad928SBarry Smith }
3544b9ad928SBarry Smith 
355*569e28a7SMatthew G. Knepley static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
356*569e28a7SMatthew G. Knepley {
357*569e28a7SMatthew G. Knepley   PC_Jacobi     *jac = (PC_Jacobi *) pc->data;
358*569e28a7SMatthew G. Knepley   PetscBool      iascii;
359*569e28a7SMatthew G. Knepley   PetscErrorCode ierr;
360*569e28a7SMatthew G. Knepley 
361*569e28a7SMatthew G. Knepley   PetscFunctionBegin;
362*569e28a7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
363*569e28a7SMatthew G. Knepley   if (iascii) {
364*569e28a7SMatthew G. Knepley     PCJacobiType      type;
365*569e28a7SMatthew G. Knepley     PetscBool         useAbs;
366*569e28a7SMatthew G. Knepley     PetscViewerFormat format;
367*569e28a7SMatthew G. Knepley 
368*569e28a7SMatthew G. Knepley     ierr = PCJacobiGetType(pc, &type);CHKERRQ(ierr);
369*569e28a7SMatthew G. Knepley     ierr = PCJacobiGetUseAbs(pc, &useAbs);CHKERRQ(ierr);
370*569e28a7SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "  type %s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "");CHKERRQ(ierr);
371*569e28a7SMatthew G. Knepley     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
372*569e28a7SMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
373*569e28a7SMatthew G. Knepley       ierr = VecView(jac->diag, viewer);CHKERRQ(ierr);
374*569e28a7SMatthew G. Knepley     }
375*569e28a7SMatthew G. Knepley   }
376*569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
377*569e28a7SMatthew G. Knepley }
378*569e28a7SMatthew G. Knepley 
3794b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
3804b9ad928SBarry Smith /*
3814b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
3824b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
3834b9ad928SBarry Smith    context, PC, that was created within PCCreate().
3844b9ad928SBarry Smith 
3854b9ad928SBarry Smith    Input Parameter:
3864b9ad928SBarry Smith .  pc - the preconditioner context
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith    Application Interface Routine: PCCreate()
3894b9ad928SBarry Smith */
3904b9ad928SBarry Smith 
3914b9ad928SBarry Smith /*MC
3925a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith    Options Database Key:
395967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
396422a814eSBarry Smith -    -pc_jacobi_abs - use the absolute value of the diagonal entry
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith    Level: beginner
3994b9ad928SBarry Smith 
40095452b02SPatrick Sanan   Notes:
40195452b02SPatrick Sanan     By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric
4024b9ad928SBarry Smith          can scale each side of the matrix by the square root of the diagonal entries.
4034b9ad928SBarry Smith 
4044b9ad928SBarry Smith          Zero entries along the diagonal are replaced with the value 1.0
4054b9ad928SBarry Smith 
406422a814eSBarry Smith          See PCPBJACOBI for a point-block Jacobi preconditioner
407422a814eSBarry Smith 
4084b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
409422a814eSBarry Smith            PCJacobiSetType(), PCJacobiSetUseAbs(), PCJacobiGetUseAbs(), PCPBJACOBI
4104b9ad928SBarry Smith M*/
4114b9ad928SBarry Smith 
4128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
4134b9ad928SBarry Smith {
4144b9ad928SBarry Smith   PC_Jacobi      *jac;
415dfbe8321SBarry Smith   PetscErrorCode ierr;
4164b9ad928SBarry Smith 
4174b9ad928SBarry Smith   PetscFunctionBegin;
4184b9ad928SBarry Smith   /*
4194b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4204b9ad928SBarry Smith      attach it to the PC object.
4214b9ad928SBarry Smith   */
422b00a9115SJed Brown   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4234b9ad928SBarry Smith   pc->data = (void*)jac;
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith   /*
4264b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4274b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4284b9ad928SBarry Smith   */
4290a545947SLisandro Dalcin   jac->diag      = NULL;
4300a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4314b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
43286697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
433cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith   /*
4364b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4374b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4384b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4394b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4404b9ad928SBarry Smith       not needed.
4414b9ad928SBarry Smith   */
4424b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4434b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4444b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
445a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4464b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4474b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
448*569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4490a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4504b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4514b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4522fa5cd67SKarl Rupp 
453baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi);CHKERRQ(ierr);
454baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi);CHKERRQ(ierr);
455bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr);
456baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi);CHKERRQ(ierr);
4574b9ad928SBarry Smith   PetscFunctionReturn(0);
4584b9ad928SBarry Smith }
459cd47f5d9SBarry Smith 
460cd47f5d9SBarry Smith /*@
461cd47f5d9SBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the
462baa89ecbSBarry Smith       absolute values of the digonal divisors in the preconditioner
463cd47f5d9SBarry Smith 
464ad4df100SBarry Smith    Logically Collective on PC
465cd47f5d9SBarry Smith 
466cd47f5d9SBarry Smith    Input Parameters:
467baa89ecbSBarry Smith +  pc - the preconditioner context
468baa89ecbSBarry Smith -  flg - whether to use absolute values or not
469baa89ecbSBarry Smith 
470baa89ecbSBarry Smith    Options Database Key:
471baa89ecbSBarry Smith .  -pc_jacobi_abs
472baa89ecbSBarry Smith 
47395452b02SPatrick Sanan    Notes:
47495452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
475baa89ecbSBarry Smith 
476baa89ecbSBarry Smith    Level: intermediate
477baa89ecbSBarry Smith 
478baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiGetUseAbs()
479baa89ecbSBarry Smith 
480baa89ecbSBarry Smith @*/
481baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetUseAbs(PC pc,PetscBool flg)
482baa89ecbSBarry Smith {
483baa89ecbSBarry Smith   PetscErrorCode ierr;
484baa89ecbSBarry Smith 
485baa89ecbSBarry Smith   PetscFunctionBegin;
486baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
487baa89ecbSBarry Smith   ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
488baa89ecbSBarry Smith   PetscFunctionReturn(0);
489baa89ecbSBarry Smith }
490baa89ecbSBarry Smith 
491baa89ecbSBarry Smith /*@
492baa89ecbSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the
493baa89ecbSBarry Smith       absolute values of the digonal divisors in the preconditioner
494baa89ecbSBarry Smith 
495baa89ecbSBarry Smith    Logically Collective on PC
496baa89ecbSBarry Smith 
497baa89ecbSBarry Smith    Input Parameter:
498cd47f5d9SBarry Smith .  pc - the preconditioner context
499cd47f5d9SBarry Smith 
500baa89ecbSBarry Smith    Output Parameter:
501baa89ecbSBarry Smith .  flg - whether to use absolute values or not
502cd47f5d9SBarry Smith 
503cd47f5d9SBarry Smith    Options Database Key:
504cd47f5d9SBarry Smith .  -pc_jacobi_abs
505cd47f5d9SBarry Smith 
506cd47f5d9SBarry Smith    Level: intermediate
507cd47f5d9SBarry Smith 
508baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiSetUseAbs(), PCJacobiGetType()
509cd47f5d9SBarry Smith 
510cd47f5d9SBarry Smith @*/
511baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetUseAbs(PC pc,PetscBool *flg)
512cd47f5d9SBarry Smith {
5134ac538c5SBarry Smith   PetscErrorCode ierr;
514cd47f5d9SBarry Smith 
515cd47f5d9SBarry Smith   PetscFunctionBegin;
5160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
517baa89ecbSBarry Smith   ierr = PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
518cd47f5d9SBarry Smith   PetscFunctionReturn(0);
519cd47f5d9SBarry Smith }
520cd47f5d9SBarry Smith 
5214b9ad928SBarry Smith /*@
522baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
523baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
5244b9ad928SBarry Smith 
525ad4df100SBarry Smith    Logically Collective on PC
5264b9ad928SBarry Smith 
5274b9ad928SBarry Smith    Input Parameters:
528baa89ecbSBarry Smith +  pc - the preconditioner context
529baa89ecbSBarry Smith -  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
5304b9ad928SBarry Smith 
5314b9ad928SBarry Smith    Options Database Key:
532baa89ecbSBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum>
5334b9ad928SBarry Smith 
5344b9ad928SBarry Smith    Level: intermediate
5354b9ad928SBarry Smith 
536baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiGetType()
5374b9ad928SBarry Smith @*/
538baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetType(PC pc,PCJacobiType type)
5394b9ad928SBarry Smith {
5404ac538c5SBarry Smith   PetscErrorCode ierr;
5414b9ad928SBarry Smith 
5424b9ad928SBarry Smith   PetscFunctionBegin;
5430700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
544baa89ecbSBarry Smith   ierr = PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));CHKERRQ(ierr);
5454b9ad928SBarry Smith   PetscFunctionReturn(0);
5464b9ad928SBarry Smith }
5474b9ad928SBarry Smith 
54886697f06SMatthew Knepley /*@
549baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
55086697f06SMatthew Knepley 
551baa89ecbSBarry Smith    Not Collective on PC
55286697f06SMatthew Knepley 
553baa89ecbSBarry Smith    Input Parameter:
55486697f06SMatthew Knepley .  pc - the preconditioner context
55586697f06SMatthew Knepley 
556baa89ecbSBarry Smith    Output Parameter:
557baa89ecbSBarry Smith .  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
55886697f06SMatthew Knepley 
55986697f06SMatthew Knepley    Level: intermediate
56086697f06SMatthew Knepley 
561baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiSetType()
56286697f06SMatthew Knepley @*/
563baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetType(PC pc,PCJacobiType *type)
56486697f06SMatthew Knepley {
5654ac538c5SBarry Smith   PetscErrorCode ierr;
56686697f06SMatthew Knepley 
56786697f06SMatthew Knepley   PetscFunctionBegin;
5680700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
569baa89ecbSBarry Smith   ierr = PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));CHKERRQ(ierr);
57086697f06SMatthew Knepley   PetscFunctionReturn(0);
57186697f06SMatthew Knepley }
572