14b9ad928SBarry Smith 24b9ad928SBarry Smith /* -------------------------------------------------------------------- 34b9ad928SBarry Smith 45a46d3faSBarry Smith This file implements a Jacobi preconditioner in PETSc as part of PC. 55a46d3faSBarry Smith You can use this as a starting point for implementing your own 65a46d3faSBarry Smith preconditioner that is not provided with PETSc. (You might also consider 75a46d3faSBarry Smith just using PCSHELL) 84b9ad928SBarry Smith 94b9ad928SBarry Smith The following basic routines are required for each preconditioner. 104b9ad928SBarry Smith PCCreate_XXX() - Creates a preconditioner context 114b9ad928SBarry Smith PCSetFromOptions_XXX() - Sets runtime options 124b9ad928SBarry Smith PCApply_XXX() - Applies the preconditioner 134b9ad928SBarry Smith PCDestroy_XXX() - Destroys the preconditioner context 144b9ad928SBarry Smith where the suffix "_XXX" denotes a particular implementation, in 154b9ad928SBarry Smith this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi). 164b9ad928SBarry Smith These routines are actually called via the common user interface 174b9ad928SBarry Smith routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 184b9ad928SBarry Smith so the application code interface remains identical for all 194b9ad928SBarry Smith preconditioners. 204b9ad928SBarry Smith 214b9ad928SBarry Smith Another key routine is: 224b9ad928SBarry Smith PCSetUp_XXX() - Prepares for the use of a preconditioner 234b9ad928SBarry Smith by setting data structures and options. The interface routine PCSetUp() 244b9ad928SBarry Smith is not usually called directly by the user, but instead is called by 254b9ad928SBarry Smith PCApply() if necessary. 264b9ad928SBarry Smith 274b9ad928SBarry Smith Additional basic routines are: 284b9ad928SBarry Smith PCView_XXX() - Prints details of runtime options that 294b9ad928SBarry Smith have actually been used. 304b9ad928SBarry Smith These are called by application codes via the interface routines 314b9ad928SBarry Smith PCView(). 324b9ad928SBarry Smith 334b9ad928SBarry Smith The various types of solvers (preconditioners, Krylov subspace methods, 344b9ad928SBarry Smith nonlinear solvers, timesteppers) are all organized similarly, so the 354b9ad928SBarry Smith above description applies to these categories also. One exception is 364b9ad928SBarry Smith that the analogues of PCApply() for these components are KSPSolve(), 374b9ad928SBarry Smith SNESSolve(), and TSSolve(). 384b9ad928SBarry Smith 394b9ad928SBarry Smith Additional optional functionality unique to preconditioners is left and 404b9ad928SBarry Smith right symmetric preconditioner application via PCApplySymmetricLeft() 414b9ad928SBarry Smith and PCApplySymmetricRight(). The Jacobi implementation is 424b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi(). 434b9ad928SBarry Smith 444b9ad928SBarry Smith -------------------------------------------------------------------- */ 454b9ad928SBarry Smith 464b9ad928SBarry Smith /* 474b9ad928SBarry Smith Include files needed for the Jacobi preconditioner: 484b9ad928SBarry Smith pcimpl.h - private include file intended for use by all preconditioners 494b9ad928SBarry Smith */ 504b9ad928SBarry Smith 51af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 524b9ad928SBarry Smith 530a545947SLisandro Dalcin const char *const PCJacobiTypes[] = {"DIAGONAL", "ROWMAX", "ROWSUM", "PCJacobiType", "PC_JACOBI_", NULL}; 54baa89ecbSBarry Smith 554b9ad928SBarry Smith /* 564b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 574b9ad928SBarry Smith */ 584b9ad928SBarry Smith typedef struct { 592fa5cd67SKarl Rupp Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 604b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 614b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 624b9ad928SBarry Smith only for symmetric preconditioner application) */ 63baa89ecbSBarry Smith PetscBool userowmax; /* set with PCJacobiSetType() */ 64ace3abfcSBarry Smith PetscBool userowsum; 65ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */ 6667ed0f3bSStefano Zampini PetscBool fixdiag; /* fix zero diagonal terms */ 674b9ad928SBarry Smith } PC_Jacobi; 684b9ad928SBarry Smith 699371c9d4SSatish Balay static PetscErrorCode PCJacobiSetType_Jacobi(PC pc, PCJacobiType type) { 70baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data; 714b9ad928SBarry Smith 724b9ad928SBarry Smith PetscFunctionBegin; 73baa89ecbSBarry Smith j->userowmax = PETSC_FALSE; 74baa89ecbSBarry Smith j->userowsum = PETSC_FALSE; 75baa89ecbSBarry Smith if (type == PC_JACOBI_ROWMAX) { 764b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 77baa89ecbSBarry Smith } else if (type == PC_JACOBI_ROWSUM) { 78baa89ecbSBarry Smith j->userowsum = PETSC_TRUE; 79baa89ecbSBarry Smith } 804b9ad928SBarry Smith PetscFunctionReturn(0); 814b9ad928SBarry Smith } 824b9ad928SBarry Smith 839371c9d4SSatish Balay static PetscErrorCode PCJacobiGetType_Jacobi(PC pc, PCJacobiType *type) { 84baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data; 8586697f06SMatthew Knepley 8686697f06SMatthew Knepley PetscFunctionBegin; 87baa89ecbSBarry Smith if (j->userowmax) { 88baa89ecbSBarry Smith *type = PC_JACOBI_ROWMAX; 89baa89ecbSBarry Smith } else if (j->userowsum) { 90baa89ecbSBarry Smith *type = PC_JACOBI_ROWSUM; 91baa89ecbSBarry Smith } else { 92baa89ecbSBarry Smith *type = PC_JACOBI_DIAGONAL; 93baa89ecbSBarry Smith } 9486697f06SMatthew Knepley PetscFunctionReturn(0); 9586697f06SMatthew Knepley } 9686697f06SMatthew Knepley 979371c9d4SSatish Balay static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc, PetscBool flg) { 98baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data; 99cd47f5d9SBarry Smith 100cd47f5d9SBarry Smith PetscFunctionBegin; 101baa89ecbSBarry Smith j->useabs = flg; 102baa89ecbSBarry Smith PetscFunctionReturn(0); 103baa89ecbSBarry Smith } 104baa89ecbSBarry Smith 1059371c9d4SSatish Balay static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc, PetscBool *flg) { 106baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi *)pc->data; 107baa89ecbSBarry Smith 108baa89ecbSBarry Smith PetscFunctionBegin; 109baa89ecbSBarry Smith *flg = j->useabs; 110cd47f5d9SBarry Smith PetscFunctionReturn(0); 111cd47f5d9SBarry Smith } 112cd47f5d9SBarry Smith 1139371c9d4SSatish Balay static PetscErrorCode PCJacobiSetFixDiagonal_Jacobi(PC pc, PetscBool flg) { 11467ed0f3bSStefano Zampini PC_Jacobi *j = (PC_Jacobi *)pc->data; 11567ed0f3bSStefano Zampini 11667ed0f3bSStefano Zampini PetscFunctionBegin; 11767ed0f3bSStefano Zampini j->fixdiag = flg; 11867ed0f3bSStefano Zampini PetscFunctionReturn(0); 11967ed0f3bSStefano Zampini } 12067ed0f3bSStefano Zampini 1219371c9d4SSatish Balay static PetscErrorCode PCJacobiGetFixDiagonal_Jacobi(PC pc, PetscBool *flg) { 12267ed0f3bSStefano Zampini PC_Jacobi *j = (PC_Jacobi *)pc->data; 12367ed0f3bSStefano Zampini 12467ed0f3bSStefano Zampini PetscFunctionBegin; 12567ed0f3bSStefano Zampini *flg = j->fixdiag; 12667ed0f3bSStefano Zampini PetscFunctionReturn(0); 12767ed0f3bSStefano Zampini } 12867ed0f3bSStefano Zampini 1294b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1304b9ad928SBarry Smith /* 1314b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1324b9ad928SBarry Smith by setting data structures and options. 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith Input Parameter: 1354b9ad928SBarry Smith . pc - the preconditioner context 1364b9ad928SBarry Smith 1374b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith Notes: 1404b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1414b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1424b9ad928SBarry Smith */ 1439371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi(PC pc) { 1444b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 1454b9ad928SBarry Smith Vec diag, diagsqrt; 146cd47f5d9SBarry Smith PetscInt n, i; 1474b9ad928SBarry Smith PetscScalar *x; 148ace3abfcSBarry Smith PetscBool zeroflag = PETSC_FALSE; 1494b9ad928SBarry Smith 1504b9ad928SBarry Smith PetscFunctionBegin; 1514b9ad928SBarry Smith /* 1524b9ad928SBarry Smith For most preconditioners the code would begin here something like 1534b9ad928SBarry Smith 1544b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 1559566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->mat,&jac->diag)); 1563bb1ff40SBarry Smith PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag); 1574b9ad928SBarry Smith } 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1604b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1614b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1624b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1634b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1644b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1654b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1664b9ad928SBarry Smith */ 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith /* 1694b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1704b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1714b9ad928SBarry Smith */ 1724b9ad928SBarry Smith diag = jac->diag; 1734b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith if (diag) { 176b94d7dedSBarry Smith PetscBool isset, isspd; 17767ed0f3bSStefano Zampini 1784b9ad928SBarry Smith if (jac->userowmax) { 1799566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(pc->pmat, diag, NULL)); 18086697f06SMatthew Knepley } else if (jac->userowsum) { 1819566063dSJacob Faibussowitsch PetscCall(MatGetRowSum(pc->pmat, diag)); 1824b9ad928SBarry Smith } else { 1839566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diag)); 1844b9ad928SBarry Smith } 1859566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 1861baa6e33SBarry Smith if (jac->useabs) PetscCall(VecAbs(diag)); 187b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(pc->pmat, &isset, &isspd)); 188b94d7dedSBarry Smith if (jac->fixdiag && (!isset || !isspd)) { 1899566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diag, &n)); 1909566063dSJacob Faibussowitsch PetscCall(VecGetArray(diag, &x)); 1914b9ad928SBarry Smith for (i = 0; i < n; i++) { 1924b9ad928SBarry Smith if (x[i] == 0.0) { 1934b9ad928SBarry Smith x[i] = 1.0; 194cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith } 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diag, &x)); 1984b9ad928SBarry Smith } 19967ed0f3bSStefano Zampini } 2004b9ad928SBarry Smith if (diagsqrt) { 2014b9ad928SBarry Smith if (jac->userowmax) { 2029566063dSJacob Faibussowitsch PetscCall(MatGetRowMaxAbs(pc->pmat, diagsqrt, NULL)); 20386697f06SMatthew Knepley } else if (jac->userowsum) { 2049566063dSJacob Faibussowitsch PetscCall(MatGetRowSum(pc->pmat, diagsqrt)); 2054b9ad928SBarry Smith } else { 2069566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, diagsqrt)); 2074b9ad928SBarry Smith } 2089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(diagsqrt, &n)); 2099566063dSJacob Faibussowitsch PetscCall(VecGetArray(diagsqrt, &x)); 2104b9ad928SBarry Smith for (i = 0; i < n; i++) { 2118f1a2a5eSBarry Smith if (x[i] != 0.0) x[i] = 1.0 / PetscSqrtReal(PetscAbsScalar(x[i])); 2124b9ad928SBarry Smith else { 2134b9ad928SBarry Smith x[i] = 1.0; 214cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 2154b9ad928SBarry Smith } 2164b9ad928SBarry Smith } 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(diagsqrt, &x)); 2184b9ad928SBarry Smith } 219*48a46eb9SPierre Jolivet if (zeroflag) PetscCall(PetscInfo(pc, "Zero detected in diagonal of matrix, using 1 at those locations\n")); 2204b9ad928SBarry Smith PetscFunctionReturn(0); 2214b9ad928SBarry Smith } 2224b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2234b9ad928SBarry Smith /* 2244b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2254b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2264b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2274b9ad928SBarry Smith 2284b9ad928SBarry Smith Input Parameter: 2294b9ad928SBarry Smith . pc - the preconditioner context 2304b9ad928SBarry Smith */ 2319371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) { 2324b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diagsqrt, NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->diagsqrt)); 2379566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc)); 2384b9ad928SBarry Smith PetscFunctionReturn(0); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2414b9ad928SBarry Smith /* 2424b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2434b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2444b9ad928SBarry Smith right application of the Jacobi preconditioner. 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith Input Parameter: 2474b9ad928SBarry Smith . pc - the preconditioner context 2484b9ad928SBarry Smith */ 2499371c9d4SSatish Balay static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) { 2504b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 2514b9ad928SBarry Smith 2524b9ad928SBarry Smith PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &jac->diag, NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->diag)); 2559566063dSJacob Faibussowitsch PetscCall(PCSetUp_Jacobi(pc)); 2564b9ad928SBarry Smith PetscFunctionReturn(0); 2574b9ad928SBarry Smith } 2584b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2594b9ad928SBarry Smith /* 2604b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith Input Parameters: 2634b9ad928SBarry Smith . pc - the preconditioner context 2644b9ad928SBarry Smith . x - input vector 2654b9ad928SBarry Smith 2664b9ad928SBarry Smith Output Parameter: 2674b9ad928SBarry Smith . y - output vector 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith Application Interface Routine: PCApply() 2704b9ad928SBarry Smith */ 2719371c9d4SSatish Balay static PetscErrorCode PCApply_Jacobi(PC pc, Vec x, Vec y) { 2724b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith PetscFunctionBegin; 275*48a46eb9SPierre Jolivet if (!jac->diag) PetscCall(PCSetUp_Jacobi_NonSymmetric(pc)); 2769566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diag)); 2774b9ad928SBarry Smith PetscFunctionReturn(0); 2784b9ad928SBarry Smith } 2794b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2804b9ad928SBarry Smith /* 2814b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2824b9ad928SBarry Smith symmetric preconditioner to a vector. 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith Input Parameters: 2854b9ad928SBarry Smith . pc - the preconditioner context 2864b9ad928SBarry Smith . x - input vector 2874b9ad928SBarry Smith 2884b9ad928SBarry Smith Output Parameter: 2894b9ad928SBarry Smith . y - output vector 2904b9ad928SBarry Smith 2914b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2924b9ad928SBarry Smith */ 2939371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc, Vec x, Vec y) { 2944b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith PetscFunctionBegin; 297*48a46eb9SPierre Jolivet if (!jac->diagsqrt) PetscCall(PCSetUp_Jacobi_Symmetric(pc)); 2989566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, jac->diagsqrt)); 2994b9ad928SBarry Smith PetscFunctionReturn(0); 3004b9ad928SBarry Smith } 30167ed0f3bSStefano Zampini 3024b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3039371c9d4SSatish Balay static PetscErrorCode PCReset_Jacobi(PC pc) { 304a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 305a06653b4SBarry Smith 306a06653b4SBarry Smith PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diag)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->diagsqrt)); 309a06653b4SBarry Smith PetscFunctionReturn(0); 310a06653b4SBarry Smith } 311a06653b4SBarry Smith 3124b9ad928SBarry Smith /* 3134b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3144b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3154b9ad928SBarry Smith 3164b9ad928SBarry Smith Input Parameter: 3174b9ad928SBarry Smith . pc - the preconditioner context 3184b9ad928SBarry Smith 3194b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3204b9ad928SBarry Smith */ 3219371c9d4SSatish Balay static PetscErrorCode PCDestroy_Jacobi(PC pc) { 3224b9ad928SBarry Smith PetscFunctionBegin; 3239566063dSJacob Faibussowitsch PetscCall(PCReset_Jacobi(pc)); 3249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", NULL)); 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", NULL)); 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", NULL)); 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", NULL)); 3289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", NULL)); 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", NULL)); 3304b9ad928SBarry Smith 3314b9ad928SBarry Smith /* 3324b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3334b9ad928SBarry Smith */ 3349566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 3354b9ad928SBarry Smith PetscFunctionReturn(0); 3364b9ad928SBarry Smith } 3374b9ad928SBarry Smith 3389371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Jacobi(PC pc, PetscOptionItems *PetscOptionsObject) { 3394b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi *)pc->data; 340baa89ecbSBarry Smith PetscBool flg; 341baa89ecbSBarry Smith PCJacobiType deflt, type; 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith PetscFunctionBegin; 3449566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &deflt)); 345d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Jacobi options"); 3469566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_jacobi_type", "How to construct diagonal matrix", "PCJacobiSetType", PCJacobiTypes, (PetscEnum)deflt, (PetscEnum *)&type, &flg)); 3471baa6e33SBarry Smith if (flg) PetscCall(PCJacobiSetType(pc, type)); 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_abs", "Use absolute values of diagonal entries", "PCJacobiSetUseAbs", jac->useabs, &jac->useabs, NULL)); 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_jacobi_fixdiagonal", "Fix null terms on diagonal", "PCJacobiSetFixDiagonal", jac->fixdiag, &jac->fixdiag, NULL)); 350d0609cedSBarry Smith PetscOptionsHeadEnd(); 3514b9ad928SBarry Smith PetscFunctionReturn(0); 3524b9ad928SBarry Smith } 3534b9ad928SBarry Smith 3549371c9d4SSatish Balay static PetscErrorCode PCView_Jacobi(PC pc, PetscViewer viewer) { 355569e28a7SMatthew G. Knepley PC_Jacobi *jac = (PC_Jacobi *)pc->data; 356569e28a7SMatthew G. Knepley PetscBool iascii; 357569e28a7SMatthew G. Knepley 358569e28a7SMatthew G. Knepley PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 360569e28a7SMatthew G. Knepley if (iascii) { 361569e28a7SMatthew G. Knepley PCJacobiType type; 36267ed0f3bSStefano Zampini PetscBool useAbs, fixdiag; 363569e28a7SMatthew G. Knepley PetscViewerFormat format; 364569e28a7SMatthew G. Knepley 3659566063dSJacob Faibussowitsch PetscCall(PCJacobiGetType(pc, &type)); 3669566063dSJacob Faibussowitsch PetscCall(PCJacobiGetUseAbs(pc, &useAbs)); 3679566063dSJacob Faibussowitsch PetscCall(PCJacobiGetFixDiagonal(pc, &fixdiag)); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type %s%s%s\n", PCJacobiTypes[type], useAbs ? ", using absolute value of entries" : "", !fixdiag ? ", not checking null diagonal entries" : "")); 3699566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 370*48a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(VecView(jac->diag, viewer)); 371569e28a7SMatthew G. Knepley } 372569e28a7SMatthew G. Knepley PetscFunctionReturn(0); 373569e28a7SMatthew G. Knepley } 374569e28a7SMatthew G. Knepley 3754b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3764b9ad928SBarry Smith /* 3774b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3784b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3794b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3804b9ad928SBarry Smith 3814b9ad928SBarry Smith Input Parameter: 3824b9ad928SBarry Smith . pc - the preconditioner context 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith Application Interface Routine: PCCreate() 3854b9ad928SBarry Smith */ 3864b9ad928SBarry Smith 3874b9ad928SBarry Smith /*MC 3885a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith Options Database Key: 391967c93d3SBarry Smith + -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner 39267ed0f3bSStefano Zampini . -pc_jacobi_abs - use the absolute value of the diagonal entry 393147403d9SBarry Smith - -pc_jacobi_fixdiag - fix for zero diagonal terms by placing 1.0 in those locations 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith Level: beginner 3964b9ad928SBarry Smith 39795452b02SPatrick Sanan Notes: 39895452b02SPatrick Sanan By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 3994b9ad928SBarry Smith can scale each side of the matrix by the square root of the diagonal entries. 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 4024b9ad928SBarry Smith 403468ae2e8SBarry Smith See PCPBJACOBI for fixed-size point block, PCVPBJACOBI for variable-sized point block, and PCBJACOBI for large size blocks 404422a814eSBarry Smith 405db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 406db781477SPatrick Sanan `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, 407db781477SPatrick Sanan `PCJacobiSetFixDiagonal()`, `PCJacobiGetFixDiagonal()` 408db781477SPatrick Sanan `PCJacobiSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetUseAbs()`, `PCPBJACOBI`, `PCBJACOBI`, `PCVPBJACOBI` 4094b9ad928SBarry Smith M*/ 4104b9ad928SBarry Smith 4119371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) { 4124b9ad928SBarry Smith PC_Jacobi *jac; 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith PetscFunctionBegin; 4154b9ad928SBarry Smith /* 4164b9ad928SBarry Smith Creates the private data structure for this preconditioner and 4174b9ad928SBarry Smith attach it to the PC object. 4184b9ad928SBarry Smith */ 4199566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &jac)); 4204b9ad928SBarry Smith pc->data = (void *)jac; 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith /* 4234b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4244b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4254b9ad928SBarry Smith */ 4260a545947SLisandro Dalcin jac->diag = NULL; 4270a545947SLisandro Dalcin jac->diagsqrt = NULL; 4284b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 42986697f06SMatthew Knepley jac->userowsum = PETSC_FALSE; 430cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 43167ed0f3bSStefano Zampini jac->fixdiag = PETSC_TRUE; 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith /* 4344b9ad928SBarry Smith Set the pointers for the functions that are provided above. 4354b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4364b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4374b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4384b9ad928SBarry Smith not needed. 4394b9ad928SBarry Smith */ 4404b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4414b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4424b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 443a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi; 4444b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4454b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 446569e28a7SMatthew G. Knepley pc->ops->view = PCView_Jacobi; 4470a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 4484b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4494b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 4502fa5cd67SKarl Rupp 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetType_C", PCJacobiSetType_Jacobi)); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetType_C", PCJacobiGetType_Jacobi)); 4539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetUseAbs_C", PCJacobiSetUseAbs_Jacobi)); 4549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetUseAbs_C", PCJacobiGetUseAbs_Jacobi)); 4559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiSetFixDiagonal_C", PCJacobiSetFixDiagonal_Jacobi)); 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCJacobiGetFixDiagonal_C", PCJacobiGetFixDiagonal_Jacobi)); 4574b9ad928SBarry Smith PetscFunctionReturn(0); 4584b9ad928SBarry Smith } 459cd47f5d9SBarry Smith 460cd47f5d9SBarry Smith /*@ 461cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 46267ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner 463cd47f5d9SBarry Smith 464ad4df100SBarry Smith Logically Collective on PC 465cd47f5d9SBarry Smith 466cd47f5d9SBarry Smith Input Parameters: 467baa89ecbSBarry Smith + pc - the preconditioner context 468baa89ecbSBarry Smith - flg - whether to use absolute values or not 469baa89ecbSBarry Smith 470baa89ecbSBarry Smith Options Database Key: 471147403d9SBarry Smith . -pc_jacobi_abs <bool> - use absolute values 472baa89ecbSBarry Smith 47395452b02SPatrick Sanan Notes: 47495452b02SPatrick Sanan This takes affect at the next construction of the preconditioner 475baa89ecbSBarry Smith 476baa89ecbSBarry Smith Level: intermediate 477baa89ecbSBarry Smith 478db781477SPatrick Sanan .seealso: `PCJacobiaSetType()`, `PCJacobiGetUseAbs()` 479baa89ecbSBarry Smith 480baa89ecbSBarry Smith @*/ 4819371c9d4SSatish Balay PetscErrorCode PCJacobiSetUseAbs(PC pc, PetscBool flg) { 482baa89ecbSBarry Smith PetscFunctionBegin; 483baa89ecbSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 484cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetUseAbs_C", (PC, PetscBool), (pc, flg)); 485baa89ecbSBarry Smith PetscFunctionReturn(0); 486baa89ecbSBarry Smith } 487baa89ecbSBarry Smith 488baa89ecbSBarry Smith /*@ 489baa89ecbSBarry Smith PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the 49067ed0f3bSStefano Zampini absolute values of the diagonal divisors in the preconditioner 491baa89ecbSBarry Smith 492baa89ecbSBarry Smith Logically Collective on PC 493baa89ecbSBarry Smith 494baa89ecbSBarry Smith Input Parameter: 495cd47f5d9SBarry Smith . pc - the preconditioner context 496cd47f5d9SBarry Smith 497baa89ecbSBarry Smith Output Parameter: 498baa89ecbSBarry Smith . flg - whether to use absolute values or not 499cd47f5d9SBarry Smith 500cd47f5d9SBarry Smith Level: intermediate 501cd47f5d9SBarry Smith 502db781477SPatrick Sanan .seealso: `PCJacobiaSetType()`, `PCJacobiSetUseAbs()`, `PCJacobiGetType()` 503cd47f5d9SBarry Smith 504cd47f5d9SBarry Smith @*/ 5059371c9d4SSatish Balay PetscErrorCode PCJacobiGetUseAbs(PC pc, PetscBool *flg) { 506cd47f5d9SBarry Smith PetscFunctionBegin; 5070700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 508cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetUseAbs_C", (PC, PetscBool *), (pc, flg)); 509cd47f5d9SBarry Smith PetscFunctionReturn(0); 510cd47f5d9SBarry Smith } 511cd47f5d9SBarry Smith 5124b9ad928SBarry Smith /*@ 513147403d9SBarry Smith PCJacobiSetFixDiagonal - Check for zero values on the diagonal and replace them with 1.0 51467ed0f3bSStefano Zampini 51567ed0f3bSStefano Zampini Logically Collective on PC 51667ed0f3bSStefano Zampini 51767ed0f3bSStefano Zampini Input Parameters: 51867ed0f3bSStefano Zampini + pc - the preconditioner context 51967ed0f3bSStefano Zampini - flg - the boolean flag 52067ed0f3bSStefano Zampini 52167ed0f3bSStefano Zampini Options Database Key: 522147403d9SBarry Smith . -pc_jacobi_fixdiagonal <bool> - check for zero values on the diagonal 52367ed0f3bSStefano Zampini 52467ed0f3bSStefano Zampini Notes: 52567ed0f3bSStefano Zampini This takes affect at the next construction of the preconditioner 52667ed0f3bSStefano Zampini 52767ed0f3bSStefano Zampini Level: intermediate 52867ed0f3bSStefano Zampini 529db781477SPatrick Sanan .seealso: `PCJacobiSetType()`, `PCJacobiGetFixDiagonal()` 53067ed0f3bSStefano Zampini 53167ed0f3bSStefano Zampini @*/ 5329371c9d4SSatish Balay PetscErrorCode PCJacobiSetFixDiagonal(PC pc, PetscBool flg) { 53367ed0f3bSStefano Zampini PetscFunctionBegin; 53467ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 535cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetFixDiagonal_C", (PC, PetscBool), (pc, flg)); 53667ed0f3bSStefano Zampini PetscFunctionReturn(0); 53767ed0f3bSStefano Zampini } 53867ed0f3bSStefano Zampini 53967ed0f3bSStefano Zampini /*@ 54067ed0f3bSStefano Zampini PCJacobiGetFixDiagonal - Determines if the Jacobi preconditioner checks for zero diagonal terms 54167ed0f3bSStefano Zampini 54267ed0f3bSStefano Zampini Logically Collective on PC 54367ed0f3bSStefano Zampini 54467ed0f3bSStefano Zampini Input Parameter: 54567ed0f3bSStefano Zampini . pc - the preconditioner context 54667ed0f3bSStefano Zampini 54767ed0f3bSStefano Zampini Output Parameter: 54867ed0f3bSStefano Zampini . flg - the boolean flag 54967ed0f3bSStefano Zampini 55067ed0f3bSStefano Zampini Options Database Key: 55167b8a455SSatish Balay . -pc_jacobi_fixdiagonal <bool> - Fix 0 terms on diagonal by using 1 55267ed0f3bSStefano Zampini 55367ed0f3bSStefano Zampini Level: intermediate 55467ed0f3bSStefano Zampini 555db781477SPatrick Sanan .seealso: `PCJacobiSetType()`, `PCJacobiSetFixDiagonal()` 55667ed0f3bSStefano Zampini 55767ed0f3bSStefano Zampini @*/ 5589371c9d4SSatish Balay PetscErrorCode PCJacobiGetFixDiagonal(PC pc, PetscBool *flg) { 55967ed0f3bSStefano Zampini PetscFunctionBegin; 56067ed0f3bSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 561cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetFixDiagonal_C", (PC, PetscBool *), (pc, flg)); 56267ed0f3bSStefano Zampini PetscFunctionReturn(0); 56367ed0f3bSStefano Zampini } 56467ed0f3bSStefano Zampini 56567ed0f3bSStefano Zampini /*@ 566baa89ecbSBarry Smith PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row, 567baa89ecbSBarry Smith of the sum of rows entries for the diagonal preconditioner 5684b9ad928SBarry Smith 569ad4df100SBarry Smith Logically Collective on PC 5704b9ad928SBarry Smith 5714b9ad928SBarry Smith Input Parameters: 572baa89ecbSBarry Smith + pc - the preconditioner context 573baa89ecbSBarry Smith - type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 5744b9ad928SBarry Smith 5754b9ad928SBarry Smith Options Database Key: 576147403d9SBarry Smith . -pc_jacobi_type <diagonal,rowmax,rowsum> - the type of diagonal matrix to use for Jacobi 5774b9ad928SBarry Smith 5784b9ad928SBarry Smith Level: intermediate 5794b9ad928SBarry Smith 580db781477SPatrick Sanan .seealso: `PCJacobiaUseAbs()`, `PCJacobiGetType()` 5814b9ad928SBarry Smith @*/ 5829371c9d4SSatish Balay PetscErrorCode PCJacobiSetType(PC pc, PCJacobiType type) { 5834b9ad928SBarry Smith PetscFunctionBegin; 5840700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 585cac4c232SBarry Smith PetscTryMethod(pc, "PCJacobiSetType_C", (PC, PCJacobiType), (pc, type)); 5864b9ad928SBarry Smith PetscFunctionReturn(0); 5874b9ad928SBarry Smith } 5884b9ad928SBarry Smith 58986697f06SMatthew Knepley /*@ 590baa89ecbSBarry Smith PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner 59186697f06SMatthew Knepley 592baa89ecbSBarry Smith Not Collective on PC 59386697f06SMatthew Knepley 594baa89ecbSBarry Smith Input Parameter: 59586697f06SMatthew Knepley . pc - the preconditioner context 59686697f06SMatthew Knepley 597baa89ecbSBarry Smith Output Parameter: 598baa89ecbSBarry Smith . type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 59986697f06SMatthew Knepley 60086697f06SMatthew Knepley Level: intermediate 60186697f06SMatthew Knepley 602db781477SPatrick Sanan .seealso: `PCJacobiaUseAbs()`, `PCJacobiSetType()` 60386697f06SMatthew Knepley @*/ 6049371c9d4SSatish Balay PetscErrorCode PCJacobiGetType(PC pc, PCJacobiType *type) { 60586697f06SMatthew Knepley PetscFunctionBegin; 6060700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 607cac4c232SBarry Smith PetscUseMethod(pc, "PCJacobiGetType_C", (PC, PCJacobiType *), (pc, type)); 60886697f06SMatthew Knepley PetscFunctionReturn(0); 60986697f06SMatthew Knepley } 610