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 50*eede4a3fSMark 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 { 562fa5cd67SKarl Rupp Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 574b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 584b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 594b9ad928SBarry Smith only for symmetric preconditioner application) */ 60*eede4a3fSMark Adams PCJacobiType type; 61ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */ 6267ed0f3bSStefano Zampini PetscBool fixdiag; /* fix zero diagonal terms */ 63*eede4a3fSMark Adams PetscReal scale; /* for scaling rowl1 off-diagonals */ 644b9ad928SBarry Smith } PC_Jacobi; 654b9ad928SBarry Smith 667a660b1eSJed Brown static PetscErrorCode PCReset_Jacobi(PC); 677a660b1eSJed Brown 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)); 77*eede4a3fSMark Adams j->type = type; 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 794b9ad928SBarry Smith } 804b9ad928SBarry Smith 81*eede4a3fSMark 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; 86*eede4a3fSMark Adams *flg = j->useabs; 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8886697f06SMatthew Knepley } 8986697f06SMatthew Knepley 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 99*eede4a3fSMark 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; 104*eede4a3fSMark Adams *type = j->type; 105*eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS); 106*eede4a3fSMark Adams } 107*eede4a3fSMark Adams 108*eede4a3fSMark Adams static PetscErrorCode PCJacobiSetRowl1Scale_Jacobi(PC pc, PetscReal flg) 109*eede4a3fSMark Adams { 110*eede4a3fSMark Adams PC_Jacobi *j = (PC_Jacobi *)pc->data; 111*eede4a3fSMark Adams 112*eede4a3fSMark Adams PetscFunctionBegin; 113*eede4a3fSMark Adams j->scale = flg; 114*eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS); 115*eede4a3fSMark Adams } 116*eede4a3fSMark Adams 117*eede4a3fSMark Adams static PetscErrorCode PCJacobiGetRowl1Scale_Jacobi(PC pc, PetscReal *flg) 118*eede4a3fSMark Adams { 119*eede4a3fSMark Adams PC_Jacobi *j = (PC_Jacobi *)pc->data; 120*eede4a3fSMark Adams 121*eede4a3fSMark Adams PetscFunctionBegin; 122*eede4a3fSMark Adams *flg = j->scale; 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124cd47f5d9SBarry Smith } 125cd47f5d9SBarry Smith 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 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 1444b9ad928SBarry Smith /* 1454b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1464b9ad928SBarry Smith by setting data structures and options. 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith Input Parameter: 1494b9ad928SBarry Smith . pc - the preconditioner context 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1524b9ad928SBarry Smith 153f1580f4eSBarry Smith Note: 1544b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1554b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1564b9ad928SBarry Smith */ 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi(PC pc) 158d71ae5a4SJacob Faibussowitsch { 1594b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 1604b9ad928SBarry Smith Vec diag, diagsqrt; 161cd47f5d9SBarry Smith PetscInt n, i; 162*eede4a3fSMark Adams PetscBool zeroflag = PETSC_FALSE, negflag = PETSC_FALSE; 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith PetscFunctionBegin; 1654b9ad928SBarry Smith /* 1664b9ad928SBarry Smith For most preconditioners the code would begin here something like 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 1699566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->mat,&jac->diag)); 1704b9ad928SBarry Smith } 1714b9ad928SBarry Smith 1724b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1734b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1744b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1754b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1764b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1774b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1784b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1794b9ad928SBarry Smith */ 1804b9ad928SBarry Smith 1814b9ad928SBarry Smith /* 1824b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1834b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1844b9ad928SBarry Smith */ 1854b9ad928SBarry Smith diag = jac->diag; 1864b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith if (diag) { 189b94d7dedSBarry Smith PetscBool isset, isspd; 19067ed0f3bSStefano Zampini 191*eede4a3fSMark Adams switch (jac->type) { 192*eede4a3fSMark Adams case PC_JACOBI_DIAGONAL: 1939566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diag)); 194*eede4a3fSMark Adams break; 195*eede4a3fSMark Adams case PC_JACOBI_ROWMAX: 196*eede4a3fSMark Adams PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL)); 197*eede4a3fSMark Adams break; 198*eede4a3fSMark Adams case PC_JACOBI_ROWL1: 199*eede4a3fSMark Adams PetscCall(MatGetRowSumAbs(pc->pmat, diag)); 200*eede4a3fSMark Adams // fix negative rows (eg, negative definite) -- this could be done for all, not needed for userowmax 201*eede4a3fSMark Adams PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd)); 202*eede4a3fSMark Adams if (jac->fixdiag && (!isset || !isspd)) { 203*eede4a3fSMark Adams PetscScalar *x2; 204*eede4a3fSMark Adams const PetscScalar *x; 205*eede4a3fSMark Adams Vec true_diag; 206*eede4a3fSMark Adams PetscCall(VecDuplicate(diag, &true_diag)); 207*eede4a3fSMark Adams PetscCall(MatGetDiagonal(pc->pmat, true_diag)); 208*eede4a3fSMark Adams PetscCall(VecGetLocalSize(diag, &n)); 209*eede4a3fSMark Adams PetscCall(VecGetArrayWrite(diag, &x2)); 210*eede4a3fSMark Adams PetscCall(VecGetArrayRead(true_diag, &x)); // to make more general -todo 211*eede4a3fSMark Adams for (i = 0; i < n; i++) { 212*eede4a3fSMark Adams if (PetscRealPart(x[i]) < 0.0) { 213*eede4a3fSMark Adams x2[i] = -x2[i]; // flip sign to keep DA > 0 214*eede4a3fSMark Adams negflag = PETSC_TRUE; 215*eede4a3fSMark Adams } 216*eede4a3fSMark Adams } 217*eede4a3fSMark Adams PetscCall(VecRestoreArrayRead(true_diag, &x)); 218*eede4a3fSMark Adams PetscCall(VecRestoreArrayWrite(diag, &x2)); 219*eede4a3fSMark Adams PetscCheck(!jac->useabs || !negflag, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Jacobi use_abs and l1 not compatible with negative diagonal"); 220*eede4a3fSMark Adams PetscCall(VecDestroy(&true_diag)); 221*eede4a3fSMark Adams } 222*eede4a3fSMark Adams if (jac->scale != 1.0) { 223*eede4a3fSMark Adams Vec true_diag; 224*eede4a3fSMark Adams PetscCall(VecDuplicate(diag, &true_diag)); 225*eede4a3fSMark Adams PetscCall(MatGetDiagonal(pc->pmat, true_diag)); 226*eede4a3fSMark Adams PetscCall(VecAXPY(diag, -1, true_diag)); // subtract off diag 227*eede4a3fSMark Adams PetscCall(VecScale(diag, jac->scale)); // scale off-diag 228*eede4a3fSMark Adams PetscCall(VecAXPY(diag, 1, true_diag)); // add diag back in 229*eede4a3fSMark Adams PetscCall(VecDestroy(&true_diag)); 230*eede4a3fSMark Adams } 231*eede4a3fSMark Adams break; 232*eede4a3fSMark Adams case PC_JACOBI_ROWSUM: 233*eede4a3fSMark Adams PetscCall(MatGetRowSum(pc->pmat, diag)); 234*eede4a3fSMark Adams break; 2354b9ad928SBarry Smith } 2369566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 2371baa6e33SBarry Smith if (jac->useabs) PetscCall(VecAbs(diag)); 238b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd)); 239b94d7dedSBarry Smith if (jac->fixdiag && (!isset || !isspd)) { 240*eede4a3fSMark Adams PetscScalar *x; 2419566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diag, &n)); 2429566063dSJacob Faibussowitsch PetscCall(VecGetArray(diag, &x)); 2434b9ad928SBarry Smith for (i = 0; i < n; i++) { 2444b9ad928SBarry Smith if (x[i] == 0.0) { 2454b9ad928SBarry Smith x[i] = 1.0; 246cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 2474b9ad928SBarry Smith } 2484b9ad928SBarry Smith } 2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diag, &x)); 2504b9ad928SBarry Smith } 25167ed0f3bSStefano Zampini } 2524b9ad928SBarry Smith if (diagsqrt) { 253*eede4a3fSMark Adams PetscScalar *x; 254*eede4a3fSMark Adams switch (jac->type) { 255*eede4a3fSMark Adams case PC_JACOBI_DIAGONAL: 2569566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diagsqrt)); 257*eede4a3fSMark Adams break; 258*eede4a3fSMark Adams case PC_JACOBI_ROWMAX: 259*eede4a3fSMark Adams PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL)); 260*eede4a3fSMark Adams break; 261*eede4a3fSMark Adams case PC_JACOBI_ROWL1: 262*eede4a3fSMark Adams PetscCall(MatGetRowSumAbs(pc->pmat, diagsqrt)); 263*eede4a3fSMark Adams break; 264*eede4a3fSMark Adams case PC_JACOBI_ROWSUM: 265*eede4a3fSMark Adams PetscCall(MatGetRowSum(pc->pmat, diagsqrt)); 266*eede4a3fSMark Adams break; 2674b9ad928SBarry Smith } 2689566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diagsqrt, &n)); 2699566063dSJacob Faibussowitsch PetscCall(VecGetArray(diagsqrt, &x)); 2704b9ad928SBarry Smith for (i = 0; i < n; i++) { 271*eede4a3fSMark Adams if (PetscRealPart(x[i]) < 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(-x[i])); 272*eede4a3fSMark Adams else if (PetscRealPart(x[i]) > 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i])); 2734b9ad928SBarry Smith else { 2744b9ad928SBarry Smith x[i] = 1.0; 275cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 2764b9ad928SBarry Smith } 2774b9ad928SBarry Smith } 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diagsqrt, &x)); 2794b9ad928SBarry Smith } 28048a46eb9SPierre Jolivet if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n")); 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2824b9ad928SBarry Smith } 283f1580f4eSBarry Smith 2844b9ad928SBarry Smith /* 2854b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2864b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2874b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith Input Parameter: 2904b9ad928SBarry Smith . pc - the preconditioner context 2914b9ad928SBarry Smith */ 292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 293d71ae5a4SJacob Faibussowitsch { 2944b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL)); 2989566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc)); 2993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3004b9ad928SBarry Smith } 301f1580f4eSBarry Smith 3024b9ad928SBarry Smith /* 3034b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 3044b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 3054b9ad928SBarry Smith right application of the Jacobi preconditioner. 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith Input Parameter: 3084b9ad928SBarry Smith . pc - the preconditioner context 3094b9ad928SBarry Smith */ 310d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 311d71ae5a4SJacob Faibussowitsch { 3124b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith PetscFunctionBegin; 3159566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL)); 3169566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc)); 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3184b9ad928SBarry Smith } 319f1580f4eSBarry Smith 3204b9ad928SBarry Smith /* 3214b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith Input Parameters: 3244b9ad928SBarry Smith . pc - the preconditioner context 3254b9ad928SBarry Smith . x - input vector 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith Output Parameter: 3284b9ad928SBarry Smith . y - output vector 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Application Interface Routine: PCApply() 3314b9ad928SBarry Smith */ 332d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y) 333d71ae5a4SJacob Faibussowitsch { 3344b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith PetscFunctionBegin; 33748a46eb9SPierre Jolivet if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc)); 3389566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diag)); 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3404b9ad928SBarry Smith } 341f1580f4eSBarry Smith 3424b9ad928SBarry Smith /* 3434b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 3444b9ad928SBarry Smith symmetric preconditioner to a vector. 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith Input Parameters: 3474b9ad928SBarry Smith . pc - the preconditioner context 3484b9ad928SBarry Smith . x - input vector 3494b9ad928SBarry Smith 3504b9ad928SBarry Smith Output Parameter: 3514b9ad928SBarry Smith . y - output vector 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 3544b9ad928SBarry Smith */ 355d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y) 356d71ae5a4SJacob Faibussowitsch { 3574b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith PetscFunctionBegin; 36048a46eb9SPierre Jolivet if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc)); 3619566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diagsqrt)); 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3634b9ad928SBarry Smith } 36467ed0f3bSStefano Zampini 365d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Jacobi(PC pc) 366d71ae5a4SJacob Faibussowitsch { 367a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 368a06653b4SBarry Smith 369a06653b4SBarry Smith PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag)); 3719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diagsqrt)); 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373a06653b4SBarry Smith } 374a06653b4SBarry Smith 3754b9ad928SBarry Smith /* 3764b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3774b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3784b9ad928SBarry Smith 3794b9ad928SBarry Smith Input Parameter: 3804b9ad928SBarry Smith . pc - the preconditioner context 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3834b9ad928SBarry Smith */ 384d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Jacobi(PC pc) 385d71ae5a4SJacob Faibussowitsch { 3864b9ad928SBarry Smith PetscFunctionBegin; 3879566063dSJacob Faibussowitsch PetscCall(PCReset_Jacobi(pc)); 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL)); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL)); 3909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL)); 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL)); 392*eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetRowl1Scale_C", NULL)); 393*eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetRowl1Scale_C", NULL)); 3949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL)); 3959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL)); 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith /* 3984b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3994b9ad928SBarry Smith */ 4009566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4024b9ad928SBarry Smith } 4034b9ad928SBarry Smith 404d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject) 405d71ae5a4SJacob Faibussowitsch { 4064b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 407baa89ecbSBarry Smith PetscBool flg; 408baa89ecbSBarry Smith PCJacobiType deflt, type; 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &deflt)); 412d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options"); 4139566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg)); 4141baa6e33SBarry Smith if (flg) PetscCall(PCJacobiSetType(pc, type)); 4159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL)); 4169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL)); 417*eede4a3fSMark Adams PetscCall(PetscOptionsRangeReal("-pc_jacobi_rowl1_scale", "scaling of off-diagonal elements for rowl1", "PCJacobiSetRowl1Scale", jac->scale, &jac->scale, NULL, 0.0, 1.0)); 418d0609cedSBarry Smith PetscOptionsHeadEnd(); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4204b9ad928SBarry Smith } 4214b9ad928SBarry Smith 422d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer) 423d71ae5a4SJacob Faibussowitsch { 424569e28a7SMatthew G. Knepley PC_Jacobi *jac = (PC_Jacobi *)pc->data; 425569e28a7SMatthew G. Knepley PetscBool iascii; 426569e28a7SMatthew G. Knepley 427569e28a7SMatthew G. Knepley PetscFunctionBegin; 4289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 429569e28a7SMatthew G. Knepley if (iascii) { 430569e28a7SMatthew G. Knepley PCJacobiType type; 43167ed0f3bSStefano Zampini PetscBool useAbs, fixdiag; 432569e28a7SMatthew G. Knepley PetscViewerFormat format; 433*eede4a3fSMark Adams PetscReal scale; 434569e28a7SMatthew G. Knepley 4359566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &type)); 4369566063dSJacob Faibussowitsch PetscCall(PCJacobiGetUseAbs(pc, &useAbs)); 4379566063dSJacob Faibussowitsch PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag)); 438*eede4a3fSMark Adams PetscCall(PCJacobiGetRowl1Scale(pc, &scale)); 439*eede4a3fSMark Adams if (type == PC_JACOBI_ROWL1) 440*eede4a3fSMark 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)); 441*eede4a3fSMark Adams else PetscCall(PetscViewerASCIIPrintf(viewer, " type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : "")); 4429566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 443b24842d1SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO_DETAIL && jac->diag) PetscCall(VecView(jac->diag, viewer)); 444569e28a7SMatthew G. Knepley } 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 446569e28a7SMatthew G. Knepley } 447569e28a7SMatthew G. Knepley 4484b9ad928SBarry Smith /* 4494b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 4504b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 4514b9ad928SBarry Smith context, PC, that was created within PCCreate(). 4524b9ad928SBarry Smith 4534b9ad928SBarry Smith Input Parameter: 4544b9ad928SBarry Smith . pc - the preconditioner context 4554b9ad928SBarry Smith 4564b9ad928SBarry Smith Application Interface Routine: PCCreate() 4574b9ad928SBarry Smith */ 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith /*MC 4605a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 4614b9ad928SBarry Smith 462f1580f4eSBarry Smith Options Database Keys: 463*eede4a3fSMark Adams + -pc_jacobi_type <diagonal,rowl1,rowmax,rowsum> - approach for forming the preconditioner 46467ed0f3bSStefano Zampini . -pc_jacobi_abs - use the absolute value of the diagonal entry 465*eede4a3fSMark Adams . -pc_jacobi_rowl1_scale - scaling of off-diagonal terms 466147403d9SBarry Smith - -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations 4674b9ad928SBarry Smith 4684b9ad928SBarry Smith Level: beginner 4694b9ad928SBarry Smith 47095452b02SPatrick Sanan Notes: 471f1580f4eSBarry Smith By using `KSPSetPCSide`(ksp,`PC_SYMMETRIC`) or -ksp_pc_side symmetric 4724b9ad928SBarry Smith can scale each side of the matrix by the square root of the diagonal entries. 4734b9ad928SBarry Smith 4744b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 4754b9ad928SBarry Smith 476f1580f4eSBarry Smith See `PCPBJACOBI` for fixed-size point block, `PCVPBJACOBI` for variable-sized point block, and `PCBJACOBI` for large size blocks 477422a814eSBarry Smith 478db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 479f1580f4eSBarry Smith `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCASM`, 480db781477SPatrick Sanan `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()` 481db781477SPatrick Sanan `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI` 4824b9ad928SBarry Smith M*/ 4834b9ad928SBarry Smith 484d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) 485d71ae5a4SJacob Faibussowitsch { 4864b9ad928SBarry Smith PC_Jacobi *jac; 4874b9ad928SBarry Smith 4884b9ad928SBarry Smith PetscFunctionBegin; 4894b9ad928SBarry Smith /* 4904b9ad928SBarry Smith Creates the private data structure for this preconditioner and 4914b9ad928SBarry Smith attach it to the PC object. 4924b9ad928SBarry Smith */ 4934dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 4944b9ad928SBarry Smith pc->data = (void *)jac; 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith /* 4974b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4984b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4994b9ad928SBarry Smith */ 5000a545947SLisandro Dalcin jac->diag = NULL; 5010a545947SLisandro Dalcin jac->diagsqrt = NULL; 502*eede4a3fSMark Adams jac->type = PC_JACOBI_DIAGONAL; 503cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 50467ed0f3bSStefano Zampini jac->fixdiag = PETSC_TRUE; 505*eede4a3fSMark Adams jac->scale = 1.0; 5064b9ad928SBarry Smith 5074b9ad928SBarry Smith /* 5084b9ad928SBarry Smith Set the pointers for the functions that are provided above. 5094b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 5104b9ad928SBarry Smith are called, they will automatically call these functions. Note we 5114b9ad928SBarry Smith choose not to provide a couple of these functions since they are 5124b9ad928SBarry Smith not needed. 5134b9ad928SBarry Smith */ 5144b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 5154b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 5164b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 517a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi; 5184b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 5194b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 520569e28a7SMatthew G. Knepley pc->ops->view = PCView_Jacobi; 5210a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 5224b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 5234b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 5242fa5cd67SKarl Rupp 5259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi)); 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi)); 527*eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetRowl1Scale_C", PCJacobiSetRowl1Scale_Jacobi)); 528*eede4a3fSMark Adams PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetRowl1Scale_C", PCJacobiGetRowl1Scale_Jacobi)); 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi)); 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi)); 5319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi)); 5329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi)); 5333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5344b9ad928SBarry Smith } 535cd47f5d9SBarry Smith 536cd47f5d9SBarry Smith /*@ 537f1580f4eSBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner `PCJACOBI` to use the 53867ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner 539cd47f5d9SBarry Smith 540c3339decSBarry Smith Logically Collective 541cd47f5d9SBarry Smith 542cd47f5d9SBarry Smith Input Parameters: 543baa89ecbSBarry Smith + pc - the preconditioner context 544baa89ecbSBarry Smith - flg - whether to use absolute values or not 545baa89ecbSBarry Smith 546baa89ecbSBarry Smith Options Database Key: 547147403d9SBarry Smith . -pc_jacobi_abs <bool> - use absolute values 548baa89ecbSBarry Smith 549f1580f4eSBarry Smith Note: 55095452b02SPatrick Sanan This takes affect at the next construction of the preconditioner 551baa89ecbSBarry Smith 552baa89ecbSBarry Smith Level: intermediate 553baa89ecbSBarry Smith 554*eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetUseAbs()` 555baa89ecbSBarry Smith @*/ 556d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg) 557d71ae5a4SJacob Faibussowitsch { 558baa89ecbSBarry Smith PetscFunctionBegin; 559baa89ecbSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 560cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg)); 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 562baa89ecbSBarry Smith } 563baa89ecbSBarry Smith 564baa89ecbSBarry Smith /*@ 565f1580f4eSBarry Smith PCJacobiGetUseAbs - Determines if the Jacobi preconditioner `PCJACOBI` uses the 56667ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner 567baa89ecbSBarry Smith 568c3339decSBarry Smith Logically Collective 569baa89ecbSBarry Smith 570baa89ecbSBarry Smith Input Parameter: 571cd47f5d9SBarry Smith . pc - the preconditioner context 572cd47f5d9SBarry Smith 573baa89ecbSBarry Smith Output Parameter: 574baa89ecbSBarry Smith . flg - whether to use absolute values or not 575cd47f5d9SBarry Smith 576cd47f5d9SBarry Smith Level: intermediate 577cd47f5d9SBarry Smith 578*eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()` 579cd47f5d9SBarry Smith @*/ 580d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg) 581d71ae5a4SJacob Faibussowitsch { 582cd47f5d9SBarry Smith PetscFunctionBegin; 5830700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 584cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg)); 5853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 586cd47f5d9SBarry Smith } 587cd47f5d9SBarry Smith 5884b9ad928SBarry Smith /*@ 589*eede4a3fSMark Adams PCJacobiSetRowl1Scale - Set scaling of off-diagonal of operator when computing l1 row norms, eg, 590*eede4a3fSMark Adams Remark 6.1 in "Multigrid Smoothers for Ultraparallel Computing", Baker et al, with 0.5 scaling 591*eede4a3fSMark Adams 592*eede4a3fSMark Adams Logically Collective 593*eede4a3fSMark Adams 594*eede4a3fSMark Adams Input Parameters: 595*eede4a3fSMark Adams + pc - the preconditioner context 596*eede4a3fSMark Adams - scale - scaling 597*eede4a3fSMark Adams 598*eede4a3fSMark Adams Options Database Key: 599*eede4a3fSMark Adams . -pc_jacobi_rowl1_scale <real> - use absolute values 600*eede4a3fSMark Adams 601*eede4a3fSMark Adams Level: intermediate 602*eede4a3fSMark Adams 603*eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetRowl1Scale()` 604*eede4a3fSMark Adams @*/ 605*eede4a3fSMark Adams PetscErrorCode PCJacobiSetRowl1Scale(PC pc, PetscReal scale) 606*eede4a3fSMark Adams { 607*eede4a3fSMark Adams PetscFunctionBegin; 608*eede4a3fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 609*eede4a3fSMark Adams PetscTryMethod(pc, "PCJacobiSetRowl1Scale_C", (PC, PetscReal), (pc, scale)); 610*eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS); 611*eede4a3fSMark Adams } 612*eede4a3fSMark Adams 613*eede4a3fSMark Adams /*@ 614*eede4a3fSMark Adams PCJacobiGetRowl1Scale - Get scaling of off-diagonal elements summed into l1-norm diagonal 615*eede4a3fSMark Adams 616*eede4a3fSMark Adams Logically Collective 617*eede4a3fSMark Adams 618*eede4a3fSMark Adams Input Parameter: 619*eede4a3fSMark Adams . pc - the preconditioner context 620*eede4a3fSMark Adams 621*eede4a3fSMark Adams Output Parameter: 622*eede4a3fSMark Adams . scale - scaling 623*eede4a3fSMark Adams 624*eede4a3fSMark Adams Level: intermediate 625*eede4a3fSMark Adams 626*eede4a3fSMark Adams .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetRowl1Scale()`, `PCJacobiGetType()` 627*eede4a3fSMark Adams @*/ 628*eede4a3fSMark Adams PetscErrorCode PCJacobiGetRowl1Scale(PC pc, PetscReal *scale) 629*eede4a3fSMark Adams { 630*eede4a3fSMark Adams PetscFunctionBegin; 631*eede4a3fSMark Adams PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 632*eede4a3fSMark Adams PetscUseMethod(pc, "PCJacobiGetRowl1Scale_C", (PC, PetscReal *), (pc, scale)); 633*eede4a3fSMark Adams PetscFunctionReturn(PETSC_SUCCESS); 634*eede4a3fSMark Adams } 635*eede4a3fSMark Adams 636*eede4a3fSMark Adams /*@ 637147403d9SBarry Smith PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0 63867ed0f3bSStefano Zampini 639c3339decSBarry Smith Logically Collective 64067ed0f3bSStefano Zampini 64167ed0f3bSStefano Zampini Input Parameters: 64267ed0f3bSStefano Zampini + pc - the preconditioner context 64367ed0f3bSStefano Zampini - flg - the boolean flag 64467ed0f3bSStefano Zampini 64567ed0f3bSStefano Zampini Options Database Key: 646147403d9SBarry Smith . -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal 64767ed0f3bSStefano Zampini 648f1580f4eSBarry Smith Note: 64967ed0f3bSStefano Zampini This takes affect at the next construction of the preconditioner 65067ed0f3bSStefano Zampini 65167ed0f3bSStefano Zampini Level: intermediate 65267ed0f3bSStefano Zampini 653562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()`, `PCJacobiSetUseAbs()` 65467ed0f3bSStefano Zampini @*/ 655d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg) 656d71ae5a4SJacob Faibussowitsch { 65767ed0f3bSStefano Zampini PetscFunctionBegin; 65867ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 659cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg)); 6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66167ed0f3bSStefano Zampini } 66267ed0f3bSStefano Zampini 66367ed0f3bSStefano Zampini /*@ 664f1580f4eSBarry Smith PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner `PCJACOBI` checks for zero diagonal terms 66567ed0f3bSStefano Zampini 666c3339decSBarry Smith Logically Collective 66767ed0f3bSStefano Zampini 66867ed0f3bSStefano Zampini Input Parameter: 66967ed0f3bSStefano Zampini . pc - the preconditioner context 67067ed0f3bSStefano Zampini 67167ed0f3bSStefano Zampini Output Parameter: 67267ed0f3bSStefano Zampini . flg - the boolean flag 67367ed0f3bSStefano Zampini 67467ed0f3bSStefano Zampini Options Database Key: 67567b8a455SSatish Balay . -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1 67667ed0f3bSStefano Zampini 67767ed0f3bSStefano Zampini Level: intermediate 67867ed0f3bSStefano Zampini 679562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()` 68067ed0f3bSStefano Zampini @*/ 681d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg) 682d71ae5a4SJacob Faibussowitsch { 68367ed0f3bSStefano Zampini PetscFunctionBegin; 68467ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 685cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg)); 6863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68767ed0f3bSStefano Zampini } 68867ed0f3bSStefano Zampini 68967ed0f3bSStefano Zampini /*@ 690baa89ecbSBarry Smith PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row, 691baa89ecbSBarry Smith of the sum of rows entries for the diagonal preconditioner 6924b9ad928SBarry Smith 693c3339decSBarry Smith Logically Collective 6944b9ad928SBarry Smith 6954b9ad928SBarry Smith Input Parameters: 696baa89ecbSBarry Smith + pc - the preconditioner context 697*eede4a3fSMark Adams - type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWL1`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM` 6984b9ad928SBarry Smith 6994b9ad928SBarry Smith Options Database Key: 700*eede4a3fSMark Adams . -pc_jacobi_type <diagonal,rowl1,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi 7014b9ad928SBarry Smith 7024b9ad928SBarry Smith Level: intermediate 7034b9ad928SBarry Smith 704feefa0e1SJacob Faibussowitsch Developer Notes: 705f1580f4eSBarry Smith Why is there a separate function for using the absolute value? 706f1580f4eSBarry Smith 707562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()` 7084b9ad928SBarry Smith @*/ 709d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type) 710d71ae5a4SJacob Faibussowitsch { 7114b9ad928SBarry Smith PetscFunctionBegin; 7120700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 713cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type)); 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7154b9ad928SBarry Smith } 7164b9ad928SBarry Smith 71786697f06SMatthew Knepley /*@ 718baa89ecbSBarry Smith PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner 71986697f06SMatthew Knepley 720c3339decSBarry Smith Not Collective 72186697f06SMatthew Knepley 722baa89ecbSBarry Smith Input Parameter: 72386697f06SMatthew Knepley . pc - the preconditioner context 72486697f06SMatthew Knepley 725baa89ecbSBarry Smith Output Parameter: 726*eede4a3fSMark Adams . type - `PC_JACOBI_DIAGONAL`, `PC_JACOBI_ROWL1`, `PC_JACOBI_ROWMAX`, `PC_JACOBI_ROWSUM` 72786697f06SMatthew Knepley 72886697f06SMatthew Knepley Level: intermediate 72986697f06SMatthew Knepley 730562efe2eSBarry Smith .seealso: [](ch_ksp), `PCJACOBI`, `PCJacobiaUseAbs()`, `PCJacobiSetType()` 73186697f06SMatthew Knepley @*/ 732d71ae5a4SJacob Faibussowitsch PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type) 733d71ae5a4SJacob Faibussowitsch { 73486697f06SMatthew Knepley PetscFunctionBegin; 7350700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 736cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type)); 7373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73886697f06SMatthew Knepley } 739