xref: /petsc/src/ksp/pc/impls/jacobi/jacobi.c (revision 67ed0f3b56ac23720fdaae488c8b128091e5c6cc)
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 */
66*67ed0f3bSStefano 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 
117*67ed0f3bSStefano Zampini static PetscErrorCode  PCJacobiSetFixDiagonal_Jacobi(PC pc,PetscBool flg)
118*67ed0f3bSStefano Zampini {
119*67ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi*)pc->data;
120*67ed0f3bSStefano Zampini 
121*67ed0f3bSStefano Zampini   PetscFunctionBegin;
122*67ed0f3bSStefano Zampini   j->fixdiag = flg;
123*67ed0f3bSStefano Zampini   PetscFunctionReturn(0);
124*67ed0f3bSStefano Zampini }
125*67ed0f3bSStefano Zampini 
126*67ed0f3bSStefano Zampini static PetscErrorCode  PCJacobiGetFixDiagonal_Jacobi(PC pc,PetscBool *flg)
127*67ed0f3bSStefano Zampini {
128*67ed0f3bSStefano Zampini   PC_Jacobi *j = (PC_Jacobi*)pc->data;
129*67ed0f3bSStefano Zampini 
130*67ed0f3bSStefano Zampini   PetscFunctionBegin;
131*67ed0f3bSStefano Zampini   *flg = j->fixdiag;
132*67ed0f3bSStefano Zampini   PetscFunctionReturn(0);
133*67ed0f3bSStefano Zampini }
134*67ed0f3bSStefano 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;
153dfbe8321SBarry Smith   PetscErrorCode ierr;
154cd47f5d9SBarry Smith   PetscInt       n,i;
1554b9ad928SBarry Smith   PetscScalar    *x;
156ace3abfcSBarry Smith   PetscBool      zeroflag = PETSC_FALSE;
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   PetscFunctionBegin;
1594b9ad928SBarry Smith   /*
1604b9ad928SBarry Smith        For most preconditioners the code would begin here something like
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith   if (pc->setupcalled == 0) { allocate space the first time this is ever called
1632a7a6963SBarry Smith     ierr = MatCreateVecs(pc->mat,&jac->diag);CHKERRQ(ierr);
1643bb1ff40SBarry Smith     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
1654b9ad928SBarry Smith   }
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith     But for this preconditioner we want to support use of both the matrix' diagonal
1684b9ad928SBarry Smith     elements (for left or right preconditioning) and square root of diagonal elements
1694b9ad928SBarry Smith     (for symmetric preconditioning).  Hence we do not allocate space here, since we
1704b9ad928SBarry Smith     don't know at this point which will be needed (diag and/or diagsqrt) until the user
1714b9ad928SBarry Smith     applies the preconditioner, and we don't want to allocate BOTH unless we need
1724b9ad928SBarry Smith     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1734b9ad928SBarry Smith     and PCSetUp_Jacobi_Symmetric(), respectively.
1744b9ad928SBarry Smith   */
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith   /*
1774b9ad928SBarry Smith     Here we set up the preconditioner; that is, we copy the diagonal values from
1784b9ad928SBarry Smith     the matrix and put them into a format to make them quick to apply as a preconditioner.
1794b9ad928SBarry Smith   */
1804b9ad928SBarry Smith   diag     = jac->diag;
1814b9ad928SBarry Smith   diagsqrt = jac->diagsqrt;
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith   if (diag) {
184*67ed0f3bSStefano Zampini     PetscBool isspd;
185*67ed0f3bSStefano Zampini 
1864b9ad928SBarry Smith     if (jac->userowmax) {
1870298fd71SBarry Smith       ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr);
18886697f06SMatthew Knepley     } else if (jac->userowsum) {
18986697f06SMatthew Knepley       ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr);
1904b9ad928SBarry Smith     } else {
1914b9ad928SBarry Smith       ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr);
1924b9ad928SBarry Smith     }
1934b9ad928SBarry Smith     ierr = VecReciprocal(diag);CHKERRQ(ierr);
194cd47f5d9SBarry Smith     if (jac->useabs) {
195c2208464SStefano Zampini       ierr = VecAbs(diag);CHKERRQ(ierr);
196cd47f5d9SBarry Smith     }
197*67ed0f3bSStefano Zampini     ierr = MatGetOption(pc->pmat,MAT_SPD,&isspd);CHKERRQ(ierr);
198*67ed0f3bSStefano Zampini     if (jac->fixdiag && !isspd) {
199*67ed0f3bSStefano Zampini       ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr);
200c2208464SStefano Zampini       ierr = VecGetArray(diag,&x);CHKERRQ(ierr);
2014b9ad928SBarry Smith       for (i=0; i<n; i++) {
2024b9ad928SBarry Smith         if (x[i] == 0.0) {
2034b9ad928SBarry Smith           x[i]     = 1.0;
204cd47f5d9SBarry Smith           zeroflag = PETSC_TRUE;
2054b9ad928SBarry Smith         }
2064b9ad928SBarry Smith       }
2074b9ad928SBarry Smith       ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr);
2084b9ad928SBarry Smith     }
209*67ed0f3bSStefano Zampini   }
2104b9ad928SBarry Smith   if (diagsqrt) {
2114b9ad928SBarry Smith     if (jac->userowmax) {
2120298fd71SBarry Smith       ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr);
21386697f06SMatthew Knepley     } else if (jac->userowsum) {
21486697f06SMatthew Knepley       ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr);
2154b9ad928SBarry Smith     } else {
2164b9ad928SBarry Smith       ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr);
2174b9ad928SBarry Smith     }
2184b9ad928SBarry Smith     ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr);
2194b9ad928SBarry Smith     ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr);
2204b9ad928SBarry Smith     for (i=0; i<n; i++) {
2218f1a2a5eSBarry Smith       if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i]));
2224b9ad928SBarry Smith       else {
2234b9ad928SBarry Smith         x[i]     = 1.0;
224cd47f5d9SBarry Smith         zeroflag = PETSC_TRUE;
2254b9ad928SBarry Smith       }
2264b9ad928SBarry Smith     }
2274b9ad928SBarry Smith     ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr);
2284b9ad928SBarry Smith   }
2294b9ad928SBarry Smith   if (zeroflag) {
230ae15b995SBarry Smith     ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr);
2314b9ad928SBarry Smith   }
2324b9ad928SBarry Smith   PetscFunctionReturn(0);
2334b9ad928SBarry Smith }
2344b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2354b9ad928SBarry Smith /*
2364b9ad928SBarry Smith    PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
2374b9ad928SBarry Smith    inverse of the square root of the diagonal entries of the matrix.  This
2384b9ad928SBarry Smith    is used for symmetric application of the Jacobi preconditioner.
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith    Input Parameter:
2414b9ad928SBarry Smith .  pc - the preconditioner context
2424b9ad928SBarry Smith */
2436849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
2444b9ad928SBarry Smith {
245dfbe8321SBarry Smith   PetscErrorCode ierr;
2464b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2474b9ad928SBarry Smith 
2484b9ad928SBarry Smith   PetscFunctionBegin;
2490a545947SLisandro Dalcin   ierr = MatCreateVecs(pc->pmat,&jac->diagsqrt,NULL);CHKERRQ(ierr);
2503bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt);CHKERRQ(ierr);
2514b9ad928SBarry Smith   ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr);
2524b9ad928SBarry Smith   PetscFunctionReturn(0);
2534b9ad928SBarry Smith }
2544b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2554b9ad928SBarry Smith /*
2564b9ad928SBarry Smith    PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
2574b9ad928SBarry Smith    inverse of the diagonal entries of the matrix.  This is used for left of
2584b9ad928SBarry Smith    right application of the Jacobi preconditioner.
2594b9ad928SBarry Smith 
2604b9ad928SBarry Smith    Input Parameter:
2614b9ad928SBarry Smith .  pc - the preconditioner context
2624b9ad928SBarry Smith */
2636849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
2644b9ad928SBarry Smith {
265dfbe8321SBarry Smith   PetscErrorCode ierr;
2664b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
2674b9ad928SBarry Smith 
2684b9ad928SBarry Smith   PetscFunctionBegin;
2690a545947SLisandro Dalcin   ierr = MatCreateVecs(pc->pmat,&jac->diag,NULL);CHKERRQ(ierr);
2703bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);CHKERRQ(ierr);
2714b9ad928SBarry Smith   ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr);
2724b9ad928SBarry Smith   PetscFunctionReturn(0);
2734b9ad928SBarry Smith }
2744b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2754b9ad928SBarry Smith /*
2764b9ad928SBarry Smith    PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith    Input Parameters:
2794b9ad928SBarry Smith .  pc - the preconditioner context
2804b9ad928SBarry Smith .  x - input vector
2814b9ad928SBarry Smith 
2824b9ad928SBarry Smith    Output Parameter:
2834b9ad928SBarry Smith .  y - output vector
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith    Application Interface Routine: PCApply()
2864b9ad928SBarry Smith  */
2876849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y)
2884b9ad928SBarry Smith {
2894b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
290dfbe8321SBarry Smith   PetscErrorCode ierr;
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith   PetscFunctionBegin;
2934b9ad928SBarry Smith   if (!jac->diag) {
2944b9ad928SBarry Smith     ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr);
2954b9ad928SBarry Smith   }
2962dcb1b2aSMatthew Knepley   ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr);
2974b9ad928SBarry Smith   PetscFunctionReturn(0);
2984b9ad928SBarry Smith }
2994b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
3004b9ad928SBarry Smith /*
3014b9ad928SBarry Smith    PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
3024b9ad928SBarry Smith    symmetric preconditioner to a vector.
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith    Input Parameters:
3054b9ad928SBarry Smith .  pc - the preconditioner context
3064b9ad928SBarry Smith .  x - input vector
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith    Output Parameter:
3094b9ad928SBarry Smith .  y - output vector
3104b9ad928SBarry Smith 
3114b9ad928SBarry Smith    Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
3124b9ad928SBarry Smith */
3136849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)
3144b9ad928SBarry Smith {
315dfbe8321SBarry Smith   PetscErrorCode ierr;
3164b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
3174b9ad928SBarry Smith 
3184b9ad928SBarry Smith   PetscFunctionBegin;
3194b9ad928SBarry Smith   if (!jac->diagsqrt) {
3204b9ad928SBarry Smith     ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr);
3214b9ad928SBarry Smith   }
322*67ed0f3bSStefano Zampini   ierr = VecPointwiseMult(y,x,jac->diagsqrt);CHKERRQ(ierr);
3234b9ad928SBarry Smith   PetscFunctionReturn(0);
3244b9ad928SBarry Smith }
325*67ed0f3bSStefano Zampini 
3264b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
327a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc)
328a06653b4SBarry Smith {
329a06653b4SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
330a06653b4SBarry Smith   PetscErrorCode ierr;
331a06653b4SBarry Smith 
332a06653b4SBarry Smith   PetscFunctionBegin;
3336bf464f9SBarry Smith   ierr = VecDestroy(&jac->diag);CHKERRQ(ierr);
3346bf464f9SBarry Smith   ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr);
335a06653b4SBarry Smith   PetscFunctionReturn(0);
336a06653b4SBarry Smith }
337a06653b4SBarry Smith 
3384b9ad928SBarry Smith /*
3394b9ad928SBarry Smith    PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
3404b9ad928SBarry Smith    that was created with PCCreate_Jacobi().
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith    Input Parameter:
3434b9ad928SBarry Smith .  pc - the preconditioner context
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith    Application Interface Routine: PCDestroy()
3464b9ad928SBarry Smith */
3476849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc)
3484b9ad928SBarry Smith {
349dfbe8321SBarry Smith   PetscErrorCode ierr;
3504b9ad928SBarry Smith 
3514b9ad928SBarry Smith   PetscFunctionBegin;
352a06653b4SBarry Smith   ierr = PCReset_Jacobi(pc);CHKERRQ(ierr);
353*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",NULL);CHKERRQ(ierr);
354*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",NULL);CHKERRQ(ierr);
355*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",NULL);CHKERRQ(ierr);
356*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",NULL);CHKERRQ(ierr);
357*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",NULL);CHKERRQ(ierr);
358*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",NULL);CHKERRQ(ierr);
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith   /*
3614b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
3624b9ad928SBarry Smith   */
363c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
3644b9ad928SBarry Smith   PetscFunctionReturn(0);
3654b9ad928SBarry Smith }
3664b9ad928SBarry Smith 
3674416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PetscOptionItems *PetscOptionsObject,PC pc)
3684b9ad928SBarry Smith {
3694b9ad928SBarry Smith   PC_Jacobi      *jac = (PC_Jacobi*)pc->data;
370dfbe8321SBarry Smith   PetscErrorCode ierr;
371baa89ecbSBarry Smith   PetscBool      flg;
372baa89ecbSBarry Smith   PCJacobiType   deflt,type;
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith   PetscFunctionBegin;
375baa89ecbSBarry Smith   ierr = PCJacobiGetType(pc,&deflt);CHKERRQ(ierr);
376e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Jacobi options");CHKERRQ(ierr);
377baa89ecbSBarry Smith   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
378baa89ecbSBarry Smith   if (flg) {
379baa89ecbSBarry Smith     ierr = PCJacobiSetType(pc,type);CHKERRQ(ierr);
380baa89ecbSBarry Smith   }
381*67ed0f3bSStefano Zampini   ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL);CHKERRQ(ierr);
382*67ed0f3bSStefano Zampini   ierr = PetscOptionsBool("-pc_jacobi_fixdiagonal","Fix null terms on diagonal","PCJacobiSetFixDiagonal",jac->fixdiag,&jac->fixdiag,NULL);CHKERRQ(ierr);
3834b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3844b9ad928SBarry Smith   PetscFunctionReturn(0);
3854b9ad928SBarry Smith }
3864b9ad928SBarry Smith 
387569e28a7SMatthew G. Knepley static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
388569e28a7SMatthew G. Knepley {
389569e28a7SMatthew G. Knepley   PC_Jacobi     *jac = (PC_Jacobi *) pc->data;
390569e28a7SMatthew G. Knepley   PetscBool      iascii;
391569e28a7SMatthew G. Knepley   PetscErrorCode ierr;
392569e28a7SMatthew G. Knepley 
393569e28a7SMatthew G. Knepley   PetscFunctionBegin;
394569e28a7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
395569e28a7SMatthew G. Knepley   if (iascii) {
396569e28a7SMatthew G. Knepley     PCJacobiType      type;
397*67ed0f3bSStefano Zampini     PetscBool         useAbs,fixdiag;
398569e28a7SMatthew G. Knepley     PetscViewerFormat format;
399569e28a7SMatthew G. Knepley 
400569e28a7SMatthew G. Knepley     ierr = PCJacobiGetType(pc, &type);CHKERRQ(ierr);
401569e28a7SMatthew G. Knepley     ierr = PCJacobiGetUseAbs(pc, &useAbs);CHKERRQ(ierr);
402*67ed0f3bSStefano Zampini     ierr = PCJacobiGetFixDiagonal(pc, &fixdiag);CHKERRQ(ierr);
403*67ed0f3bSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer, "  type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : "");CHKERRQ(ierr);
404569e28a7SMatthew G. Knepley     ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
405569e28a7SMatthew G. Knepley     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
406569e28a7SMatthew G. Knepley       ierr = VecView(jac->diag, viewer);CHKERRQ(ierr);
407569e28a7SMatthew G. Knepley     }
408569e28a7SMatthew G. Knepley   }
409569e28a7SMatthew G. Knepley   PetscFunctionReturn(0);
410569e28a7SMatthew G. Knepley }
411569e28a7SMatthew G. Knepley 
4124b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
4134b9ad928SBarry Smith /*
4144b9ad928SBarry Smith    PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
4154b9ad928SBarry Smith    and sets this as the private data within the generic preconditioning
4164b9ad928SBarry Smith    context, PC, that was created within PCCreate().
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith    Input Parameter:
4194b9ad928SBarry Smith .  pc - the preconditioner context
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith    Application Interface Routine: PCCreate()
4224b9ad928SBarry Smith */
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith /*MC
4255a46d3faSBarry Smith      PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Options Database Key:
428967c93d3SBarry Smith +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
429*67ed0f3bSStefano Zampini .    -pc_jacobi_abs - use the absolute value of the diagonal entry
430*67ed0f3bSStefano Zampini -    -pc_jacobi_fixdiag - fix for zero diagonal terms
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith    Level: beginner
4334b9ad928SBarry Smith 
43495452b02SPatrick Sanan   Notes:
43595452b02SPatrick Sanan     By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric
4364b9ad928SBarry Smith          can scale each side of the matrix by the square root of the diagonal entries.
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith          Zero entries along the diagonal are replaced with the value 1.0
4394b9ad928SBarry Smith 
440422a814eSBarry Smith          See PCPBJACOBI for a point-block Jacobi preconditioner
441422a814eSBarry Smith 
4424b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
443*67ed0f3bSStefano Zampini            PCJacobiSetType(), PCJacobiSetUseAbs(), PCJacobiGetUseAbs(),
444*67ed0f3bSStefano Zampini            PCJacobiSetFixDiagonal(), PCJacobiGetFixDiagonal(), PCPBJACOBI
4454b9ad928SBarry Smith M*/
4464b9ad928SBarry Smith 
4478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
4484b9ad928SBarry Smith {
4494b9ad928SBarry Smith   PC_Jacobi      *jac;
450dfbe8321SBarry Smith   PetscErrorCode ierr;
4514b9ad928SBarry Smith 
4524b9ad928SBarry Smith   PetscFunctionBegin;
4534b9ad928SBarry Smith   /*
4544b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
4554b9ad928SBarry Smith      attach it to the PC object.
4564b9ad928SBarry Smith   */
457b00a9115SJed Brown   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4584b9ad928SBarry Smith   pc->data = (void*)jac;
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith   /*
4614b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
4624b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
4634b9ad928SBarry Smith   */
4640a545947SLisandro Dalcin   jac->diag      = NULL;
4650a545947SLisandro Dalcin   jac->diagsqrt  = NULL;
4664b9ad928SBarry Smith   jac->userowmax = PETSC_FALSE;
46786697f06SMatthew Knepley   jac->userowsum = PETSC_FALSE;
468cd47f5d9SBarry Smith   jac->useabs    = PETSC_FALSE;
469*67ed0f3bSStefano Zampini   jac->fixdiag   = PETSC_TRUE;
4704b9ad928SBarry Smith 
4714b9ad928SBarry Smith   /*
4724b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
4734b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
4744b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
4754b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
4764b9ad928SBarry Smith       not needed.
4774b9ad928SBarry Smith   */
4784b9ad928SBarry Smith   pc->ops->apply               = PCApply_Jacobi;
4794b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_Jacobi;
4804b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_Jacobi;
481a06653b4SBarry Smith   pc->ops->reset               = PCReset_Jacobi;
4824b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Jacobi;
4834b9ad928SBarry Smith   pc->ops->setfromoptions      = PCSetFromOptions_Jacobi;
484569e28a7SMatthew G. Knepley   pc->ops->view                = PCView_Jacobi;
4850a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
4864b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeftOrRight_Jacobi;
4874b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
4882fa5cd67SKarl Rupp 
489baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi);CHKERRQ(ierr);
490baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi);CHKERRQ(ierr);
491bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr);
492baa89ecbSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi);CHKERRQ(ierr);
493*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetFixDiagonal_C",PCJacobiSetFixDiagonal_Jacobi);CHKERRQ(ierr);
494*67ed0f3bSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetFixDiagonal_C",PCJacobiGetFixDiagonal_Jacobi);CHKERRQ(ierr);
4954b9ad928SBarry Smith   PetscFunctionReturn(0);
4964b9ad928SBarry Smith }
497cd47f5d9SBarry Smith 
498cd47f5d9SBarry Smith /*@
499cd47f5d9SBarry Smith    PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the
500*67ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
501cd47f5d9SBarry Smith 
502ad4df100SBarry Smith    Logically Collective on PC
503cd47f5d9SBarry Smith 
504cd47f5d9SBarry Smith    Input Parameters:
505baa89ecbSBarry Smith +  pc - the preconditioner context
506baa89ecbSBarry Smith -  flg - whether to use absolute values or not
507baa89ecbSBarry Smith 
508baa89ecbSBarry Smith    Options Database Key:
509baa89ecbSBarry Smith .  -pc_jacobi_abs
510baa89ecbSBarry Smith 
51195452b02SPatrick Sanan    Notes:
51295452b02SPatrick Sanan     This takes affect at the next construction of the preconditioner
513baa89ecbSBarry Smith 
514baa89ecbSBarry Smith    Level: intermediate
515baa89ecbSBarry Smith 
516baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiGetUseAbs()
517baa89ecbSBarry Smith 
518baa89ecbSBarry Smith @*/
519baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetUseAbs(PC pc,PetscBool flg)
520baa89ecbSBarry Smith {
521baa89ecbSBarry Smith   PetscErrorCode ierr;
522baa89ecbSBarry Smith 
523baa89ecbSBarry Smith   PetscFunctionBegin;
524baa89ecbSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
525baa89ecbSBarry Smith   ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
526baa89ecbSBarry Smith   PetscFunctionReturn(0);
527baa89ecbSBarry Smith }
528baa89ecbSBarry Smith 
529baa89ecbSBarry Smith /*@
530baa89ecbSBarry Smith    PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the
531*67ed0f3bSStefano Zampini       absolute values of the diagonal divisors in the preconditioner
532baa89ecbSBarry Smith 
533baa89ecbSBarry Smith    Logically Collective on PC
534baa89ecbSBarry Smith 
535baa89ecbSBarry Smith    Input Parameter:
536cd47f5d9SBarry Smith .  pc - the preconditioner context
537cd47f5d9SBarry Smith 
538baa89ecbSBarry Smith    Output Parameter:
539baa89ecbSBarry Smith .  flg - whether to use absolute values or not
540cd47f5d9SBarry Smith 
541cd47f5d9SBarry Smith    Options Database Key:
542cd47f5d9SBarry Smith .  -pc_jacobi_abs
543cd47f5d9SBarry Smith 
544cd47f5d9SBarry Smith    Level: intermediate
545cd47f5d9SBarry Smith 
546baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiSetUseAbs(), PCJacobiGetType()
547cd47f5d9SBarry Smith 
548cd47f5d9SBarry Smith @*/
549baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetUseAbs(PC pc,PetscBool *flg)
550cd47f5d9SBarry Smith {
5514ac538c5SBarry Smith   PetscErrorCode ierr;
552cd47f5d9SBarry Smith 
553cd47f5d9SBarry Smith   PetscFunctionBegin;
5540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
555baa89ecbSBarry Smith   ierr = PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
556cd47f5d9SBarry Smith   PetscFunctionReturn(0);
557cd47f5d9SBarry Smith }
558cd47f5d9SBarry Smith 
5594b9ad928SBarry Smith /*@
560*67ed0f3bSStefano Zampini    PCJacobiSetFixDiagonal - Do not check for zero values on diagonal
561*67ed0f3bSStefano Zampini 
562*67ed0f3bSStefano Zampini    Logically Collective on PC
563*67ed0f3bSStefano Zampini 
564*67ed0f3bSStefano Zampini    Input Parameters:
565*67ed0f3bSStefano Zampini +  pc - the preconditioner context
566*67ed0f3bSStefano Zampini -  flg - the boolean flag
567*67ed0f3bSStefano Zampini 
568*67ed0f3bSStefano Zampini    Options Database Key:
569*67ed0f3bSStefano Zampini .  -pc_jacobi_fixdiagonal
570*67ed0f3bSStefano Zampini 
571*67ed0f3bSStefano Zampini    Notes:
572*67ed0f3bSStefano Zampini     This takes affect at the next construction of the preconditioner
573*67ed0f3bSStefano Zampini 
574*67ed0f3bSStefano Zampini    Level: intermediate
575*67ed0f3bSStefano Zampini 
576*67ed0f3bSStefano Zampini .seealso: PCJacobiSetType(), PCJacobiGetFixDiagonal()
577*67ed0f3bSStefano Zampini 
578*67ed0f3bSStefano Zampini @*/
579*67ed0f3bSStefano Zampini PetscErrorCode  PCJacobiSetFixDiagonal(PC pc,PetscBool flg)
580*67ed0f3bSStefano Zampini {
581*67ed0f3bSStefano Zampini   PetscErrorCode ierr;
582*67ed0f3bSStefano Zampini 
583*67ed0f3bSStefano Zampini   PetscFunctionBegin;
584*67ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
585*67ed0f3bSStefano Zampini   ierr = PetscTryMethod(pc,"PCJacobiSetFixDiagonal_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
586*67ed0f3bSStefano Zampini   PetscFunctionReturn(0);
587*67ed0f3bSStefano Zampini }
588*67ed0f3bSStefano Zampini 
589*67ed0f3bSStefano Zampini /*@
590*67ed0f3bSStefano Zampini    PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner checks for zero diagonal terms
591*67ed0f3bSStefano Zampini 
592*67ed0f3bSStefano Zampini    Logically Collective on PC
593*67ed0f3bSStefano Zampini 
594*67ed0f3bSStefano Zampini    Input Parameter:
595*67ed0f3bSStefano Zampini .  pc - the preconditioner context
596*67ed0f3bSStefano Zampini 
597*67ed0f3bSStefano Zampini    Output Parameter:
598*67ed0f3bSStefano Zampini .  flg - the boolean flag
599*67ed0f3bSStefano Zampini 
600*67ed0f3bSStefano Zampini    Options Database Key:
601*67ed0f3bSStefano Zampini .  -pc_jacobi_fixdiagonal
602*67ed0f3bSStefano Zampini 
603*67ed0f3bSStefano Zampini    Level: intermediate
604*67ed0f3bSStefano Zampini 
605*67ed0f3bSStefano Zampini .seealso: PCJacobiSetType(), PCJacobiSetFixDiagonal()
606*67ed0f3bSStefano Zampini 
607*67ed0f3bSStefano Zampini @*/
608*67ed0f3bSStefano Zampini PetscErrorCode  PCJacobiGetFixDiagonal(PC pc,PetscBool *flg)
609*67ed0f3bSStefano Zampini {
610*67ed0f3bSStefano Zampini   PetscErrorCode ierr;
611*67ed0f3bSStefano Zampini 
612*67ed0f3bSStefano Zampini   PetscFunctionBegin;
613*67ed0f3bSStefano Zampini   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
614*67ed0f3bSStefano Zampini   ierr = PetscUseMethod(pc,"PCJacobiGetFixDiagonal_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
615*67ed0f3bSStefano Zampini   PetscFunctionReturn(0);
616*67ed0f3bSStefano Zampini }
617*67ed0f3bSStefano Zampini 
618*67ed0f3bSStefano Zampini /*@
619baa89ecbSBarry Smith    PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
620baa89ecbSBarry Smith       of the sum of rows entries for the diagonal preconditioner
6214b9ad928SBarry Smith 
622ad4df100SBarry Smith    Logically Collective on PC
6234b9ad928SBarry Smith 
6244b9ad928SBarry Smith    Input Parameters:
625baa89ecbSBarry Smith +  pc - the preconditioner context
626baa89ecbSBarry Smith -  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
6274b9ad928SBarry Smith 
6284b9ad928SBarry Smith    Options Database Key:
629baa89ecbSBarry Smith .  -pc_jacobi_type <diagonal,rowmax,rowsum>
6304b9ad928SBarry Smith 
6314b9ad928SBarry Smith    Level: intermediate
6324b9ad928SBarry Smith 
633baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiGetType()
6344b9ad928SBarry Smith @*/
635baa89ecbSBarry Smith PetscErrorCode  PCJacobiSetType(PC pc,PCJacobiType type)
6364b9ad928SBarry Smith {
6374ac538c5SBarry Smith   PetscErrorCode ierr;
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith   PetscFunctionBegin;
6400700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
641baa89ecbSBarry Smith   ierr = PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));CHKERRQ(ierr);
6424b9ad928SBarry Smith   PetscFunctionReturn(0);
6434b9ad928SBarry Smith }
6444b9ad928SBarry Smith 
64586697f06SMatthew Knepley /*@
646baa89ecbSBarry Smith    PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
64786697f06SMatthew Knepley 
648baa89ecbSBarry Smith    Not Collective on PC
64986697f06SMatthew Knepley 
650baa89ecbSBarry Smith    Input Parameter:
65186697f06SMatthew Knepley .  pc - the preconditioner context
65286697f06SMatthew Knepley 
653baa89ecbSBarry Smith    Output Parameter:
654baa89ecbSBarry Smith .  type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM
65586697f06SMatthew Knepley 
65686697f06SMatthew Knepley    Level: intermediate
65786697f06SMatthew Knepley 
658baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiSetType()
65986697f06SMatthew Knepley @*/
660baa89ecbSBarry Smith PetscErrorCode  PCJacobiGetType(PC pc,PCJacobiType *type)
66186697f06SMatthew Knepley {
6624ac538c5SBarry Smith   PetscErrorCode ierr;
66386697f06SMatthew Knepley 
66486697f06SMatthew Knepley   PetscFunctionBegin;
6650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
666baa89ecbSBarry Smith   ierr = PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));CHKERRQ(ierr);
66786697f06SMatthew Knepley   PetscFunctionReturn(0);
66886697f06SMatthew Knepley }
669