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 534b9ad928SBarry Smith /* 544b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 554b9ad928SBarry Smith */ 564b9ad928SBarry Smith typedef struct { 57*2fa5cd67SKarl Rupp Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 584b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 594b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 604b9ad928SBarry Smith only for symmetric preconditioner application) */ 61ace3abfcSBarry Smith PetscBool userowmax; 62ace3abfcSBarry Smith PetscBool userowsum; 63ace3abfcSBarry Smith PetscBool useabs; /* use the absolute values of the diagonal entries */ 644b9ad928SBarry Smith } PC_Jacobi; 654b9ad928SBarry Smith 664b9ad928SBarry Smith EXTERN_C_BEGIN 674b9ad928SBarry Smith #undef __FUNCT__ 684b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 697087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 704b9ad928SBarry Smith { 714b9ad928SBarry Smith PC_Jacobi *j; 724b9ad928SBarry Smith 734b9ad928SBarry Smith PetscFunctionBegin; 744b9ad928SBarry Smith j = (PC_Jacobi*)pc->data; 754b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith EXTERN_C_END 794b9ad928SBarry Smith 80cd47f5d9SBarry Smith EXTERN_C_BEGIN 81cd47f5d9SBarry Smith #undef __FUNCT__ 8286697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 837087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowSum_Jacobi(PC pc) 8486697f06SMatthew Knepley { 8586697f06SMatthew Knepley PC_Jacobi *j; 8686697f06SMatthew Knepley 8786697f06SMatthew Knepley PetscFunctionBegin; 8886697f06SMatthew Knepley j = (PC_Jacobi*)pc->data; 8986697f06SMatthew Knepley j->userowsum = PETSC_TRUE; 9086697f06SMatthew Knepley PetscFunctionReturn(0); 9186697f06SMatthew Knepley } 9286697f06SMatthew Knepley EXTERN_C_END 9386697f06SMatthew Knepley 9486697f06SMatthew Knepley EXTERN_C_BEGIN 9586697f06SMatthew Knepley #undef __FUNCT__ 96cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 977087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc) 98cd47f5d9SBarry Smith { 99cd47f5d9SBarry Smith PC_Jacobi *j; 100cd47f5d9SBarry Smith 101cd47f5d9SBarry Smith PetscFunctionBegin; 102cd47f5d9SBarry Smith j = (PC_Jacobi*)pc->data; 103cd47f5d9SBarry Smith j->useabs = PETSC_TRUE; 104cd47f5d9SBarry Smith PetscFunctionReturn(0); 105cd47f5d9SBarry Smith } 106cd47f5d9SBarry Smith EXTERN_C_END 107cd47f5d9SBarry Smith 1084b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1094b9ad928SBarry Smith /* 1104b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1114b9ad928SBarry Smith by setting data structures and options. 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith Input Parameter: 1144b9ad928SBarry Smith . pc - the preconditioner context 1154b9ad928SBarry Smith 1164b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1174b9ad928SBarry Smith 1184b9ad928SBarry Smith Notes: 1194b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1204b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1214b9ad928SBarry Smith */ 1224b9ad928SBarry Smith #undef __FUNCT__ 1234b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1246849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1254b9ad928SBarry Smith { 1264b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1274b9ad928SBarry Smith Vec diag,diagsqrt; 128dfbe8321SBarry Smith PetscErrorCode ierr; 129cd47f5d9SBarry Smith PetscInt n,i; 1304b9ad928SBarry Smith PetscScalar *x; 131ace3abfcSBarry Smith PetscBool zeroflag = PETSC_FALSE; 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith PetscFunctionBegin; 1344b9ad928SBarry Smith /* 1354b9ad928SBarry Smith For most preconditioners the code would begin here something like 1364b9ad928SBarry Smith 1374b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 13823ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1394b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1404b9ad928SBarry Smith } 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1434b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1444b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1454b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1464b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1474b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1484b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1494b9ad928SBarry Smith */ 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith /* 1524b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1534b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1544b9ad928SBarry Smith */ 1554b9ad928SBarry Smith diag = jac->diag; 1564b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith if (diag) { 1594b9ad928SBarry Smith if (jac->userowmax) { 160985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diag,PETSC_NULL);CHKERRQ(ierr); 16186697f06SMatthew Knepley } else if (jac->userowsum) { 16286697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 1634b9ad928SBarry Smith } else { 1644b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1654b9ad928SBarry Smith } 1664b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1674b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1684b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 169cd47f5d9SBarry Smith if (jac->useabs) { 170*2fa5cd67SKarl Rupp for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 171cd47f5d9SBarry Smith } 1724b9ad928SBarry Smith for (i=0; i<n; i++) { 1734b9ad928SBarry Smith if (x[i] == 0.0) { 1744b9ad928SBarry Smith x[i] = 1.0; 175cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1764b9ad928SBarry Smith } 1774b9ad928SBarry Smith } 1784b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1794b9ad928SBarry Smith } 1804b9ad928SBarry Smith if (diagsqrt) { 1814b9ad928SBarry Smith if (jac->userowmax) { 182985db425SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,PETSC_NULL);CHKERRQ(ierr); 18386697f06SMatthew Knepley } else if (jac->userowsum) { 18486697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 1854b9ad928SBarry Smith } else { 1864b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1874b9ad928SBarry Smith } 1884b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1894b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1904b9ad928SBarry Smith for (i=0; i<n; i++) { 1918f1a2a5eSBarry Smith if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 1924b9ad928SBarry Smith else { 1934b9ad928SBarry Smith x[i] = 1.0; 194cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith } 1974b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 1984b9ad928SBarry Smith } 1994b9ad928SBarry Smith if (zeroflag) { 200ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 2014b9ad928SBarry Smith } 2024b9ad928SBarry Smith PetscFunctionReturn(0); 2034b9ad928SBarry Smith } 2044b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2054b9ad928SBarry Smith /* 2064b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2074b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2084b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith Input Parameter: 2114b9ad928SBarry Smith . pc - the preconditioner context 2124b9ad928SBarry Smith */ 2134b9ad928SBarry Smith #undef __FUNCT__ 2144b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2156849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2164b9ad928SBarry Smith { 217dfbe8321SBarry Smith PetscErrorCode ierr; 2184b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2194b9ad928SBarry Smith 2204b9ad928SBarry Smith PetscFunctionBegin; 22123ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 22252e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 2234b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2244b9ad928SBarry Smith PetscFunctionReturn(0); 2254b9ad928SBarry Smith } 2264b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2274b9ad928SBarry Smith /* 2284b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2294b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2304b9ad928SBarry Smith right application of the Jacobi preconditioner. 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith Input Parameter: 2334b9ad928SBarry Smith . pc - the preconditioner context 2344b9ad928SBarry Smith */ 2354b9ad928SBarry Smith #undef __FUNCT__ 2364b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 2376849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2384b9ad928SBarry Smith { 239dfbe8321SBarry Smith PetscErrorCode ierr; 2404b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith PetscFunctionBegin; 24323ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 24452e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 2454b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2464b9ad928SBarry Smith PetscFunctionReturn(0); 2474b9ad928SBarry Smith } 2484b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2494b9ad928SBarry Smith /* 2504b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2514b9ad928SBarry Smith 2524b9ad928SBarry Smith Input Parameters: 2534b9ad928SBarry Smith . pc - the preconditioner context 2544b9ad928SBarry Smith . x - input vector 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith Output Parameter: 2574b9ad928SBarry Smith . y - output vector 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith Application Interface Routine: PCApply() 2604b9ad928SBarry Smith */ 2614b9ad928SBarry Smith #undef __FUNCT__ 2624b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2636849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2644b9ad928SBarry Smith { 2654b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 266dfbe8321SBarry Smith PetscErrorCode ierr; 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith PetscFunctionBegin; 2694b9ad928SBarry Smith if (!jac->diag) { 2704b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2714b9ad928SBarry Smith } 2722dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 2734b9ad928SBarry Smith PetscFunctionReturn(0); 2744b9ad928SBarry Smith } 2754b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2764b9ad928SBarry Smith /* 2774b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2784b9ad928SBarry Smith symmetric preconditioner to a vector. 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith Input Parameters: 2814b9ad928SBarry Smith . pc - the preconditioner context 2824b9ad928SBarry Smith . x - input vector 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith Output Parameter: 2854b9ad928SBarry Smith . y - output vector 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2884b9ad928SBarry Smith */ 2894b9ad928SBarry Smith #undef __FUNCT__ 2904b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 2916849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2924b9ad928SBarry Smith { 293dfbe8321SBarry Smith PetscErrorCode ierr; 2944b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith PetscFunctionBegin; 2974b9ad928SBarry Smith if (!jac->diagsqrt) { 2984b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 2994b9ad928SBarry Smith } 3002dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 3014b9ad928SBarry Smith PetscFunctionReturn(0); 3024b9ad928SBarry Smith } 3034b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 304a06653b4SBarry Smith #undef __FUNCT__ 305a06653b4SBarry Smith #define __FUNCT__ "PCReset_Jacobi" 306a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc) 307a06653b4SBarry Smith { 308a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 309a06653b4SBarry Smith PetscErrorCode ierr; 310a06653b4SBarry Smith 311a06653b4SBarry Smith PetscFunctionBegin; 3126bf464f9SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 3136bf464f9SBarry Smith ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 314a06653b4SBarry Smith PetscFunctionReturn(0); 315a06653b4SBarry Smith } 316a06653b4SBarry Smith 3174b9ad928SBarry Smith /* 3184b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3194b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3204b9ad928SBarry Smith 3214b9ad928SBarry Smith Input Parameter: 3224b9ad928SBarry Smith . pc - the preconditioner context 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3254b9ad928SBarry Smith */ 3264b9ad928SBarry Smith #undef __FUNCT__ 3274b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3286849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3294b9ad928SBarry Smith { 330dfbe8321SBarry Smith PetscErrorCode ierr; 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith PetscFunctionBegin; 333a06653b4SBarry Smith ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith /* 3364b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3374b9ad928SBarry Smith */ 338c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 3394b9ad928SBarry Smith PetscFunctionReturn(0); 3404b9ad928SBarry Smith } 3414b9ad928SBarry Smith 3424b9ad928SBarry Smith #undef __FUNCT__ 3434b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3446849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 3454b9ad928SBarry Smith { 3464b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 347dfbe8321SBarry Smith PetscErrorCode ierr; 3484b9ad928SBarry Smith 3494b9ad928SBarry Smith PetscFunctionBegin; 3504b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 351acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3524b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 353acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 35486697f06SMatthew Knepley &jac->userowsum,PETSC_NULL);CHKERRQ(ierr); 355acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 356cd47f5d9SBarry Smith &jac->useabs,PETSC_NULL);CHKERRQ(ierr); 3574b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3584b9ad928SBarry Smith PetscFunctionReturn(0); 3594b9ad928SBarry Smith } 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3624b9ad928SBarry Smith /* 3634b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3644b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3654b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith Input Parameter: 3684b9ad928SBarry Smith . pc - the preconditioner context 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Application Interface Routine: PCCreate() 3714b9ad928SBarry Smith */ 3724b9ad928SBarry Smith 3734b9ad928SBarry Smith /*MC 3745a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith Options Database Key: 377cd47f5d9SBarry Smith + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3784b9ad928SBarry Smith rather than the diagonal 379079a340bSHong Zhang . -pc_jacobi_rowsum - use the row sums as the scaling factor, 38086697f06SMatthew Knepley rather than the diagonal 381cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Level: beginner 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3864b9ad928SBarry Smith 387b037da10SBarry Smith Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 3884b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 39386697f06SMatthew Knepley PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 3944b9ad928SBarry Smith M*/ 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith EXTERN_C_BEGIN 3974b9ad928SBarry Smith #undef __FUNCT__ 3984b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 3997087cfbeSBarry Smith PetscErrorCode PCCreate_Jacobi(PC pc) 4004b9ad928SBarry Smith { 4014b9ad928SBarry Smith PC_Jacobi *jac; 402dfbe8321SBarry Smith PetscErrorCode ierr; 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith PetscFunctionBegin; 4054b9ad928SBarry Smith /* 4064b9ad928SBarry Smith Creates the private data structure for this preconditioner and 4074b9ad928SBarry Smith attach it to the PC object. 4084b9ad928SBarry Smith */ 40938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Jacobi,&jac);CHKERRQ(ierr); 4104b9ad928SBarry Smith pc->data = (void*)jac; 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith /* 4134b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 4144b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 4154b9ad928SBarry Smith */ 4164b9ad928SBarry Smith jac->diag = 0; 4174b9ad928SBarry Smith jac->diagsqrt = 0; 4184b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 41986697f06SMatthew Knepley jac->userowsum = PETSC_FALSE; 420cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith /* 4234b9ad928SBarry Smith Set the pointers for the functions that are provided above. 4244b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4254b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4264b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4274b9ad928SBarry Smith not needed. 4284b9ad928SBarry Smith */ 4294b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4304b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4314b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 432a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi; 4334b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4344b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4354b9ad928SBarry Smith pc->ops->view = 0; 4364b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4374b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4384b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 439*2fa5cd67SKarl Rupp 440cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 44186697f06SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 442cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 4434b9ad928SBarry Smith PetscFunctionReturn(0); 4444b9ad928SBarry Smith } 4454b9ad928SBarry Smith EXTERN_C_END 4464b9ad928SBarry Smith 447cd47f5d9SBarry Smith 448cd47f5d9SBarry Smith #undef __FUNCT__ 449cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs" 450cd47f5d9SBarry Smith /*@ 451cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 452cd47f5d9SBarry Smith absolute value of the diagonal to for the preconditioner 453cd47f5d9SBarry Smith 454ad4df100SBarry Smith Logically Collective on PC 455cd47f5d9SBarry Smith 456cd47f5d9SBarry Smith Input Parameters: 457cd47f5d9SBarry Smith . pc - the preconditioner context 458cd47f5d9SBarry Smith 459cd47f5d9SBarry Smith 460cd47f5d9SBarry Smith Options Database Key: 461cd47f5d9SBarry Smith . -pc_jacobi_abs 462cd47f5d9SBarry Smith 463cd47f5d9SBarry Smith Level: intermediate 464cd47f5d9SBarry Smith 465cd47f5d9SBarry Smith Concepts: Jacobi preconditioner 466cd47f5d9SBarry Smith 46786697f06SMatthew Knepley .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum() 468cd47f5d9SBarry Smith 469cd47f5d9SBarry Smith @*/ 4707087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseAbs(PC pc) 471cd47f5d9SBarry Smith { 4724ac538c5SBarry Smith PetscErrorCode ierr; 473cd47f5d9SBarry Smith 474cd47f5d9SBarry Smith PetscFunctionBegin; 4750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4764ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC),(pc));CHKERRQ(ierr); 477cd47f5d9SBarry Smith PetscFunctionReturn(0); 478cd47f5d9SBarry Smith } 479cd47f5d9SBarry Smith 4804b9ad928SBarry Smith #undef __FUNCT__ 4814b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4824b9ad928SBarry Smith /*@ 4834b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4844b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4854b9ad928SBarry Smith the diagonal entry 4864b9ad928SBarry Smith 487ad4df100SBarry Smith Logically Collective on PC 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith Input Parameters: 4904b9ad928SBarry Smith . pc - the preconditioner context 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith 4934b9ad928SBarry Smith Options Database Key: 4944b9ad928SBarry Smith . -pc_jacobi_rowmax 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith Level: intermediate 4974b9ad928SBarry Smith 4984b9ad928SBarry Smith Concepts: Jacobi preconditioner 4994b9ad928SBarry Smith 500cd47f5d9SBarry Smith .seealso: PCJacobiaUseAbs() 5014b9ad928SBarry Smith @*/ 5027087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowMax(PC pc) 5034b9ad928SBarry Smith { 5044ac538c5SBarry Smith PetscErrorCode ierr; 5054b9ad928SBarry Smith 5064b9ad928SBarry Smith PetscFunctionBegin; 5070700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5084ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowMax_C",(PC),(pc));CHKERRQ(ierr); 5094b9ad928SBarry Smith PetscFunctionReturn(0); 5104b9ad928SBarry Smith } 5114b9ad928SBarry Smith 51286697f06SMatthew Knepley #undef __FUNCT__ 51386697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum" 51486697f06SMatthew Knepley /*@ 51586697f06SMatthew Knepley PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 51686697f06SMatthew Knepley sum of each row as the diagonal preconditioner, instead of 51786697f06SMatthew Knepley the diagonal entry 51886697f06SMatthew Knepley 519ad4df100SBarry Smith Logical Collective on PC 52086697f06SMatthew Knepley 52186697f06SMatthew Knepley Input Parameters: 52286697f06SMatthew Knepley . pc - the preconditioner context 52386697f06SMatthew Knepley 52486697f06SMatthew Knepley 52586697f06SMatthew Knepley Options Database Key: 52686697f06SMatthew Knepley . -pc_jacobi_rowsum 52786697f06SMatthew Knepley 52886697f06SMatthew Knepley Level: intermediate 52986697f06SMatthew Knepley 53086697f06SMatthew Knepley Concepts: Jacobi preconditioner 53186697f06SMatthew Knepley 53286697f06SMatthew Knepley .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 53386697f06SMatthew Knepley @*/ 5347087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowSum(PC pc) 53586697f06SMatthew Knepley { 5364ac538c5SBarry Smith PetscErrorCode ierr; 53786697f06SMatthew Knepley 53886697f06SMatthew Knepley PetscFunctionBegin; 5390700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5404ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowSum_C",(PC),(pc));CHKERRQ(ierr); 54186697f06SMatthew Knepley PetscFunctionReturn(0); 54286697f06SMatthew Knepley } 54386697f06SMatthew Knepley 544