104c3f3b8SBarry Smith /*
25a46d3faSBarry Smith This file implements a Jacobi preconditioner in PETSc as part of PC.
35a46d3faSBarry Smith You can use this as a starting point for implementing your own
45a46d3faSBarry Smith preconditioner that is not provided with PETSc. (You might also consider
55a46d3faSBarry Smith just using PCSHELL)
64b9ad928SBarry Smith
74b9ad928SBarry Smith The following basic routines are required for each preconditioner.
84b9ad928SBarry Smith PCCreate_XXX() - Creates a preconditioner context
94b9ad928SBarry Smith PCSetFromOptions_XXX() - Sets runtime options
104b9ad928SBarry Smith PCApply_XXX() - Applies the preconditioner
114b9ad928SBarry Smith PCDestroy_XXX() - Destroys the preconditioner context
124b9ad928SBarry Smith where the suffix "_XXX" denotes a particular implementation, in
134b9ad928SBarry Smith this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi).
144b9ad928SBarry Smith These routines are actually called via the common user interface
154b9ad928SBarry Smith routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
164b9ad928SBarry Smith so the application code interface remains identical for all
174b9ad928SBarry Smith preconditioners.
184b9ad928SBarry Smith
194b9ad928SBarry Smith Another key routine is:
204b9ad928SBarry Smith PCSetUp_XXX() - Prepares for the use of a preconditioner
214b9ad928SBarry Smith by setting data structures and options. The interface routine PCSetUp()
224b9ad928SBarry Smith is not usually called directly by the user, but instead is called by
234b9ad928SBarry Smith PCApply() if necessary.
244b9ad928SBarry Smith
254b9ad928SBarry Smith Additional basic routines are:
264b9ad928SBarry Smith PCView_XXX() - Prints details of runtime options that
274b9ad928SBarry Smith have actually been used.
284b9ad928SBarry Smith These are called by application codes via the interface routines
294b9ad928SBarry Smith PCView().
304b9ad928SBarry Smith
314b9ad928SBarry Smith The various types of solvers (preconditioners, Krylov subspace methods,
324b9ad928SBarry Smith nonlinear solvers, timesteppers) are all organized similarly, so the
334b9ad928SBarry Smith above description applies to these categories also. One exception is
344b9ad928SBarry Smith that the analogues of PCApply() for these components are KSPSolve(),
354b9ad928SBarry Smith SNESSolve(), and TSSolve().
364b9ad928SBarry Smith
374b9ad928SBarry Smith Additional optional functionality unique to preconditioners is left and
384b9ad928SBarry Smith right symmetric preconditioner application via PCApplySymmetricLeft()
394b9ad928SBarry Smith and PCApplySymmetricRight(). The Jacobi implementation is
404b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi().
4104c3f3b8SBarry Smith */
424b9ad928SBarry Smith
434b9ad928SBarry Smith /*
444b9ad928SBarry Smith Include files needed for the Jacobi preconditioner:
454b9ad928SBarry Smith pcimpl.h - private include file intended for use by all preconditioners
464b9ad928SBarry Smith */
474b9ad928SBarry Smith
48af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
494b9ad928SBarry Smith
50eede4a3fSMark Adams const char *const PCJacobiTypes[] = {"DIAGONAL", "ROWL1", "ROWMAX", "ROWSUM", "PCJacobiType", "PC_JACOBI_", NULL};
51baa89ecbSBarry Smith
524b9ad928SBarry Smith /*
534b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner.
544b9ad928SBarry Smith */
554b9ad928SBarry Smith typedef struct {
567addb90fSBarry Smith Vec diag; /* vector containing the reciprocals of the diagonal elements of the matrix used to construct the preconditioner */
574b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of
587addb90fSBarry Smith the diagonal elements of the matrix used to compute the preconditioner (used
594b9ad928SBarry Smith only for symmetric preconditioner application) */
60eede4a3fSMark Adams PCJacobiType type;
61ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */
6267ed0f3bSStefano Zampini PetscBool fixdiag; /* fix zero diagonal terms */
63eede4a3fSMark Adams PetscReal scale; /* for scaling rowl1 off-diagonals */
644b9ad928SBarry Smith } PC_Jacobi;
654b9ad928SBarry Smith
667a660b1eSJed Brown static PetscErrorCode PCReset_Jacobi(PC);
677a660b1eSJed Brown
PCJacobiSetType_Jacobi(PC pc,PCJacobiType type)68d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetType_Jacobi(PC pc, PCJacobiType type)
69d71ae5a4SJacob Faibussowitsch {
70baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data;
717a660b1eSJed Brown PCJacobiType old_type;
724b9ad928SBarry Smith
734b9ad928SBarry Smith PetscFunctionBegin;
747a660b1eSJed Brown PetscCall(PCJacobiGetType(pc, &old_type));
753ba16761SJacob Faibussowitsch if (old_type == type) PetscFunctionReturn(PETSC_SUCCESS);
767a660b1eSJed Brown PetscCall(PCReset_Jacobi(pc));
77eede4a3fSMark Adams j->type = type;
783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
794b9ad928SBarry Smith }
804b9ad928SBarry Smith
PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool * flg)81eede4a3fSMark Adams static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc, PetscBool *flg)
82d71ae5a4SJacob Faibussowitsch {
83baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data;
8486697f06SMatthew Knepley
8586697f06SMatthew Knepley PetscFunctionBegin;
86eede4a3fSMark Adams *flg = j->useabs;
873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8886697f06SMatthew Knepley }
8986697f06SMatthew Knepley
PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg)90d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc, PetscBool flg)
91d71ae5a4SJacob Faibussowitsch {
92baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data;
93cd47f5d9SBarry Smith
94cd47f5d9SBarry Smith PetscFunctionBegin;
95baa89ecbSBarry Smith j->useabs = flg;
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97baa89ecbSBarry Smith }
98baa89ecbSBarry Smith
PCJacobiGetType_Jacobi(PC pc,PCJacobiType * type)99eede4a3fSMark Adams static PetscErrorCode PCJacobiGetType_Jacobi(PC pc, PCJacobiType *type)
100d71ae5a4SJacob Faibussowitsch {
101baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data;
102baa89ecbSBarry Smith
103baa89ecbSBarry Smith PetscFunctionBegin;
104eede4a3fSMark Adams *type = j->type;
105eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
106eede4a3fSMark Adams }
107eede4a3fSMark Adams
PCJacobiSetRowl1Scale_Jacobi(PC pc,PetscReal flg)108eede4a3fSMark Adams static PetscErrorCode PCJacobiSetRowl1Scale_Jacobi(PC pc, PetscReal flg)
109eede4a3fSMark Adams {
110eede4a3fSMark Adams PC_Jacobi *j = (PC_Jacobi *)pc->data;
111eede4a3fSMark Adams
112eede4a3fSMark Adams PetscFunctionBegin;
113eede4a3fSMark Adams j->scale = flg;
114eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
115eede4a3fSMark Adams }
116eede4a3fSMark Adams
PCJacobiGetRowl1Scale_Jacobi(PC pc,PetscReal * flg)117eede4a3fSMark Adams static PetscErrorCode PCJacobiGetRowl1Scale_Jacobi(PC pc, PetscReal *flg)
118eede4a3fSMark Adams {
119eede4a3fSMark Adams PC_Jacobi *j = (PC_Jacobi *)pc->data;
120eede4a3fSMark Adams
121eede4a3fSMark Adams PetscFunctionBegin;
122eede4a3fSMark Adams *flg = j->scale;
1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
124cd47f5d9SBarry Smith }
125cd47f5d9SBarry Smith
PCJacobiSetFixDiagonal_Jacobi(PC pc,PetscBool flg)126d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiSetFixDiagonal_Jacobi(PC pc, PetscBool flg)
127d71ae5a4SJacob Faibussowitsch {
12867ed0f3bSStefano Zampini PC_Jacobi *j = (PC_Jacobi *)pc->data;
12967ed0f3bSStefano Zampini
13067ed0f3bSStefano Zampini PetscFunctionBegin;
13167ed0f3bSStefano Zampini j->fixdiag = flg;
1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13367ed0f3bSStefano Zampini }
13467ed0f3bSStefano Zampini
PCJacobiGetFixDiagonal_Jacobi(PC pc,PetscBool * flg)135d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCJacobiGetFixDiagonal_Jacobi(PC pc, PetscBool *flg)
136d71ae5a4SJacob Faibussowitsch {
13767ed0f3bSStefano Zampini PC_Jacobi *j = (PC_Jacobi *)pc->data;
13867ed0f3bSStefano Zampini
13967ed0f3bSStefano Zampini PetscFunctionBegin;
14067ed0f3bSStefano Zampini *flg = j->fixdiag;
1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14267ed0f3bSStefano Zampini }
14367ed0f3bSStefano Zampini
PCJacobiGetDiagonal_Jacobi(PC pc,Vec diag,Vec diagsqrt)144e03abaaeSJames Wright static PetscErrorCode PCJacobiGetDiagonal_Jacobi(PC pc, Vec diag, Vec diagsqrt)
145e03abaaeSJames Wright {
146e03abaaeSJames Wright PC_Jacobi *j = (PC_Jacobi *)pc->data;
147e03abaaeSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)pc);
148e03abaaeSJames Wright
149e03abaaeSJames Wright PetscFunctionBegin;
150e03abaaeSJames Wright PetscCheck(j->diag || j->diagsqrt, comm, PETSC_ERR_ARG_WRONGSTATE, "Jacobi diagonal has not been created yet. Use PCApply to force creation");
151e03abaaeSJames Wright PetscCheck(!diag || (diag && j->diag), comm, PETSC_ERR_ARG_WRONGSTATE, "Jacobi diagonal not available. Check if PC is non-symmetric");
152e03abaaeSJames Wright PetscCheck(!diagsqrt || (diagsqrt && j->diagsqrt), comm, PETSC_ERR_ARG_WRONGSTATE, "Jacobi diagonal squareroot not available. Check if PC is symmetric");
153e03abaaeSJames Wright
154e03abaaeSJames Wright if (diag) PetscCall(VecCopy(j->diag, diag));
155e03abaaeSJames Wright if (diagsqrt) PetscCall(VecCopy(j->diagsqrt, diagsqrt));
156e03abaaeSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
157e03abaaeSJames Wright }
158e03abaaeSJames Wright
1594b9ad928SBarry Smith /*
1604b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner
1614b9ad928SBarry Smith by setting data structures and options.
1624b9ad928SBarry Smith
1634b9ad928SBarry Smith Input Parameter:
1644b9ad928SBarry Smith . pc - the preconditioner context
1654b9ad928SBarry Smith
1664b9ad928SBarry Smith Application Interface Routine: PCSetUp()
1674b9ad928SBarry Smith
168f1580f4eSBarry Smith Note:
1694b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by
1704b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary.
1714b9ad928SBarry Smith */
PCSetUp_Jacobi(PC pc)172d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi(PC pc)
173d71ae5a4SJacob Faibussowitsch {
1744b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
1754b9ad928SBarry Smith Vec diag, diagsqrt;
176cd47f5d9SBarry Smith PetscInt n, i;
177eede4a3fSMark Adams PetscBool zeroflag = PETSC_FALSE, negflag = PETSC_FALSE;
1784b9ad928SBarry Smith
1794b9ad928SBarry Smith PetscFunctionBegin;
1804b9ad928SBarry Smith /*
1814b9ad928SBarry Smith For most preconditioners the code would begin here something like
1824b9ad928SBarry Smith
183371d2eb7SMartin Diehl if (!pc->setupcalled) { allocate space the first time this is ever called
1849566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->mat,&jac->diag));
1854b9ad928SBarry Smith }
1864b9ad928SBarry Smith
1874b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal
1884b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements
1894b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we
1904b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user
1914b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need
1924b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric()
1934b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively.
1944b9ad928SBarry Smith */
1954b9ad928SBarry Smith
1964b9ad928SBarry Smith /*
1974b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from
1984b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner.
1994b9ad928SBarry Smith */
2004b9ad928SBarry Smith diag = jac->diag;
2014b9ad928SBarry Smith diagsqrt = jac->diagsqrt;
2024b9ad928SBarry Smith
2034b9ad928SBarry Smith if (diag) {
204b94d7dedSBarry Smith PetscBool isset, isspd;
20567ed0f3bSStefano Zampini
206e03abaaeSJames Wright PetscCall(VecLockReadPop(diag));
207eede4a3fSMark Adams switch (jac->type) {
208eede4a3fSMark Adams case PC_JACOBI_DIAGONAL:
2099566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diag));
210eede4a3fSMark Adams break;
211eede4a3fSMark Adams case PC_JACOBI_ROWMAX:
212eede4a3fSMark Adams PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL));
213eede4a3fSMark Adams break;
214eede4a3fSMark Adams case PC_JACOBI_ROWL1:
215eede4a3fSMark Adams PetscCall(MatGetRowSumAbs(pc->pmat, diag));
216eede4a3fSMark Adams // fix negative rows (eg, negative definite) -- this could be done for all, not needed for userowmax
217eede4a3fSMark Adams PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
218eede4a3fSMark Adams if (jac->fixdiag && (!isset || !isspd)) {
219eede4a3fSMark Adams PetscScalar *x2;
220eede4a3fSMark Adams const PetscScalar *x;
221eede4a3fSMark Adams Vec true_diag;
222eede4a3fSMark Adams PetscCall(VecDuplicate(diag, &true_diag));
223eede4a3fSMark Adams PetscCall(MatGetDiagonal(pc->pmat, true_diag));
224eede4a3fSMark Adams PetscCall(VecGetLocalSize(diag, &n));
225eede4a3fSMark Adams PetscCall(VecGetArrayWrite(diag, &x2));
226eede4a3fSMark Adams PetscCall(VecGetArrayRead(true_diag, &x)); // to make more general -todo
227eede4a3fSMark Adams for (i = 0; i < n; i++) {
228eede4a3fSMark Adams if (PetscRealPart(x[i]) < 0.0) {
229eede4a3fSMark Adams x2[i] = -x2[i]; // flip sign to keep DA > 0
230eede4a3fSMark Adams negflag = PETSC_TRUE;
231eede4a3fSMark Adams }
232eede4a3fSMark Adams }
233eede4a3fSMark Adams PetscCall(VecRestoreArrayRead(true_diag, &x));
234eede4a3fSMark Adams PetscCall(VecRestoreArrayWrite(diag, &x2));
235eede4a3fSMark Adams PetscCheck(!jac->useabs || !negflag, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Jacobi use_abs and l1 not compatible with negative diagonal");
236eede4a3fSMark Adams PetscCall(VecDestroy(&true_diag));
237eede4a3fSMark Adams }
238eede4a3fSMark Adams if (jac->scale != 1.0) {
239eede4a3fSMark Adams Vec true_diag;
240eede4a3fSMark Adams PetscCall(VecDuplicate(diag, &true_diag));
241eede4a3fSMark Adams PetscCall(MatGetDiagonal(pc->pmat, true_diag));
242eede4a3fSMark Adams PetscCall(VecAXPY(diag, -1, true_diag)); // subtract off diag
243eede4a3fSMark Adams PetscCall(VecScale(diag, jac->scale)); // scale off-diag
244eede4a3fSMark Adams PetscCall(VecAXPY(diag, 1, true_diag)); // add diag back in
245eede4a3fSMark Adams PetscCall(VecDestroy(&true_diag));
246eede4a3fSMark Adams }
247eede4a3fSMark Adams break;
248eede4a3fSMark Adams case PC_JACOBI_ROWSUM:
249eede4a3fSMark Adams PetscCall(MatGetRowSum(pc->pmat, diag));
250eede4a3fSMark Adams break;
2514b9ad928SBarry Smith }
2529566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag));
2531baa6e33SBarry Smith if (jac->useabs) PetscCall(VecAbs(diag));
254b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd));
255b94d7dedSBarry Smith if (jac->fixdiag && (!isset || !isspd)) {
256eede4a3fSMark Adams PetscScalar *x;
2579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diag, &n));
2589566063dSJacob Faibussowitsch PetscCall(VecGetArray(diag, &x));
2594b9ad928SBarry Smith for (i = 0; i < n; i++) {
2604b9ad928SBarry Smith if (x[i] == 0.0) {
2614b9ad928SBarry Smith x[i] = 1.0;
262cd47f5d9SBarry Smith zeroflag = PETSC_TRUE;
2634b9ad928SBarry Smith }
2644b9ad928SBarry Smith }
2659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diag, &x));
2664b9ad928SBarry Smith }
267e03abaaeSJames Wright PetscCall(VecLockReadPush(diag));
26867ed0f3bSStefano Zampini }
2694b9ad928SBarry Smith if (diagsqrt) {
270eede4a3fSMark Adams PetscScalar *x;
271e03abaaeSJames Wright
272e03abaaeSJames Wright PetscCall(VecLockReadPop(diagsqrt));
273eede4a3fSMark Adams switch (jac->type) {
274eede4a3fSMark Adams case PC_JACOBI_DIAGONAL:
2759566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diagsqrt));
276eede4a3fSMark Adams break;
277eede4a3fSMark Adams case PC_JACOBI_ROWMAX:
278eede4a3fSMark Adams PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL));
279eede4a3fSMark Adams break;
280eede4a3fSMark Adams case PC_JACOBI_ROWL1:
281eede4a3fSMark Adams PetscCall(MatGetRowSumAbs(pc->pmat, diagsqrt));
282eede4a3fSMark Adams break;
283eede4a3fSMark Adams case PC_JACOBI_ROWSUM:
284eede4a3fSMark Adams PetscCall(MatGetRowSum(pc->pmat, diagsqrt));
285eede4a3fSMark Adams break;
2864b9ad928SBarry Smith }
2879566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diagsqrt, &n));
2889566063dSJacob Faibussowitsch PetscCall(VecGetArray(diagsqrt, &x));
2894b9ad928SBarry Smith for (i = 0; i < n; i++) {
290eede4a3fSMark Adams if (PetscRealPart(x[i]) < 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(-x[i]));
291eede4a3fSMark Adams else if (PetscRealPart(x[i]) > 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i]));
2924b9ad928SBarry Smith else {
2934b9ad928SBarry Smith x[i] = 1.0;
294cd47f5d9SBarry Smith zeroflag = PETSC_TRUE;
2954b9ad928SBarry Smith }
2964b9ad928SBarry Smith }
2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diagsqrt, &x));
298e03abaaeSJames Wright PetscCall(VecLockReadPush(diagsqrt));
2994b9ad928SBarry Smith }
30048a46eb9SPierre Jolivet if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n"));
3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3024b9ad928SBarry Smith }
303f1580f4eSBarry Smith
3044b9ad928SBarry Smith /*
3054b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the
3064b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This
3074b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner.
3084b9ad928SBarry Smith
3094b9ad928SBarry Smith Input Parameter:
3104b9ad928SBarry Smith . pc - the preconditioner context
3114b9ad928SBarry Smith */
PCSetUp_Jacobi_Symmetric(PC pc)312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc)
313d71ae5a4SJacob Faibussowitsch {
3144b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3154b9ad928SBarry Smith
3164b9ad928SBarry Smith PetscFunctionBegin;
3179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL));
318e03abaaeSJames Wright PetscCall(VecLockReadPush(jac->diagsqrt));
3199566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc));
3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3214b9ad928SBarry Smith }
322f1580f4eSBarry Smith
3234b9ad928SBarry Smith /*
3244b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the
3254b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of
3264b9ad928SBarry Smith right application of the Jacobi preconditioner.
3274b9ad928SBarry Smith
3284b9ad928SBarry Smith Input Parameter:
3294b9ad928SBarry Smith . pc - the preconditioner context
3304b9ad928SBarry Smith */
PCSetUp_Jacobi_NonSymmetric(PC pc)331d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc)
332d71ae5a4SJacob Faibussowitsch {
3334b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3344b9ad928SBarry Smith
3354b9ad928SBarry Smith PetscFunctionBegin;
3369566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL));
337e03abaaeSJames Wright PetscCall(VecLockReadPush(jac->diag));
3389566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc));
3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3404b9ad928SBarry Smith }
341f1580f4eSBarry Smith
3424b9ad928SBarry Smith /*
3434b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector.
3444b9ad928SBarry Smith
3454b9ad928SBarry Smith Input Parameters:
3464b9ad928SBarry Smith . pc - the preconditioner context
3474b9ad928SBarry Smith . x - input vector
3484b9ad928SBarry Smith
3494b9ad928SBarry Smith Output Parameter:
3504b9ad928SBarry Smith . y - output vector
3514b9ad928SBarry Smith
3524b9ad928SBarry Smith Application Interface Routine: PCApply()
3534b9ad928SBarry Smith */
PCApply_Jacobi(PC pc,Vec x,Vec y)354d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y)
355d71ae5a4SJacob Faibussowitsch {
3564b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3574b9ad928SBarry Smith
3584b9ad928SBarry Smith PetscFunctionBegin;
35948a46eb9SPierre Jolivet if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc));
3609566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diag));
3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3624b9ad928SBarry Smith }
363f1580f4eSBarry Smith
3644b9ad928SBarry Smith /*
3654b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a
3664b9ad928SBarry Smith symmetric preconditioner to a vector.
3674b9ad928SBarry Smith
3684b9ad928SBarry Smith Input Parameters:
3694b9ad928SBarry Smith . pc - the preconditioner context
3704b9ad928SBarry Smith . x - input vector
3714b9ad928SBarry Smith
3724b9ad928SBarry Smith Output Parameter:
3734b9ad928SBarry Smith . y - output vector
3744b9ad928SBarry Smith
3754b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight()
3764b9ad928SBarry Smith */
PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y)377d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y)
378d71ae5a4SJacob Faibussowitsch {
3794b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
3804b9ad928SBarry Smith
3814b9ad928SBarry Smith PetscFunctionBegin;
38248a46eb9SPierre Jolivet if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc));
3839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diagsqrt));
3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3854b9ad928SBarry Smith }
38667ed0f3bSStefano Zampini
PCReset_Jacobi(PC pc)387d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Jacobi(PC pc)
388d71ae5a4SJacob Faibussowitsch {
389a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
390a06653b4SBarry Smith
391a06653b4SBarry Smith PetscFunctionBegin;
392e03abaaeSJames Wright if (jac->diag) PetscCall(VecLockReadPop(jac->diag));
393e03abaaeSJames Wright if (jac->diagsqrt) PetscCall(VecLockReadPop(jac->diagsqrt));
3949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag));
3959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diagsqrt));
3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
397a06653b4SBarry Smith }
398a06653b4SBarry Smith
3994b9ad928SBarry Smith /*
4004b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner
4014b9ad928SBarry Smith that was created with PCCreate_Jacobi().
4024b9ad928SBarry Smith
4034b9ad928SBarry Smith Input Parameter:
4044b9ad928SBarry Smith . pc - the preconditioner context
4054b9ad928SBarry Smith
4064b9ad928SBarry Smith Application Interface Routine: PCDestroy()
4074b9ad928SBarry Smith */
PCDestroy_Jacobi(PC pc)408d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Jacobi(PC pc)
409d71ae5a4SJacob Faibussowitsch {
4104b9ad928SBarry Smith PetscFunctionBegin;
4119566063dSJacob Faibussowitsch PetscCall(PCReset_Jacobi(pc));
4129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL));
4139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL));
4149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL));
4159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL));
416eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetRowl1Scale_C", NULL));
417eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetRowl1Scale_C", NULL));
4189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL));
4199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL));
420e03abaaeSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetDiagonal_C", NULL));
4214b9ad928SBarry Smith
4224b9ad928SBarry Smith /*
4234b9ad928SBarry Smith Free the private data structure that was hanging off the PC
4244b9ad928SBarry Smith */
4259566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4274b9ad928SBarry Smith }
4284b9ad928SBarry Smith
PCSetFromOptions_Jacobi(PC pc,PetscOptionItems PetscOptionsObject)429ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems PetscOptionsObject)
430d71ae5a4SJacob Faibussowitsch {
4314b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data;
432baa89ecbSBarry Smith PetscBool flg;
433baa89ecbSBarry Smith PCJacobiType deflt, type;
4344b9ad928SBarry Smith
4354b9ad928SBarry Smith PetscFunctionBegin;
4369566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &deflt));
437d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options");
4389566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg));
4391baa6e33SBarry Smith if (flg) PetscCall(PCJacobiSetType(pc, type));
4409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL));
4419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL));
442eede4a3fSMark Adams PetscCall(PetscOptionsRangeReal("-pc_jacobi_rowl1_scale", "scaling of off-diagonal elements for rowl1", "PCJacobiSetRowl1Scale", jac->scale, &jac->scale, NULL, 0.0, 1.0));
443d0609cedSBarry Smith PetscOptionsHeadEnd();
4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4454b9ad928SBarry Smith }
4464b9ad928SBarry Smith
PCView_Jacobi(PC pc,PetscViewer viewer)447d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer)
448d71ae5a4SJacob Faibussowitsch {
449569e28a7SMatthew G. Knepley PC_Jacobi *jac = (PC_Jacobi *)pc->data;
450*9f196a02SMartin Diehl PetscBool isascii;
451569e28a7SMatthew G. Knepley
452569e28a7SMatthew G. Knepley PetscFunctionBegin;
453*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
454*9f196a02SMartin Diehl if (isascii) {
455569e28a7SMatthew G. Knepley PCJacobiType type;
45667ed0f3bSStefano Zampini PetscBool useAbs, fixdiag;
457569e28a7SMatthew G. Knepley PetscViewerFormat format;
458eede4a3fSMark Adams PetscReal scale;
459569e28a7SMatthew G. Knepley
4609566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &type));
4619566063dSJacob Faibussowitsch PetscCall(PCJacobiGetUseAbs(pc, &useAbs));
4629566063dSJacob Faibussowitsch PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag));
463eede4a3fSMark Adams PetscCall(PCJacobiGetRowl1Scale(pc, &scale));
464eede4a3fSMark Adams if (type == PC_JACOBI_ROWL1)
465eede4a3fSMark Adams PetscCall(PetscViewerASCIIPrintf(viewer, " type %s%s%s (l1-norm off-diagonal scaling %e)\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : "", (double)scale));
466eede4a3fSMark Adams else PetscCall(PetscViewerASCIIPrintf(viewer, " type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : ""));
4679566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format));
468b24842d1SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO_DETAIL && jac->diag) PetscCall(VecView(jac->diag, viewer));
469569e28a7SMatthew G. Knepley }
4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
471569e28a7SMatthew G. Knepley }
472569e28a7SMatthew G. Knepley
4734b9ad928SBarry Smith /*
4744b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi,
4754b9ad928SBarry Smith and sets this as the private data within the generic preconditioning
4764b9ad928SBarry Smith context, PC, that was created within PCCreate().
4774b9ad928SBarry Smith
4784b9ad928SBarry Smith Input Parameter:
4794b9ad928SBarry Smith . pc - the preconditioner context
4804b9ad928SBarry Smith
4814b9ad928SBarry Smith Application Interface Routine: PCCreate()
4824b9ad928SBarry Smith */
4834b9ad928SBarry Smith
4844b9ad928SBarry Smith /*MC
4855a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning)
4864b9ad928SBarry Smith
487f1580f4eSBarry Smith Options Database Keys:
488eede4a3fSMark Adams + -pc_jacobi_type <diagonal,rowl1,rowmax,rowsum> - approach for forming the preconditioner
48967ed0f3bSStefano Zampini . -pc_jacobi_abs - use the absolute value of the diagonal entry
490eede4a3fSMark Adams . -pc_jacobi_rowl1_scale - scaling of off-diagonal terms
491147403d9SBarry Smith - -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations
4924b9ad928SBarry Smith
4934b9ad928SBarry Smith Level: beginner
4944b9ad928SBarry Smith
49595452b02SPatrick Sanan Notes:
496f1580f4eSBarry Smith By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric
4974b9ad928SBarry Smith can scale each side of the matrix by the square root of the diagonal entries.
4984b9ad928SBarry Smith
4994b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0
5004b9ad928SBarry Smith
501f1580f4eSBarry Smith See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks
502422a814eSBarry Smith
503db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
504f1580f4eSBarry Smith `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`,
505db781477SPatrick Sanan `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()`
506db781477SPatrick Sanan `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI`
5074b9ad928SBarry Smith M*/
5084b9ad928SBarry Smith
PCCreate_Jacobi(PC pc)509d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc)
510d71ae5a4SJacob Faibussowitsch {
5114b9ad928SBarry Smith PC_Jacobi *jac;
5124b9ad928SBarry Smith
5134b9ad928SBarry Smith PetscFunctionBegin;
5144b9ad928SBarry Smith /*
5154b9ad928SBarry Smith Creates the private data structure for this preconditioner and
5164b9ad928SBarry Smith attach it to the PC object.
5174b9ad928SBarry Smith */
5184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac));
5194b9ad928SBarry Smith pc->data = (void *)jac;
5204b9ad928SBarry Smith
5214b9ad928SBarry Smith /*
5224b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store
5234b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application.
5244b9ad928SBarry Smith */
5250a545947SLisandro Dalcin jac->diag = NULL;
5260a545947SLisandro Dalcin jac->diagsqrt = NULL;
527eede4a3fSMark Adams jac->type = PC_JACOBI_DIAGONAL;
528cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE;
52967ed0f3bSStefano Zampini jac->fixdiag = PETSC_TRUE;
530eede4a3fSMark Adams jac->scale = 1.0;
5314b9ad928SBarry Smith
5324b9ad928SBarry Smith /*
5334b9ad928SBarry Smith Set the pointers for the functions that are provided above.
5344b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
5354b9ad928SBarry Smith are called, they will automatically call these functions. Note we
5364b9ad928SBarry Smith choose not to provide a couple of these functions since they are
5374b9ad928SBarry Smith not needed.
5384b9ad928SBarry Smith */
5394b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi;
5404b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi;
5414b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi;
542a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi;
5434b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi;
5444b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi;
545569e28a7SMatthew G. Knepley pc->ops->view = PCView_Jacobi;
5460a545947SLisandro Dalcin pc->ops->applyrichardson = NULL;
5474b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi;
5484b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi;
5492fa5cd67SKarl Rupp
5509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi));
5519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi));
552eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetRowl1Scale_C", PCJacobiSetRowl1Scale_Jacobi));
553eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetRowl1Scale_C", PCJacobiGetRowl1Scale_Jacobi));
5549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi));
5559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi));
5569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi));
5579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi));
558e03abaaeSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetDiagonal_C", PCJacobiGetDiagonal_Jacobi));
5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5604b9ad928SBarry Smith }
561cd47f5d9SBarry Smith
562cd47f5d9SBarry Smith /*@
563f1580f4eSBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the
56467ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner
565cd47f5d9SBarry Smith
566c3339decSBarry Smith Logically Collective
567cd47f5d9SBarry Smith
568cd47f5d9SBarry Smith Input Parameters:
569baa89ecbSBarry Smith + pc - the preconditioner context
570baa89ecbSBarry Smith - flg - whether to use absolute values or not
571baa89ecbSBarry Smith
572baa89ecbSBarry Smith Options Database Key:
573147403d9SBarry Smith . -pc_jacobi_abs <bool> - use absolute values
574baa89ecbSBarry Smith
575f1580f4eSBarry Smith Note:
57695452b02SPatrick Sanan This takes affect at the next construction of the preconditioner
577baa89ecbSBarry Smith
578baa89ecbSBarry Smith Level: intermediate
579baa89ecbSBarry Smith
580eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetUseAbs()`
581baa89ecbSBarry Smith @*/
PCJacobiSetUseAbs(PC pc,PetscBool flg)582d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg)
583d71ae5a4SJacob Faibussowitsch {
584baa89ecbSBarry Smith PetscFunctionBegin;
585baa89ecbSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
586cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg));
5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
588baa89ecbSBarry Smith }
589baa89ecbSBarry Smith
590baa89ecbSBarry Smith /*@
591f1580f4eSBarry Smith PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the
59267ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner
593baa89ecbSBarry Smith
594c3339decSBarry Smith Logically Collective
595baa89ecbSBarry Smith
596baa89ecbSBarry Smith Input Parameter:
597cd47f5d9SBarry Smith . pc - the preconditioner context
598cd47f5d9SBarry Smith
599baa89ecbSBarry Smith Output Parameter:
600baa89ecbSBarry Smith . flg - whether to use absolute values or not
601cd47f5d9SBarry Smith
602cd47f5d9SBarry Smith Level: intermediate
603cd47f5d9SBarry Smith
604eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
605cd47f5d9SBarry Smith @*/
PCJacobiGetUseAbs(PC pc,PetscBool * flg)606d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg)
607d71ae5a4SJacob Faibussowitsch {
608cd47f5d9SBarry Smith PetscFunctionBegin;
6090700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
610cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg));
6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
612cd47f5d9SBarry Smith }
613cd47f5d9SBarry Smith
6144b9ad928SBarry Smith /*@
615eede4a3fSMark Adams PCJacobiSetRowl1Scale - Set scaling of off-diagonal of operator when computing l1 row norms, eg,
616eede4a3fSMark Adams Remark 6.1 in "Multigrid Smoothers for Ultraparallel Computing", Baker et al, with 0.5 scaling
617eede4a3fSMark Adams
618eede4a3fSMark Adams Logically Collective
619eede4a3fSMark Adams
620eede4a3fSMark Adams Input Parameters:
621eede4a3fSMark Adams + pc - the preconditioner context
622eede4a3fSMark Adams - scale - scaling
623eede4a3fSMark Adams
624eede4a3fSMark Adams Options Database Key:
625eede4a3fSMark Adams . -pc_jacobi_rowl1_scale <real> - use absolute values
626eede4a3fSMark Adams
627eede4a3fSMark Adams Level: intermediate
628eede4a3fSMark Adams
629eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetRowl1Scale()`
630eede4a3fSMark Adams @*/
PCJacobiSetRowl1Scale(PC pc,PetscReal scale)631eede4a3fSMark Adams PetscErrorCode PCJacobiSetRowl1Scale(PC pc, PetscReal scale)
632eede4a3fSMark Adams {
633eede4a3fSMark Adams PetscFunctionBegin;
634eede4a3fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
635eede4a3fSMark Adams PetscTryMethod(pc, "PCJacobiSetRowl1Scale_C", (PC, PetscReal), (pc, scale));
636eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
637eede4a3fSMark Adams }
638eede4a3fSMark Adams
639eede4a3fSMark Adams /*@
640eede4a3fSMark Adams PCJacobiGetRowl1Scale - Get scaling of off-diagonal elements summed into l1-norm diagonal
641eede4a3fSMark Adams
642eede4a3fSMark Adams Logically Collective
643eede4a3fSMark Adams
644eede4a3fSMark Adams Input Parameter:
645eede4a3fSMark Adams . pc - the preconditioner context
646eede4a3fSMark Adams
647eede4a3fSMark Adams Output Parameter:
648eede4a3fSMark Adams . scale - scaling
649eede4a3fSMark Adams
650eede4a3fSMark Adams Level: intermediate
651eede4a3fSMark Adams
652eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetRowl1Scale()`, `PCJacobiGetType()`
653eede4a3fSMark Adams @*/
PCJacobiGetRowl1Scale(PC pc,PetscReal * scale)654eede4a3fSMark Adams PetscErrorCode PCJacobiGetRowl1Scale(PC pc, PetscReal *scale)
655eede4a3fSMark Adams {
656eede4a3fSMark Adams PetscFunctionBegin;
657eede4a3fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
658eede4a3fSMark Adams PetscUseMethod(pc, "PCJacobiGetRowl1Scale_C", (PC, PetscReal *), (pc, scale));
659eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
660eede4a3fSMark Adams }
661eede4a3fSMark Adams
662eede4a3fSMark Adams /*@
663147403d9SBarry Smith PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0
66467ed0f3bSStefano Zampini
665c3339decSBarry Smith Logically Collective
66667ed0f3bSStefano Zampini
66767ed0f3bSStefano Zampini Input Parameters:
66867ed0f3bSStefano Zampini + pc - the preconditioner context
66967ed0f3bSStefano Zampini - flg - the boolean flag
67067ed0f3bSStefano Zampini
67167ed0f3bSStefano Zampini Options Database Key:
672147403d9SBarry Smith . -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal
67367ed0f3bSStefano Zampini
674f1580f4eSBarry Smith Note:
67567ed0f3bSStefano Zampini This takes affect at the next construction of the preconditioner
67667ed0f3bSStefano Zampini
67767ed0f3bSStefano Zampini Level: intermediate
67867ed0f3bSStefano Zampini
679562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()`
68067ed0f3bSStefano Zampini @*/
PCJacobiSetFixDiagonal(PC pc,PetscBool flg)681d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg)
682d71ae5a4SJacob Faibussowitsch {
68367ed0f3bSStefano Zampini PetscFunctionBegin;
68467ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
685cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg));
6863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
68767ed0f3bSStefano Zampini }
68867ed0f3bSStefano Zampini
68967ed0f3bSStefano Zampini /*@
690f1580f4eSBarry Smith PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms
69167ed0f3bSStefano Zampini
692c3339decSBarry Smith Logically Collective
69367ed0f3bSStefano Zampini
69467ed0f3bSStefano Zampini Input Parameter:
69567ed0f3bSStefano Zampini . pc - the preconditioner context
69667ed0f3bSStefano Zampini
69767ed0f3bSStefano Zampini Output Parameter:
69867ed0f3bSStefano Zampini . flg - the boolean flag
69967ed0f3bSStefano Zampini
70067ed0f3bSStefano Zampini Options Database Key:
70167b8a455SSatish Balay . -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1
70267ed0f3bSStefano Zampini
70367ed0f3bSStefano Zampini Level: intermediate
70467ed0f3bSStefano Zampini
705562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()`
70667ed0f3bSStefano Zampini @*/
PCJacobiGetFixDiagonal(PC pc,PetscBool * flg)707d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg)
708d71ae5a4SJacob Faibussowitsch {
70967ed0f3bSStefano Zampini PetscFunctionBegin;
71067ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
711cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg));
7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
71367ed0f3bSStefano Zampini }
71467ed0f3bSStefano Zampini
71567ed0f3bSStefano Zampini /*@
716e03abaaeSJames Wright PCJacobiGetDiagonal - Returns copy of the diagonal and/or diagonal squareroot `Vec`
717e03abaaeSJames Wright
718e03abaaeSJames Wright Logically Collective
719e03abaaeSJames Wright
720e03abaaeSJames Wright Input Parameter:
721e03abaaeSJames Wright . pc - the preconditioner context
722e03abaaeSJames Wright
723e03abaaeSJames Wright Output Parameters:
724e03abaaeSJames Wright + diagonal - Copy of `Vec` of the inverted diagonal
725e03abaaeSJames Wright - diagonal_sqrt - Copy of `Vec` of the inverted square root diagonal
726e03abaaeSJames Wright
727e03abaaeSJames Wright Level: developer
728e03abaaeSJames Wright
729e03abaaeSJames Wright .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`
730e03abaaeSJames Wright @*/
PCJacobiGetDiagonal(PC pc,Vec diagonal,Vec diagonal_sqrt)731e03abaaeSJames Wright PetscErrorCode PCJacobiGetDiagonal(PC pc, Vec diagonal, Vec diagonal_sqrt)
732e03abaaeSJames Wright {
733e03abaaeSJames Wright PetscFunctionBegin;
734e03abaaeSJames Wright PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
735e03abaaeSJames Wright PetscUseMethod(pc, "PCJacobiGetDiagonal_C", (PC, Vec, Vec), (pc, diagonal, diagonal_sqrt));
736e03abaaeSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
737e03abaaeSJames Wright }
738e03abaaeSJames Wright
739e03abaaeSJames Wright /*@
740baa89ecbSBarry Smith PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row,
741baa89ecbSBarry Smith of the sum of rows entries for the diagonal preconditioner
7424b9ad928SBarry Smith
743c3339decSBarry Smith Logically Collective
7444b9ad928SBarry Smith
7454b9ad928SBarry Smith Input Parameters:
746baa89ecbSBarry Smith + pc - the preconditioner context
747eede4a3fSMark Adams - type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWL1`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
7484b9ad928SBarry Smith
7494b9ad928SBarry Smith Options Database Key:
750eede4a3fSMark Adams . -pc_jacobi_type <diagonal,rowl1,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi
7514b9ad928SBarry Smith
7524b9ad928SBarry Smith Level: intermediate
7534b9ad928SBarry Smith
754feefa0e1SJacob Faibussowitsch Developer Notes:
755f1580f4eSBarry Smith Why is there a separate function for using the absolute value?
756f1580f4eSBarry Smith
757562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()`
7584b9ad928SBarry Smith @*/
PCJacobiSetType(PC pc,PCJacobiType type)759d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type)
760d71ae5a4SJacob Faibussowitsch {
7614b9ad928SBarry Smith PetscFunctionBegin;
7620700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
763cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type));
7643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7654b9ad928SBarry Smith }
7664b9ad928SBarry Smith
76786697f06SMatthew Knepley /*@
768baa89ecbSBarry Smith PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner
76986697f06SMatthew Knepley
770c3339decSBarry Smith Not Collective
77186697f06SMatthew Knepley
772baa89ecbSBarry Smith Input Parameter:
77386697f06SMatthew Knepley . pc - the preconditioner context
77486697f06SMatthew Knepley
775baa89ecbSBarry Smith Output Parameter:
776eede4a3fSMark Adams . type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWL1`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM`
77786697f06SMatthew Knepley
77886697f06SMatthew Knepley Level: intermediate
77986697f06SMatthew Knepley
78054c05997SPierre Jolivet .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiSetType()`
78186697f06SMatthew Knepley @*/
PCJacobiGetType(PC pc,PCJacobiType * type)782d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type)
783d71ae5a4SJacob Faibussowitsch {
78486697f06SMatthew Knepley PetscFunctionBegin;
7850700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
786cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type));
7873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
78886697f06SMatthew Knepley }
789