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 51*c6db04a5SJed Brown #include <private/pcimpl.h> /*I "petscpc.h" I*/ 524b9ad928SBarry Smith 534b9ad928SBarry Smith /* 544b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 554b9ad928SBarry Smith */ 564b9ad928SBarry Smith typedef struct { 574b9ad928SBarry Smith Vec diag; /* vector containing the reciprocals of the diagonal elements 584b9ad928SBarry Smith of the preconditioner matrix */ 594b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 604b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 614b9ad928SBarry Smith only for symmetric preconditioner application) */ 62ace3abfcSBarry Smith PetscBool userowmax; 63ace3abfcSBarry Smith PetscBool userowsum; 64ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */ 654b9ad928SBarry Smith } PC_Jacobi; 664b9ad928SBarry Smith 674b9ad928SBarry Smith EXTERN_C_BEGIN 684b9ad928SBarry Smith #undef __FUNCT__ 694b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 707087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 714b9ad928SBarry Smith { 724b9ad928SBarry Smith PC_Jacobi *j; 734b9ad928SBarry Smith 744b9ad928SBarry Smith PetscFunctionBegin; 754b9ad928SBarry Smith j = (PC_Jacobi*)pc->data; 764b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 774b9ad928SBarry Smith PetscFunctionReturn(0); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith EXTERN_C_END 804b9ad928SBarry Smith 81cd47f5d9SBarry Smith EXTERN_C_BEGIN 82cd47f5d9SBarry Smith #undef __FUNCT__ 8386697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 847087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowSum_Jacobi(PC pc) 8586697f06SMatthew Knepley { 8686697f06SMatthew Knepley PC_Jacobi *j; 8786697f06SMatthew Knepley 8886697f06SMatthew Knepley PetscFunctionBegin; 8986697f06SMatthew Knepley j = (PC_Jacobi*)pc->data; 9086697f06SMatthew Knepley j->userowsum = PETSC_TRUE; 9186697f06SMatthew Knepley PetscFunctionReturn(0); 9286697f06SMatthew Knepley } 9386697f06SMatthew Knepley EXTERN_C_END 9486697f06SMatthew Knepley 9586697f06SMatthew Knepley EXTERN_C_BEGIN 9686697f06SMatthew Knepley #undef __FUNCT__ 97cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 987087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc) 99cd47f5d9SBarry Smith { 100cd47f5d9SBarry Smith PC_Jacobi *j; 101cd47f5d9SBarry Smith 102cd47f5d9SBarry Smith PetscFunctionBegin; 103cd47f5d9SBarry Smith j = (PC_Jacobi*)pc->data; 104cd47f5d9SBarry Smith j->useabs = PETSC_TRUE; 105cd47f5d9SBarry Smith PetscFunctionReturn(0); 106cd47f5d9SBarry Smith } 107cd47f5d9SBarry Smith EXTERN_C_END 108cd47f5d9SBarry Smith 1094b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1104b9ad928SBarry Smith /* 1114b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1124b9ad928SBarry Smith by setting data structures and options. 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith Input Parameter: 1154b9ad928SBarry Smith . pc - the preconditioner context 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith Notes: 1204b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1214b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1224b9ad928SBarry Smith */ 1234b9ad928SBarry Smith #undef __FUNCT__ 1244b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1256849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1264b9ad928SBarry Smith { 1274b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1284b9ad928SBarry Smith Vec diag,diagsqrt; 129dfbe8321SBarry Smith PetscErrorCode ierr; 130cd47f5d9SBarry Smith PetscInt n,i; 1314b9ad928SBarry Smith PetscScalar *x; 132ace3abfcSBarry Smith PetscBool zeroflag = PETSC_FALSE; 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith PetscFunctionBegin; 1354b9ad928SBarry Smith /* 1364b9ad928SBarry Smith For most preconditioners the code would begin here something like 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 13923ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1404b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1414b9ad928SBarry Smith } 1424b9ad928SBarry Smith 1434b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1444b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1454b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1464b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1474b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1484b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1494b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1504b9ad928SBarry Smith */ 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith /* 1534b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1544b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1554b9ad928SBarry Smith */ 1564b9ad928SBarry Smith diag = jac->diag; 1574b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith if (diag) { 1604b9ad928SBarry Smith if (jac->userowmax) { 161985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diag,PETSC_NULL);CHKERRQ(ierr); 16286697f06SMatthew Knepley } else if (jac->userowsum) { 16386697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 1644b9ad928SBarry Smith } else { 1654b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1664b9ad928SBarry Smith } 1674b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1684b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1694b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 170cd47f5d9SBarry Smith if (jac->useabs) { 171cd47f5d9SBarry Smith for (i=0; i<n; i++) { 172cd47f5d9SBarry Smith x[i] = PetscAbsScalar(x[i]); 173cd47f5d9SBarry Smith } 174cd47f5d9SBarry Smith } 1754b9ad928SBarry Smith for (i=0; i<n; i++) { 1764b9ad928SBarry Smith if (x[i] == 0.0) { 1774b9ad928SBarry Smith x[i] = 1.0; 178cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1794b9ad928SBarry Smith } 1804b9ad928SBarry Smith } 1814b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1824b9ad928SBarry Smith } 1834b9ad928SBarry Smith if (diagsqrt) { 1844b9ad928SBarry Smith if (jac->userowmax) { 185985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,PETSC_NULL);CHKERRQ(ierr); 18686697f06SMatthew Knepley } else if (jac->userowsum) { 18786697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 1884b9ad928SBarry Smith } else { 1894b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1904b9ad928SBarry Smith } 1914b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1924b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1934b9ad928SBarry Smith for (i=0; i<n; i++) { 1944b9ad928SBarry Smith if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 1954b9ad928SBarry Smith else { 1964b9ad928SBarry Smith x[i] = 1.0; 197cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1984b9ad928SBarry Smith } 1994b9ad928SBarry Smith } 2004b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 2014b9ad928SBarry Smith } 2024b9ad928SBarry Smith if (zeroflag) { 203ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 2044b9ad928SBarry Smith } 2054b9ad928SBarry Smith PetscFunctionReturn(0); 2064b9ad928SBarry Smith } 2074b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2084b9ad928SBarry Smith /* 2094b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2104b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2114b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2124b9ad928SBarry Smith 2134b9ad928SBarry Smith Input Parameter: 2144b9ad928SBarry Smith . pc - the preconditioner context 2154b9ad928SBarry Smith */ 2164b9ad928SBarry Smith #undef __FUNCT__ 2174b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2186849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2194b9ad928SBarry Smith { 220dfbe8321SBarry Smith PetscErrorCode ierr; 2214b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith PetscFunctionBegin; 22423ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 22552e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 2264b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2274b9ad928SBarry Smith PetscFunctionReturn(0); 2284b9ad928SBarry Smith } 2294b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2304b9ad928SBarry Smith /* 2314b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2324b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2334b9ad928SBarry Smith right application of the Jacobi preconditioner. 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith Input Parameter: 2364b9ad928SBarry Smith . pc - the preconditioner context 2374b9ad928SBarry Smith */ 2384b9ad928SBarry Smith #undef __FUNCT__ 2394b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 2406849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2414b9ad928SBarry Smith { 242dfbe8321SBarry Smith PetscErrorCode ierr; 2434b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2444b9ad928SBarry Smith 2454b9ad928SBarry Smith PetscFunctionBegin; 24623ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 24752e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 2484b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2494b9ad928SBarry Smith PetscFunctionReturn(0); 2504b9ad928SBarry Smith } 2514b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2524b9ad928SBarry Smith /* 2534b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2544b9ad928SBarry Smith 2554b9ad928SBarry Smith Input Parameters: 2564b9ad928SBarry Smith . pc - the preconditioner context 2574b9ad928SBarry Smith . x - input vector 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith Output Parameter: 2604b9ad928SBarry Smith . y - output vector 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith Application Interface Routine: PCApply() 2634b9ad928SBarry Smith */ 2644b9ad928SBarry Smith #undef __FUNCT__ 2654b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2666849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2674b9ad928SBarry Smith { 2684b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 269dfbe8321SBarry Smith PetscErrorCode ierr; 2704b9ad928SBarry Smith 2714b9ad928SBarry Smith PetscFunctionBegin; 2724b9ad928SBarry Smith if (!jac->diag) { 2734b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2744b9ad928SBarry Smith } 2752dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 2764b9ad928SBarry Smith PetscFunctionReturn(0); 2774b9ad928SBarry Smith } 2784b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2794b9ad928SBarry Smith /* 2804b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2814b9ad928SBarry Smith symmetric preconditioner to a vector. 2824b9ad928SBarry Smith 2834b9ad928SBarry Smith Input Parameters: 2844b9ad928SBarry Smith . pc - the preconditioner context 2854b9ad928SBarry Smith . x - input vector 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith Output Parameter: 2884b9ad928SBarry Smith . y - output vector 2894b9ad928SBarry Smith 2904b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2914b9ad928SBarry Smith */ 2924b9ad928SBarry Smith #undef __FUNCT__ 2934b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 2946849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2954b9ad928SBarry Smith { 296dfbe8321SBarry Smith PetscErrorCode ierr; 2974b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith PetscFunctionBegin; 3004b9ad928SBarry Smith if (!jac->diagsqrt) { 3014b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 3024b9ad928SBarry Smith } 3032dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 3044b9ad928SBarry Smith PetscFunctionReturn(0); 3054b9ad928SBarry Smith } 3064b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3074b9ad928SBarry Smith /* 3084b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3094b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3104b9ad928SBarry Smith 3114b9ad928SBarry Smith Input Parameter: 3124b9ad928SBarry Smith . pc - the preconditioner context 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3154b9ad928SBarry Smith */ 3164b9ad928SBarry Smith #undef __FUNCT__ 3174b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3186849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3194b9ad928SBarry Smith { 3204b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 321dfbe8321SBarry Smith PetscErrorCode ierr; 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith PetscFunctionBegin; 3244b9ad928SBarry Smith if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 3254b9ad928SBarry Smith if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith /* 3284b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3294b9ad928SBarry Smith */ 330c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 3314b9ad928SBarry Smith PetscFunctionReturn(0); 3324b9ad928SBarry Smith } 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith #undef __FUNCT__ 3354b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3366849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 3374b9ad928SBarry Smith { 3384b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 339dfbe8321SBarry Smith PetscErrorCode ierr; 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith PetscFunctionBegin; 3424b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 343acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3444b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 345acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 34686697f06SMatthew Knepley &jac->userowsum,PETSC_NULL);CHKERRQ(ierr); 347acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 348cd47f5d9SBarry Smith &jac->useabs,PETSC_NULL);CHKERRQ(ierr); 3494b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3504b9ad928SBarry Smith PetscFunctionReturn(0); 3514b9ad928SBarry Smith } 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3544b9ad928SBarry Smith /* 3554b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3564b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3574b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith Input Parameter: 3604b9ad928SBarry Smith . pc - the preconditioner context 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith Application Interface Routine: PCCreate() 3634b9ad928SBarry Smith */ 3644b9ad928SBarry Smith 3654b9ad928SBarry Smith /*MC 3665a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3674b9ad928SBarry Smith 3684b9ad928SBarry Smith Options Database Key: 369cd47f5d9SBarry Smith + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3704b9ad928SBarry Smith rather than the diagonal 37186697f06SMatthew Knepley . -pc_jacobi_rowsum - use the maximum absolute value in each row as the scaling factor, 37286697f06SMatthew Knepley rather than the diagonal 373cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith Level: beginner 3764b9ad928SBarry Smith 3774b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3784b9ad928SBarry Smith 379b037da10SBarry Smith Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 3804b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 38586697f06SMatthew Knepley PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 3864b9ad928SBarry Smith M*/ 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith EXTERN_C_BEGIN 3894b9ad928SBarry Smith #undef __FUNCT__ 3904b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 3917087cfbeSBarry Smith PetscErrorCode PCCreate_Jacobi(PC pc) 3924b9ad928SBarry Smith { 3934b9ad928SBarry Smith PC_Jacobi *jac; 394dfbe8321SBarry Smith PetscErrorCode ierr; 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith PetscFunctionBegin; 3974b9ad928SBarry Smith /* 3984b9ad928SBarry Smith Creates the private data structure for this preconditioner and 3994b9ad928SBarry Smith attach it to the PC object. 4004b9ad928SBarry Smith */ 40138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Jacobi,&jac);CHKERRQ(ierr); 4024b9ad928SBarry Smith pc->data = (void*)jac; 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith /* 4054b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4064b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4074b9ad928SBarry Smith */ 4084b9ad928SBarry Smith jac->diag = 0; 4094b9ad928SBarry Smith jac->diagsqrt = 0; 4104b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 41186697f06SMatthew Knepley jac->userowsum = PETSC_FALSE; 412cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith /* 4154b9ad928SBarry Smith Set the pointers for the functions that are provided above. 4164b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4174b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4184b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4194b9ad928SBarry Smith not needed. 4204b9ad928SBarry Smith */ 4214b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4224b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4234b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 4244b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4254b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4264b9ad928SBarry Smith pc->ops->view = 0; 4274b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4284b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4294b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 430cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 43186697f06SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 432cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 4334b9ad928SBarry Smith PetscFunctionReturn(0); 4344b9ad928SBarry Smith } 4354b9ad928SBarry Smith EXTERN_C_END 4364b9ad928SBarry Smith 437cd47f5d9SBarry Smith 438cd47f5d9SBarry Smith #undef __FUNCT__ 439cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs" 440cd47f5d9SBarry Smith /*@ 441cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 442cd47f5d9SBarry Smith absolute value of the diagonal to for the preconditioner 443cd47f5d9SBarry Smith 444ad4df100SBarry Smith Logically Collective on PC 445cd47f5d9SBarry Smith 446cd47f5d9SBarry Smith Input Parameters: 447cd47f5d9SBarry Smith . pc - the preconditioner context 448cd47f5d9SBarry Smith 449cd47f5d9SBarry Smith 450cd47f5d9SBarry Smith Options Database Key: 451cd47f5d9SBarry Smith . -pc_jacobi_abs 452cd47f5d9SBarry Smith 453cd47f5d9SBarry Smith Level: intermediate 454cd47f5d9SBarry Smith 455cd47f5d9SBarry Smith Concepts: Jacobi preconditioner 456cd47f5d9SBarry Smith 45786697f06SMatthew Knepley .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum() 458cd47f5d9SBarry Smith 459cd47f5d9SBarry Smith @*/ 4607087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseAbs(PC pc) 461cd47f5d9SBarry Smith { 4624ac538c5SBarry Smith PetscErrorCode ierr; 463cd47f5d9SBarry Smith 464cd47f5d9SBarry Smith PetscFunctionBegin; 4650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4664ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC),(pc));CHKERRQ(ierr); 467cd47f5d9SBarry Smith PetscFunctionReturn(0); 468cd47f5d9SBarry Smith } 469cd47f5d9SBarry Smith 4704b9ad928SBarry Smith #undef __FUNCT__ 4714b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4724b9ad928SBarry Smith /*@ 4734b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4744b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4754b9ad928SBarry Smith the diagonal entry 4764b9ad928SBarry Smith 477ad4df100SBarry Smith Logically Collective on PC 4784b9ad928SBarry Smith 4794b9ad928SBarry Smith Input Parameters: 4804b9ad928SBarry Smith . pc - the preconditioner context 4814b9ad928SBarry Smith 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith Options Database Key: 4844b9ad928SBarry Smith . -pc_jacobi_rowmax 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith Level: intermediate 4874b9ad928SBarry Smith 4884b9ad928SBarry Smith Concepts: Jacobi preconditioner 4894b9ad928SBarry Smith 490cd47f5d9SBarry Smith .seealso: PCJacobiaUseAbs() 4914b9ad928SBarry Smith @*/ 4927087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowMax(PC pc) 4934b9ad928SBarry Smith { 4944ac538c5SBarry Smith PetscErrorCode ierr; 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith PetscFunctionBegin; 4970700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4984ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowMax_C",(PC),(pc));CHKERRQ(ierr); 4994b9ad928SBarry Smith PetscFunctionReturn(0); 5004b9ad928SBarry Smith } 5014b9ad928SBarry Smith 50286697f06SMatthew Knepley #undef __FUNCT__ 50386697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum" 50486697f06SMatthew Knepley /*@ 50586697f06SMatthew Knepley PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 50686697f06SMatthew Knepley sum of each row as the diagonal preconditioner, instead of 50786697f06SMatthew Knepley the diagonal entry 50886697f06SMatthew Knepley 509ad4df100SBarry Smith Logical Collective on PC 51086697f06SMatthew Knepley 51186697f06SMatthew Knepley Input Parameters: 51286697f06SMatthew Knepley . pc - the preconditioner context 51386697f06SMatthew Knepley 51486697f06SMatthew Knepley 51586697f06SMatthew Knepley Options Database Key: 51686697f06SMatthew Knepley . -pc_jacobi_rowsum 51786697f06SMatthew Knepley 51886697f06SMatthew Knepley Level: intermediate 51986697f06SMatthew Knepley 52086697f06SMatthew Knepley Concepts: Jacobi preconditioner 52186697f06SMatthew Knepley 52286697f06SMatthew Knepley .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 52386697f06SMatthew Knepley @*/ 5247087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowSum(PC pc) 52586697f06SMatthew Knepley { 5264ac538c5SBarry Smith PetscErrorCode ierr; 52786697f06SMatthew Knepley 52886697f06SMatthew Knepley PetscFunctionBegin; 5290700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5304ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowSum_C",(PC),(pc));CHKERRQ(ierr); 53186697f06SMatthew Knepley PetscFunctionReturn(0); 53286697f06SMatthew Knepley } 53386697f06SMatthew Knepley 534