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 51b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 524b9ad928SBarry Smith 53*baa89ecbSBarry Smith const char *const PCJacobiTypes[] = {"DIAGONAL","ROWMAX","ROWSUM","PCJacobiType","PC_JACOBI_",0}; 54*baa89ecbSBarry 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) */ 63*baa89ecbSBarry Smith PetscBool userowmax; /* set with PCJacobiSetType() */ 64ace3abfcSBarry Smith PetscBool userowsum; 65ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */ 664b9ad928SBarry Smith } PC_Jacobi; 674b9ad928SBarry Smith 684b9ad928SBarry Smith #undef __FUNCT__ 69*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiSetType_Jacobi" 70*baa89ecbSBarry Smith static PetscErrorCode PCJacobiSetType_Jacobi(PC pc,PCJacobiType type) 714b9ad928SBarry Smith { 72*baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi*)pc->data; 734b9ad928SBarry Smith 744b9ad928SBarry Smith PetscFunctionBegin; 75*baa89ecbSBarry Smith j->userowmax = PETSC_FALSE; 76*baa89ecbSBarry Smith j->userowsum = PETSC_FALSE; 77*baa89ecbSBarry Smith if (type == PC_JACOBI_ROWMAX) { 784b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 79*baa89ecbSBarry Smith } else if (type == PC_JACOBI_ROWSUM) { 80*baa89ecbSBarry Smith j->userowsum = PETSC_TRUE; 81*baa89ecbSBarry Smith } 824b9ad928SBarry Smith PetscFunctionReturn(0); 834b9ad928SBarry Smith } 844b9ad928SBarry Smith 85cd47f5d9SBarry Smith #undef __FUNCT__ 86*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiGetType_Jacobi" 87*baa89ecbSBarry Smith static PetscErrorCode PCJacobiGetType_Jacobi(PC pc,PCJacobiType *type) 8886697f06SMatthew Knepley { 89*baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi*)pc->data; 9086697f06SMatthew Knepley 9186697f06SMatthew Knepley PetscFunctionBegin; 92*baa89ecbSBarry Smith if (j->userowmax) { 93*baa89ecbSBarry Smith *type = PC_JACOBI_ROWMAX; 94*baa89ecbSBarry Smith } else if (j->userowsum) { 95*baa89ecbSBarry Smith *type = PC_JACOBI_ROWSUM; 96*baa89ecbSBarry Smith } else { 97*baa89ecbSBarry Smith *type = PC_JACOBI_DIAGONAL; 98*baa89ecbSBarry Smith } 9986697f06SMatthew Knepley PetscFunctionReturn(0); 10086697f06SMatthew Knepley } 10186697f06SMatthew Knepley 10286697f06SMatthew Knepley #undef __FUNCT__ 103cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 104*baa89ecbSBarry Smith static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg) 105cd47f5d9SBarry Smith { 106*baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi*)pc->data; 107cd47f5d9SBarry Smith 108cd47f5d9SBarry Smith PetscFunctionBegin; 109*baa89ecbSBarry Smith j->useabs = flg; 110*baa89ecbSBarry Smith PetscFunctionReturn(0); 111*baa89ecbSBarry Smith } 112*baa89ecbSBarry Smith 113*baa89ecbSBarry Smith #undef __FUNCT__ 114*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiGetUseAbs_Jacobi" 115*baa89ecbSBarry Smith static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool *flg) 116*baa89ecbSBarry Smith { 117*baa89ecbSBarry Smith PC_Jacobi *j = (PC_Jacobi*)pc->data; 118*baa89ecbSBarry Smith 119*baa89ecbSBarry Smith PetscFunctionBegin; 120*baa89ecbSBarry Smith *flg = j->useabs; 121cd47f5d9SBarry Smith PetscFunctionReturn(0); 122cd47f5d9SBarry Smith } 123cd47f5d9SBarry Smith 1244b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1254b9ad928SBarry Smith /* 1264b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1274b9ad928SBarry Smith by setting data structures and options. 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith Input Parameter: 1304b9ad928SBarry Smith . pc - the preconditioner context 1314b9ad928SBarry Smith 1324b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith Notes: 1354b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1364b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1374b9ad928SBarry Smith */ 1384b9ad928SBarry Smith #undef __FUNCT__ 1394b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1406849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1414b9ad928SBarry Smith { 1424b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1434b9ad928SBarry Smith Vec diag,diagsqrt; 144dfbe8321SBarry Smith PetscErrorCode ierr; 145cd47f5d9SBarry Smith PetscInt n,i; 1464b9ad928SBarry Smith PetscScalar *x; 147ace3abfcSBarry Smith PetscBool zeroflag = PETSC_FALSE; 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith PetscFunctionBegin; 1504b9ad928SBarry Smith /* 1514b9ad928SBarry Smith For most preconditioners the code would begin here something like 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 1542a7a6963SBarry Smith ierr = MatCreateVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1553bb1ff40SBarry Smith PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag); 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1594b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1604b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1614b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1624b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1634b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1644b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1654b9ad928SBarry Smith */ 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith /* 1684b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1694b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1704b9ad928SBarry Smith */ 1714b9ad928SBarry Smith diag = jac->diag; 1724b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith if (diag) { 1754b9ad928SBarry Smith if (jac->userowmax) { 1760298fd71SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr); 17786697f06SMatthew Knepley } else if (jac->userowsum) { 17886697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 1794b9ad928SBarry Smith } else { 1804b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1834b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1844b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 185cd47f5d9SBarry Smith if (jac->useabs) { 1862fa5cd67SKarl Rupp for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 187cd47f5d9SBarry Smith } 1884b9ad928SBarry Smith for (i=0; i<n; i++) { 1894b9ad928SBarry Smith if (x[i] == 0.0) { 1904b9ad928SBarry Smith x[i] = 1.0; 191cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1924b9ad928SBarry Smith } 1934b9ad928SBarry Smith } 1944b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith if (diagsqrt) { 1974b9ad928SBarry Smith if (jac->userowmax) { 1980298fd71SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr); 19986697f06SMatthew Knepley } else if (jac->userowsum) { 20086697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 2014b9ad928SBarry Smith } else { 2024b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 2034b9ad928SBarry Smith } 2044b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 2054b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 2064b9ad928SBarry Smith for (i=0; i<n; i++) { 2078f1a2a5eSBarry Smith if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 2084b9ad928SBarry Smith else { 2094b9ad928SBarry Smith x[i] = 1.0; 210cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 2114b9ad928SBarry Smith } 2124b9ad928SBarry Smith } 2134b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 2144b9ad928SBarry Smith } 2154b9ad928SBarry Smith if (zeroflag) { 216ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 2174b9ad928SBarry Smith } 2184b9ad928SBarry Smith PetscFunctionReturn(0); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2214b9ad928SBarry Smith /* 2224b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2234b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2244b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Input Parameter: 2274b9ad928SBarry Smith . pc - the preconditioner context 2284b9ad928SBarry Smith */ 2294b9ad928SBarry Smith #undef __FUNCT__ 2304b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2316849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2324b9ad928SBarry Smith { 233dfbe8321SBarry Smith PetscErrorCode ierr; 2344b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith PetscFunctionBegin; 2372a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 2383bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt);CHKERRQ(ierr); 2394b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2404b9ad928SBarry Smith PetscFunctionReturn(0); 2414b9ad928SBarry Smith } 2424b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2434b9ad928SBarry Smith /* 2444b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2454b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2464b9ad928SBarry Smith right application of the Jacobi preconditioner. 2474b9ad928SBarry Smith 2484b9ad928SBarry Smith Input Parameter: 2494b9ad928SBarry Smith . pc - the preconditioner context 2504b9ad928SBarry Smith */ 2514b9ad928SBarry Smith #undef __FUNCT__ 2524b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 2536849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2544b9ad928SBarry Smith { 255dfbe8321SBarry Smith PetscErrorCode ierr; 2564b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith PetscFunctionBegin; 2592a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 2603bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);CHKERRQ(ierr); 2614b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2624b9ad928SBarry Smith PetscFunctionReturn(0); 2634b9ad928SBarry Smith } 2644b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2654b9ad928SBarry Smith /* 2664b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith Input Parameters: 2694b9ad928SBarry Smith . pc - the preconditioner context 2704b9ad928SBarry Smith . x - input vector 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith Output Parameter: 2734b9ad928SBarry Smith . y - output vector 2744b9ad928SBarry Smith 2754b9ad928SBarry Smith Application Interface Routine: PCApply() 2764b9ad928SBarry Smith */ 2774b9ad928SBarry Smith #undef __FUNCT__ 2784b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2796849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2804b9ad928SBarry Smith { 2814b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 282dfbe8321SBarry Smith PetscErrorCode ierr; 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith PetscFunctionBegin; 2854b9ad928SBarry Smith if (!jac->diag) { 2864b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2874b9ad928SBarry Smith } 2882dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 2894b9ad928SBarry Smith PetscFunctionReturn(0); 2904b9ad928SBarry Smith } 2914b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2924b9ad928SBarry Smith /* 2934b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2944b9ad928SBarry Smith symmetric preconditioner to a vector. 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith Input Parameters: 2974b9ad928SBarry Smith . pc - the preconditioner context 2984b9ad928SBarry Smith . x - input vector 2994b9ad928SBarry Smith 3004b9ad928SBarry Smith Output Parameter: 3014b9ad928SBarry Smith . y - output vector 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 3044b9ad928SBarry Smith */ 3054b9ad928SBarry Smith #undef __FUNCT__ 3064b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 3076849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 3084b9ad928SBarry Smith { 309dfbe8321SBarry Smith PetscErrorCode ierr; 3104b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 3114b9ad928SBarry Smith 3124b9ad928SBarry Smith PetscFunctionBegin; 3134b9ad928SBarry Smith if (!jac->diagsqrt) { 3144b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 3154b9ad928SBarry Smith } 3162dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 3174b9ad928SBarry Smith PetscFunctionReturn(0); 3184b9ad928SBarry Smith } 3194b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 320a06653b4SBarry Smith #undef __FUNCT__ 321a06653b4SBarry Smith #define __FUNCT__ "PCReset_Jacobi" 322a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc) 323a06653b4SBarry Smith { 324a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 325a06653b4SBarry Smith PetscErrorCode ierr; 326a06653b4SBarry Smith 327a06653b4SBarry Smith PetscFunctionBegin; 3286bf464f9SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 3296bf464f9SBarry Smith ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 330a06653b4SBarry Smith PetscFunctionReturn(0); 331a06653b4SBarry Smith } 332a06653b4SBarry Smith 3334b9ad928SBarry Smith /* 3344b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3354b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3364b9ad928SBarry Smith 3374b9ad928SBarry Smith Input Parameter: 3384b9ad928SBarry Smith . pc - the preconditioner context 3394b9ad928SBarry Smith 3404b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3414b9ad928SBarry Smith */ 3424b9ad928SBarry Smith #undef __FUNCT__ 3434b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3446849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3454b9ad928SBarry Smith { 346dfbe8321SBarry Smith PetscErrorCode ierr; 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith PetscFunctionBegin; 349a06653b4SBarry Smith ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith /* 3524b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3534b9ad928SBarry Smith */ 354c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 3554b9ad928SBarry Smith PetscFunctionReturn(0); 3564b9ad928SBarry Smith } 3574b9ad928SBarry Smith 3584b9ad928SBarry Smith #undef __FUNCT__ 3594b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3606849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 3614b9ad928SBarry Smith { 3624b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 363dfbe8321SBarry Smith PetscErrorCode ierr; 364*baa89ecbSBarry Smith PetscBool flg; 365*baa89ecbSBarry Smith PCJacobiType deflt,type; 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith PetscFunctionBegin; 368*baa89ecbSBarry Smith ierr = PCJacobiGetType(pc,&deflt);CHKERRQ(ierr); 3694b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 370*baa89ecbSBarry Smith ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 371*baa89ecbSBarry Smith if (flg) { 372*baa89ecbSBarry Smith ierr = PCJacobiSetType(pc,type);CHKERRQ(ierr); 373*baa89ecbSBarry Smith } 3748afaa268SBarry Smith ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL);CHKERRQ(ierr); 3754b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3764b9ad928SBarry Smith PetscFunctionReturn(0); 3774b9ad928SBarry Smith } 3784b9ad928SBarry Smith 3794b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3804b9ad928SBarry Smith /* 3814b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3824b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3834b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Input Parameter: 3864b9ad928SBarry Smith . pc - the preconditioner context 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith Application Interface Routine: PCCreate() 3894b9ad928SBarry Smith */ 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith /*MC 3925a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith Options Database Key: 395*baa89ecbSBarry Smith + -pc_jacobi_type <diagonal,rowmax,rowsum> 396cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith Level: beginner 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 4014b9ad928SBarry Smith 402b037da10SBarry Smith Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 4034b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 4044b9ad928SBarry Smith 4054b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 408*baa89ecbSBarry Smith PCJacobiSetType(), PCJacobiSetUseAbs(), PCJacobiGetUseAbs() 4094b9ad928SBarry Smith M*/ 4104b9ad928SBarry Smith 4114b9ad928SBarry Smith #undef __FUNCT__ 4124b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 4138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) 4144b9ad928SBarry Smith { 4154b9ad928SBarry Smith PC_Jacobi *jac; 416dfbe8321SBarry Smith PetscErrorCode ierr; 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith PetscFunctionBegin; 4194b9ad928SBarry Smith /* 4204b9ad928SBarry Smith Creates the private data structure for this preconditioner and 4214b9ad928SBarry Smith attach it to the PC object. 4224b9ad928SBarry Smith */ 423b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 4244b9ad928SBarry Smith pc->data = (void*)jac; 4254b9ad928SBarry Smith 4264b9ad928SBarry Smith /* 4274b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4284b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4294b9ad928SBarry Smith */ 4304b9ad928SBarry Smith jac->diag = 0; 4314b9ad928SBarry Smith jac->diagsqrt = 0; 4324b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 43386697f06SMatthew Knepley jac->userowsum = PETSC_FALSE; 434cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith /* 4374b9ad928SBarry Smith Set the pointers for the functions that are provided above. 4384b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4394b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4404b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4414b9ad928SBarry Smith not needed. 4424b9ad928SBarry Smith */ 4434b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4444b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4454b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 446a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi; 4474b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4484b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4494b9ad928SBarry Smith pc->ops->view = 0; 4504b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4514b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4524b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 4532fa5cd67SKarl Rupp 454*baa89ecbSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi);CHKERRQ(ierr); 455*baa89ecbSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi);CHKERRQ(ierr); 456bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 457*baa89ecbSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi);CHKERRQ(ierr); 4584b9ad928SBarry Smith PetscFunctionReturn(0); 4594b9ad928SBarry Smith } 460cd47f5d9SBarry Smith 461cd47f5d9SBarry Smith #undef __FUNCT__ 462cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs" 463cd47f5d9SBarry Smith /*@ 464cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 465*baa89ecbSBarry Smith absolute values of the digonal divisors in the preconditioner 466cd47f5d9SBarry Smith 467ad4df100SBarry Smith Logically Collective on PC 468cd47f5d9SBarry Smith 469cd47f5d9SBarry Smith Input Parameters: 470*baa89ecbSBarry Smith + pc - the preconditioner context 471*baa89ecbSBarry Smith - flg - whether to use absolute values or not 472*baa89ecbSBarry Smith 473*baa89ecbSBarry Smith Options Database Key: 474*baa89ecbSBarry Smith . -pc_jacobi_abs 475*baa89ecbSBarry Smith 476*baa89ecbSBarry Smith Notes: This takes affect at the next construction of the preconditioner 477*baa89ecbSBarry Smith 478*baa89ecbSBarry Smith Level: intermediate 479*baa89ecbSBarry Smith 480*baa89ecbSBarry Smith Concepts: Jacobi preconditioner 481*baa89ecbSBarry Smith 482*baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiGetUseAbs() 483*baa89ecbSBarry Smith 484*baa89ecbSBarry Smith @*/ 485*baa89ecbSBarry Smith PetscErrorCode PCJacobiSetUseAbs(PC pc,PetscBool flg) 486*baa89ecbSBarry Smith { 487*baa89ecbSBarry Smith PetscErrorCode ierr; 488*baa89ecbSBarry Smith 489*baa89ecbSBarry Smith PetscFunctionBegin; 490*baa89ecbSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 491*baa89ecbSBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 492*baa89ecbSBarry Smith PetscFunctionReturn(0); 493*baa89ecbSBarry Smith } 494*baa89ecbSBarry Smith 495*baa89ecbSBarry Smith #undef __FUNCT__ 496*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiGetUseAbs" 497*baa89ecbSBarry Smith /*@ 498*baa89ecbSBarry Smith PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the 499*baa89ecbSBarry Smith absolute values of the digonal divisors in the preconditioner 500*baa89ecbSBarry Smith 501*baa89ecbSBarry Smith Logically Collective on PC 502*baa89ecbSBarry Smith 503*baa89ecbSBarry Smith Input Parameter: 504cd47f5d9SBarry Smith . pc - the preconditioner context 505cd47f5d9SBarry Smith 506*baa89ecbSBarry Smith Output Parameter: 507*baa89ecbSBarry Smith . flg - whether to use absolute values or not 508cd47f5d9SBarry Smith 509cd47f5d9SBarry Smith Options Database Key: 510cd47f5d9SBarry Smith . -pc_jacobi_abs 511cd47f5d9SBarry Smith 512cd47f5d9SBarry Smith Level: intermediate 513cd47f5d9SBarry Smith 514cd47f5d9SBarry Smith Concepts: Jacobi preconditioner 515cd47f5d9SBarry Smith 516*baa89ecbSBarry Smith .seealso: PCJacobiaSetType(), PCJacobiSetUseAbs(), PCJacobiGetType() 517cd47f5d9SBarry Smith 518cd47f5d9SBarry Smith @*/ 519*baa89ecbSBarry Smith PetscErrorCode PCJacobiGetUseAbs(PC pc,PetscBool *flg) 520cd47f5d9SBarry Smith { 5214ac538c5SBarry Smith PetscErrorCode ierr; 522cd47f5d9SBarry Smith 523cd47f5d9SBarry Smith PetscFunctionBegin; 5240700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 525*baa89ecbSBarry Smith ierr = PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr); 526cd47f5d9SBarry Smith PetscFunctionReturn(0); 527cd47f5d9SBarry Smith } 528cd47f5d9SBarry Smith 5294b9ad928SBarry Smith #undef __FUNCT__ 530*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiSetType" 5314b9ad928SBarry Smith /*@ 532*baa89ecbSBarry Smith PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row, 533*baa89ecbSBarry Smith of the sum of rows entries for the diagonal preconditioner 5344b9ad928SBarry Smith 535ad4df100SBarry Smith Logically Collective on PC 5364b9ad928SBarry Smith 5374b9ad928SBarry Smith Input Parameters: 538*baa89ecbSBarry Smith + pc - the preconditioner context 539*baa89ecbSBarry Smith - type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 5404b9ad928SBarry Smith 5414b9ad928SBarry Smith Options Database Key: 542*baa89ecbSBarry Smith . -pc_jacobi_type <diagonal,rowmax,rowsum> 5434b9ad928SBarry Smith 5444b9ad928SBarry Smith Level: intermediate 5454b9ad928SBarry Smith 5464b9ad928SBarry Smith Concepts: Jacobi preconditioner 5474b9ad928SBarry Smith 548*baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiGetType() 5494b9ad928SBarry Smith @*/ 550*baa89ecbSBarry Smith PetscErrorCode PCJacobiSetType(PC pc,PCJacobiType type) 5514b9ad928SBarry Smith { 5524ac538c5SBarry Smith PetscErrorCode ierr; 5534b9ad928SBarry Smith 5544b9ad928SBarry Smith PetscFunctionBegin; 5550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 556*baa89ecbSBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));CHKERRQ(ierr); 5574b9ad928SBarry Smith PetscFunctionReturn(0); 5584b9ad928SBarry Smith } 5594b9ad928SBarry Smith 56086697f06SMatthew Knepley #undef __FUNCT__ 561*baa89ecbSBarry Smith #define __FUNCT__ "PCJacobiGetType" 56286697f06SMatthew Knepley /*@ 563*baa89ecbSBarry Smith PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner 56486697f06SMatthew Knepley 565*baa89ecbSBarry Smith Not Collective on PC 56686697f06SMatthew Knepley 567*baa89ecbSBarry Smith Input Parameter: 56886697f06SMatthew Knepley . pc - the preconditioner context 56986697f06SMatthew Knepley 570*baa89ecbSBarry Smith Output Parameter: 571*baa89ecbSBarry Smith . type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 57286697f06SMatthew Knepley 57386697f06SMatthew Knepley Level: intermediate 57486697f06SMatthew Knepley 57586697f06SMatthew Knepley Concepts: Jacobi preconditioner 57686697f06SMatthew Knepley 577*baa89ecbSBarry Smith .seealso: PCJacobiaUseAbs(), PCJacobiSetType() 57886697f06SMatthew Knepley @*/ 579*baa89ecbSBarry Smith PetscErrorCode PCJacobiGetType(PC pc,PCJacobiType *type) 58086697f06SMatthew Knepley { 5814ac538c5SBarry Smith PetscErrorCode ierr; 58286697f06SMatthew Knepley 58386697f06SMatthew Knepley PetscFunctionBegin; 5840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 585*baa89ecbSBarry Smith ierr = PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));CHKERRQ(ierr); 58686697f06SMatthew Knepley PetscFunctionReturn(0); 58786697f06SMatthew Knepley } 588