1dba47a55SKris Buschelman #define PETSCKSP_DLL 24b9ad928SBarry Smith 34b9ad928SBarry Smith /* -------------------------------------------------------------------- 44b9ad928SBarry Smith 55a46d3faSBarry Smith This file implements a Jacobi preconditioner in PETSc as part of PC. 65a46d3faSBarry Smith You can use this as a starting point for implementing your own 75a46d3faSBarry Smith preconditioner that is not provided with PETSc. (You might also consider 85a46d3faSBarry Smith just using PCSHELL) 94b9ad928SBarry Smith 104b9ad928SBarry Smith The following basic routines are required for each preconditioner. 114b9ad928SBarry Smith PCCreate_XXX() - Creates a preconditioner context 124b9ad928SBarry Smith PCSetFromOptions_XXX() - Sets runtime options 134b9ad928SBarry Smith PCApply_XXX() - Applies the preconditioner 144b9ad928SBarry Smith PCDestroy_XXX() - Destroys the preconditioner context 154b9ad928SBarry Smith where the suffix "_XXX" denotes a particular implementation, in 164b9ad928SBarry Smith this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi). 174b9ad928SBarry Smith These routines are actually called via the common user interface 184b9ad928SBarry Smith routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 194b9ad928SBarry Smith so the application code interface remains identical for all 204b9ad928SBarry Smith preconditioners. 214b9ad928SBarry Smith 224b9ad928SBarry Smith Another key routine is: 234b9ad928SBarry Smith PCSetUp_XXX() - Prepares for the use of a preconditioner 244b9ad928SBarry Smith by setting data structures and options. The interface routine PCSetUp() 254b9ad928SBarry Smith is not usually called directly by the user, but instead is called by 264b9ad928SBarry Smith PCApply() if necessary. 274b9ad928SBarry Smith 284b9ad928SBarry Smith Additional basic routines are: 294b9ad928SBarry Smith PCView_XXX() - Prints details of runtime options that 304b9ad928SBarry Smith have actually been used. 314b9ad928SBarry Smith These are called by application codes via the interface routines 324b9ad928SBarry Smith PCView(). 334b9ad928SBarry Smith 344b9ad928SBarry Smith The various types of solvers (preconditioners, Krylov subspace methods, 354b9ad928SBarry Smith nonlinear solvers, timesteppers) are all organized similarly, so the 364b9ad928SBarry Smith above description applies to these categories also. One exception is 374b9ad928SBarry Smith that the analogues of PCApply() for these components are KSPSolve(), 384b9ad928SBarry Smith SNESSolve(), and TSSolve(). 394b9ad928SBarry Smith 404b9ad928SBarry Smith Additional optional functionality unique to preconditioners is left and 414b9ad928SBarry Smith right symmetric preconditioner application via PCApplySymmetricLeft() 424b9ad928SBarry Smith and PCApplySymmetricRight(). The Jacobi implementation is 434b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi(). 444b9ad928SBarry Smith 454b9ad928SBarry Smith -------------------------------------------------------------------- */ 464b9ad928SBarry Smith 474b9ad928SBarry Smith /* 484b9ad928SBarry Smith Include files needed for the Jacobi preconditioner: 494b9ad928SBarry Smith pcimpl.h - private include file intended for use by all preconditioners 504b9ad928SBarry Smith */ 514b9ad928SBarry Smith 526356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 534b9ad928SBarry Smith 544b9ad928SBarry Smith /* 554b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 564b9ad928SBarry Smith */ 574b9ad928SBarry Smith typedef struct { 584b9ad928SBarry Smith Vec diag; /* vector containing the reciprocals of the diagonal elements 594b9ad928SBarry Smith 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) */ 634b9ad928SBarry Smith PetscTruth userowmax; 6486697f06SMatthew Knepley PetscTruth userowsum; 65cd47f5d9SBarry Smith PetscTruth useabs; /* use the absolute values of the diagonal entries */ 664b9ad928SBarry Smith } PC_Jacobi; 674b9ad928SBarry Smith 684b9ad928SBarry Smith EXTERN_C_BEGIN 694b9ad928SBarry Smith #undef __FUNCT__ 704b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 71dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax_Jacobi(PC pc) 724b9ad928SBarry Smith { 734b9ad928SBarry Smith PC_Jacobi *j; 744b9ad928SBarry Smith 754b9ad928SBarry Smith PetscFunctionBegin; 764b9ad928SBarry Smith j = (PC_Jacobi*)pc->data; 774b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 784b9ad928SBarry Smith PetscFunctionReturn(0); 794b9ad928SBarry Smith } 804b9ad928SBarry Smith EXTERN_C_END 814b9ad928SBarry Smith 82cd47f5d9SBarry Smith EXTERN_C_BEGIN 83cd47f5d9SBarry Smith #undef __FUNCT__ 8486697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 8586697f06SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum_Jacobi(PC pc) 8686697f06SMatthew Knepley { 8786697f06SMatthew Knepley PC_Jacobi *j; 8886697f06SMatthew Knepley 8986697f06SMatthew Knepley PetscFunctionBegin; 9086697f06SMatthew Knepley j = (PC_Jacobi*)pc->data; 9186697f06SMatthew Knepley j->userowsum = PETSC_TRUE; 9286697f06SMatthew Knepley PetscFunctionReturn(0); 9386697f06SMatthew Knepley } 9486697f06SMatthew Knepley EXTERN_C_END 9586697f06SMatthew Knepley 9686697f06SMatthew Knepley EXTERN_C_BEGIN 9786697f06SMatthew Knepley #undef __FUNCT__ 98cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 99cd47f5d9SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs_Jacobi(PC pc) 100cd47f5d9SBarry Smith { 101cd47f5d9SBarry Smith PC_Jacobi *j; 102cd47f5d9SBarry Smith 103cd47f5d9SBarry Smith PetscFunctionBegin; 104cd47f5d9SBarry Smith j = (PC_Jacobi*)pc->data; 105cd47f5d9SBarry Smith j->useabs = PETSC_TRUE; 106cd47f5d9SBarry Smith PetscFunctionReturn(0); 107cd47f5d9SBarry Smith } 108cd47f5d9SBarry Smith EXTERN_C_END 109cd47f5d9SBarry Smith 1104b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1114b9ad928SBarry Smith /* 1124b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1134b9ad928SBarry Smith by setting data structures and options. 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith Input Parameter: 1164b9ad928SBarry Smith . pc - the preconditioner context 1174b9ad928SBarry Smith 1184b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith Notes: 1214b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1224b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1234b9ad928SBarry Smith */ 1244b9ad928SBarry Smith #undef __FUNCT__ 1254b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1266849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1274b9ad928SBarry Smith { 1284b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1294b9ad928SBarry Smith Vec diag,diagsqrt; 130dfbe8321SBarry Smith PetscErrorCode ierr; 131cd47f5d9SBarry Smith PetscInt n,i; 1324b9ad928SBarry Smith PetscScalar *x; 133cd47f5d9SBarry Smith PetscTruth zeroflag = PETSC_FALSE; 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith PetscFunctionBegin; 1364b9ad928SBarry Smith /* 1374b9ad928SBarry Smith For most preconditioners the code would begin here something like 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 14023ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1414b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1424b9ad928SBarry Smith } 1434b9ad928SBarry Smith 1444b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1454b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1464b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1474b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1484b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1494b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1504b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1514b9ad928SBarry Smith */ 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith /* 1544b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1554b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1564b9ad928SBarry Smith */ 1574b9ad928SBarry Smith diag = jac->diag; 1584b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith if (diag) { 1614b9ad928SBarry Smith if (jac->userowmax) { 162985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diag,PETSC_NULL);CHKERRQ(ierr); 16386697f06SMatthew Knepley } else if (jac->userowsum) { 16486697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 1654b9ad928SBarry Smith } else { 1664b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1674b9ad928SBarry Smith } 1684b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1694b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1704b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 171cd47f5d9SBarry Smith if (jac->useabs) { 172cd47f5d9SBarry Smith for (i=0; i<n; i++) { 173cd47f5d9SBarry Smith x[i] = PetscAbsScalar(x[i]); 174cd47f5d9SBarry Smith } 175cd47f5d9SBarry Smith } 1764b9ad928SBarry Smith for (i=0; i<n; i++) { 1774b9ad928SBarry Smith if (x[i] == 0.0) { 1784b9ad928SBarry Smith x[i] = 1.0; 179cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1804b9ad928SBarry Smith } 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith if (diagsqrt) { 1854b9ad928SBarry Smith if (jac->userowmax) { 186985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,PETSC_NULL);CHKERRQ(ierr); 18786697f06SMatthew Knepley } else if (jac->userowsum) { 18886697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 1894b9ad928SBarry Smith } else { 1904b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1914b9ad928SBarry Smith } 1924b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1934b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1944b9ad928SBarry Smith for (i=0; i<n; i++) { 1954b9ad928SBarry Smith if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 1964b9ad928SBarry Smith else { 1974b9ad928SBarry Smith x[i] = 1.0; 198cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1994b9ad928SBarry Smith } 2004b9ad928SBarry Smith } 2014b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 2024b9ad928SBarry Smith } 2034b9ad928SBarry Smith if (zeroflag) { 204ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 2054b9ad928SBarry Smith } 2064b9ad928SBarry Smith PetscFunctionReturn(0); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2094b9ad928SBarry Smith /* 2104b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2114b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2124b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith Input Parameter: 2154b9ad928SBarry Smith . pc - the preconditioner context 2164b9ad928SBarry Smith */ 2174b9ad928SBarry Smith #undef __FUNCT__ 2184b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2196849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2204b9ad928SBarry Smith { 221dfbe8321SBarry Smith PetscErrorCode ierr; 2224b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith PetscFunctionBegin; 22523ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 22652e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 2274b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2284b9ad928SBarry Smith PetscFunctionReturn(0); 2294b9ad928SBarry Smith } 2304b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2314b9ad928SBarry Smith /* 2324b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2334b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2344b9ad928SBarry Smith right application of the Jacobi preconditioner. 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith Input Parameter: 2374b9ad928SBarry Smith . pc - the preconditioner context 2384b9ad928SBarry Smith */ 2394b9ad928SBarry Smith #undef __FUNCT__ 2404b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 2416849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2424b9ad928SBarry Smith { 243dfbe8321SBarry Smith PetscErrorCode ierr; 2444b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith PetscFunctionBegin; 24723ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 24852e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 2494b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2504b9ad928SBarry Smith PetscFunctionReturn(0); 2514b9ad928SBarry Smith } 2524b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2534b9ad928SBarry Smith /* 2544b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith Input Parameters: 2574b9ad928SBarry Smith . pc - the preconditioner context 2584b9ad928SBarry Smith . x - input vector 2594b9ad928SBarry Smith 2604b9ad928SBarry Smith Output Parameter: 2614b9ad928SBarry Smith . y - output vector 2624b9ad928SBarry Smith 2634b9ad928SBarry Smith Application Interface Routine: PCApply() 2644b9ad928SBarry Smith */ 2654b9ad928SBarry Smith #undef __FUNCT__ 2664b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2676849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2684b9ad928SBarry Smith { 2694b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 270dfbe8321SBarry Smith PetscErrorCode ierr; 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith PetscFunctionBegin; 2734b9ad928SBarry Smith if (!jac->diag) { 2744b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2754b9ad928SBarry Smith } 2762dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 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 */ 2934b9ad928SBarry Smith #undef __FUNCT__ 2944b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 2956849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2964b9ad928SBarry Smith { 297dfbe8321SBarry Smith PetscErrorCode ierr; 2984b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2994b9ad928SBarry Smith 3004b9ad928SBarry Smith PetscFunctionBegin; 3014b9ad928SBarry Smith if (!jac->diagsqrt) { 3024b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 3034b9ad928SBarry Smith } 3042dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 3054b9ad928SBarry Smith PetscFunctionReturn(0); 3064b9ad928SBarry Smith } 3074b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3084b9ad928SBarry Smith /* 3094b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3104b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3114b9ad928SBarry Smith 3124b9ad928SBarry Smith Input Parameter: 3134b9ad928SBarry Smith . pc - the preconditioner context 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3164b9ad928SBarry Smith */ 3174b9ad928SBarry Smith #undef __FUNCT__ 3184b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3196849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3204b9ad928SBarry Smith { 3214b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 322dfbe8321SBarry Smith PetscErrorCode ierr; 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith PetscFunctionBegin; 3254b9ad928SBarry Smith if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 3264b9ad928SBarry Smith if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 3274b9ad928SBarry Smith 3284b9ad928SBarry Smith /* 3294b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3304b9ad928SBarry Smith */ 3314b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 3324b9ad928SBarry Smith PetscFunctionReturn(0); 3334b9ad928SBarry Smith } 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith #undef __FUNCT__ 3364b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3376849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 3384b9ad928SBarry Smith { 3394b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 340dfbe8321SBarry Smith PetscErrorCode ierr; 3414b9ad928SBarry Smith 3424b9ad928SBarry Smith PetscFunctionBegin; 3434b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 3444dc4c822SBarry Smith ierr = PetscOptionsTruth("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3454b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 34686697f06SMatthew Knepley ierr = PetscOptionsTruth("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 34786697f06SMatthew Knepley &jac->userowsum,PETSC_NULL);CHKERRQ(ierr); 348cd47f5d9SBarry Smith ierr = PetscOptionsTruth("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 349cd47f5d9SBarry Smith &jac->useabs,PETSC_NULL);CHKERRQ(ierr); 3504b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3514b9ad928SBarry Smith PetscFunctionReturn(0); 3524b9ad928SBarry Smith } 3534b9ad928SBarry Smith 3544b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3554b9ad928SBarry Smith /* 3564b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3574b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3584b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3594b9ad928SBarry Smith 3604b9ad928SBarry Smith Input Parameter: 3614b9ad928SBarry Smith . pc - the preconditioner context 3624b9ad928SBarry Smith 3634b9ad928SBarry Smith Application Interface Routine: PCCreate() 3644b9ad928SBarry Smith */ 3654b9ad928SBarry Smith 3664b9ad928SBarry Smith /*MC 3675a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3684b9ad928SBarry Smith 3694b9ad928SBarry Smith Options Database Key: 370cd47f5d9SBarry Smith + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3714b9ad928SBarry Smith rather than the diagonal 37286697f06SMatthew Knepley . -pc_jacobi_rowsum - use the maximum absolute value in each row as the scaling factor, 37386697f06SMatthew Knepley rather than the diagonal 374cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith Level: beginner 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 3814b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 38686697f06SMatthew Knepley PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 3874b9ad928SBarry Smith M*/ 3884b9ad928SBarry Smith 3894b9ad928SBarry Smith EXTERN_C_BEGIN 3904b9ad928SBarry Smith #undef __FUNCT__ 3914b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 392dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Jacobi(PC pc) 3934b9ad928SBarry Smith { 3944b9ad928SBarry Smith PC_Jacobi *jac; 395dfbe8321SBarry Smith PetscErrorCode ierr; 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith PetscFunctionBegin; 3984b9ad928SBarry Smith /* 3994b9ad928SBarry Smith Creates the private data structure for this preconditioner and 4004b9ad928SBarry Smith attach it to the PC object. 4014b9ad928SBarry Smith */ 40238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Jacobi,&jac);CHKERRQ(ierr); 4034b9ad928SBarry Smith pc->data = (void*)jac; 4044b9ad928SBarry Smith 4054b9ad928SBarry Smith /* 4064b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4074b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4084b9ad928SBarry Smith */ 4094b9ad928SBarry Smith jac->diag = 0; 4104b9ad928SBarry Smith jac->diagsqrt = 0; 4114b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 41286697f06SMatthew Knepley jac->userowsum = PETSC_FALSE; 413cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith /* 4164b9ad928SBarry Smith Set the pointers for the functions that are provided above. 4174b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4184b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4194b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4204b9ad928SBarry Smith not needed. 4214b9ad928SBarry Smith */ 4224b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4234b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4244b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 4254b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4264b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4274b9ad928SBarry Smith pc->ops->view = 0; 4284b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4294b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4304b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 431cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 43286697f06SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 433cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 4344b9ad928SBarry Smith PetscFunctionReturn(0); 4354b9ad928SBarry Smith } 4364b9ad928SBarry Smith EXTERN_C_END 4374b9ad928SBarry Smith 438cd47f5d9SBarry Smith 439cd47f5d9SBarry Smith #undef __FUNCT__ 440cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs" 441cd47f5d9SBarry Smith /*@ 442cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 443cd47f5d9SBarry Smith absolute value of the diagonal to for the preconditioner 444cd47f5d9SBarry Smith 445cd47f5d9SBarry Smith Collective on PC 446cd47f5d9SBarry Smith 447cd47f5d9SBarry Smith Input Parameters: 448cd47f5d9SBarry Smith . pc - the preconditioner context 449cd47f5d9SBarry Smith 450cd47f5d9SBarry Smith 451cd47f5d9SBarry Smith Options Database Key: 452cd47f5d9SBarry Smith . -pc_jacobi_abs 453cd47f5d9SBarry Smith 454cd47f5d9SBarry Smith Level: intermediate 455cd47f5d9SBarry Smith 456cd47f5d9SBarry Smith Concepts: Jacobi preconditioner 457cd47f5d9SBarry Smith 45886697f06SMatthew Knepley .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum() 459cd47f5d9SBarry Smith 460cd47f5d9SBarry Smith @*/ 461cd47f5d9SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC pc) 462cd47f5d9SBarry Smith { 463cd47f5d9SBarry Smith PetscErrorCode ierr,(*f)(PC); 464cd47f5d9SBarry Smith 465cd47f5d9SBarry Smith PetscFunctionBegin; 466*0700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 467cd47f5d9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",(void (**)(void))&f);CHKERRQ(ierr); 468cd47f5d9SBarry Smith if (f) { 469cd47f5d9SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 470cd47f5d9SBarry Smith } 471cd47f5d9SBarry Smith PetscFunctionReturn(0); 472cd47f5d9SBarry Smith } 473cd47f5d9SBarry Smith 4744b9ad928SBarry Smith #undef __FUNCT__ 4754b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4764b9ad928SBarry Smith /*@ 4774b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4784b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4794b9ad928SBarry Smith the diagonal entry 4804b9ad928SBarry Smith 4814b9ad928SBarry Smith Collective on PC 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith Input Parameters: 4844b9ad928SBarry Smith . pc - the preconditioner context 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith 4874b9ad928SBarry Smith Options Database Key: 4884b9ad928SBarry Smith . -pc_jacobi_rowmax 4894b9ad928SBarry Smith 4904b9ad928SBarry Smith Level: intermediate 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith Concepts: Jacobi preconditioner 4934b9ad928SBarry Smith 494cd47f5d9SBarry Smith .seealso: PCJacobiaUseAbs() 4954b9ad928SBarry Smith @*/ 496dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC pc) 4974b9ad928SBarry Smith { 498dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 4994b9ad928SBarry Smith 5004b9ad928SBarry Smith PetscFunctionBegin; 501*0700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 502e34fafa9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 5034b9ad928SBarry Smith if (f) { 5044b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 5054b9ad928SBarry Smith } 5064b9ad928SBarry Smith PetscFunctionReturn(0); 5074b9ad928SBarry Smith } 5084b9ad928SBarry Smith 50986697f06SMatthew Knepley #undef __FUNCT__ 51086697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum" 51186697f06SMatthew Knepley /*@ 51286697f06SMatthew Knepley PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 51386697f06SMatthew Knepley sum of each row as the diagonal preconditioner, instead of 51486697f06SMatthew Knepley the diagonal entry 51586697f06SMatthew Knepley 51686697f06SMatthew Knepley Collective on PC 51786697f06SMatthew Knepley 51886697f06SMatthew Knepley Input Parameters: 51986697f06SMatthew Knepley . pc - the preconditioner context 52086697f06SMatthew Knepley 52186697f06SMatthew Knepley 52286697f06SMatthew Knepley Options Database Key: 52386697f06SMatthew Knepley . -pc_jacobi_rowsum 52486697f06SMatthew Knepley 52586697f06SMatthew Knepley Level: intermediate 52686697f06SMatthew Knepley 52786697f06SMatthew Knepley Concepts: Jacobi preconditioner 52886697f06SMatthew Knepley 52986697f06SMatthew Knepley .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 53086697f06SMatthew Knepley @*/ 53186697f06SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC pc) 53286697f06SMatthew Knepley { 53386697f06SMatthew Knepley PetscErrorCode ierr,(*f)(PC); 53486697f06SMatthew Knepley 53586697f06SMatthew Knepley PetscFunctionBegin; 536*0700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 53786697f06SMatthew Knepley ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseRowSum_C",(void (**)(void))&f);CHKERRQ(ierr); 53886697f06SMatthew Knepley if (f) { 53986697f06SMatthew Knepley ierr = (*f)(pc);CHKERRQ(ierr); 54086697f06SMatthew Knepley } 54186697f06SMatthew Knepley PetscFunctionReturn(0); 54286697f06SMatthew Knepley } 54386697f06SMatthew Knepley 544