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 { 572fa5cd67SKarl 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 #undef __FUNCT__ 674b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 68f7a08781SBarry Smith static PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 694b9ad928SBarry Smith { 704b9ad928SBarry Smith PC_Jacobi *j; 714b9ad928SBarry Smith 724b9ad928SBarry Smith PetscFunctionBegin; 734b9ad928SBarry Smith j = (PC_Jacobi*)pc->data; 744b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 754b9ad928SBarry Smith PetscFunctionReturn(0); 764b9ad928SBarry Smith } 774b9ad928SBarry Smith 78cd47f5d9SBarry Smith #undef __FUNCT__ 7986697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 80f7a08781SBarry Smith static PetscErrorCode PCJacobiSetUseRowSum_Jacobi(PC pc) 8186697f06SMatthew Knepley { 8286697f06SMatthew Knepley PC_Jacobi *j; 8386697f06SMatthew Knepley 8486697f06SMatthew Knepley PetscFunctionBegin; 8586697f06SMatthew Knepley j = (PC_Jacobi*)pc->data; 8686697f06SMatthew Knepley j->userowsum = PETSC_TRUE; 8786697f06SMatthew Knepley PetscFunctionReturn(0); 8886697f06SMatthew Knepley } 8986697f06SMatthew Knepley 9086697f06SMatthew Knepley #undef __FUNCT__ 91cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 92f7a08781SBarry Smith static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc) 93cd47f5d9SBarry Smith { 94cd47f5d9SBarry Smith PC_Jacobi *j; 95cd47f5d9SBarry Smith 96cd47f5d9SBarry Smith PetscFunctionBegin; 97cd47f5d9SBarry Smith j = (PC_Jacobi*)pc->data; 98cd47f5d9SBarry Smith j->useabs = PETSC_TRUE; 99cd47f5d9SBarry Smith PetscFunctionReturn(0); 100cd47f5d9SBarry Smith } 101cd47f5d9SBarry Smith 1024b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1034b9ad928SBarry Smith /* 1044b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 1054b9ad928SBarry Smith by setting data structures and options. 1064b9ad928SBarry Smith 1074b9ad928SBarry Smith Input Parameter: 1084b9ad928SBarry Smith . pc - the preconditioner context 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith Notes: 1134b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1144b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1154b9ad928SBarry Smith */ 1164b9ad928SBarry Smith #undef __FUNCT__ 1174b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1186849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1194b9ad928SBarry Smith { 1204b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1214b9ad928SBarry Smith Vec diag,diagsqrt; 122dfbe8321SBarry Smith PetscErrorCode ierr; 123cd47f5d9SBarry Smith PetscInt n,i; 1244b9ad928SBarry Smith PetscScalar *x; 125ace3abfcSBarry Smith PetscBool zeroflag = PETSC_FALSE; 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith PetscFunctionBegin; 1284b9ad928SBarry Smith /* 1294b9ad928SBarry Smith For most preconditioners the code would begin here something like 1304b9ad928SBarry Smith 1314b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 13223ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1334b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1344b9ad928SBarry Smith } 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1374b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1384b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1394b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1404b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1414b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1424b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1434b9ad928SBarry Smith */ 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith /* 1464b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1474b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1484b9ad928SBarry Smith */ 1494b9ad928SBarry Smith diag = jac->diag; 1504b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith if (diag) { 1534b9ad928SBarry Smith if (jac->userowmax) { 1540298fd71SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr); 15586697f06SMatthew Knepley } else if (jac->userowsum) { 15686697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 1574b9ad928SBarry Smith } else { 1584b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1614b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1624b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 163cd47f5d9SBarry Smith if (jac->useabs) { 1642fa5cd67SKarl Rupp for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 165cd47f5d9SBarry Smith } 1664b9ad928SBarry Smith for (i=0; i<n; i++) { 1674b9ad928SBarry Smith if (x[i] == 0.0) { 1684b9ad928SBarry Smith x[i] = 1.0; 169cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1704b9ad928SBarry Smith } 1714b9ad928SBarry Smith } 1724b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1734b9ad928SBarry Smith } 1744b9ad928SBarry Smith if (diagsqrt) { 1754b9ad928SBarry Smith if (jac->userowmax) { 1760298fd71SBarry Smith ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr); 17786697f06SMatthew Knepley } else if (jac->userowsum) { 17886697f06SMatthew Knepley ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 1794b9ad928SBarry Smith } else { 1804b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1834b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1844b9ad928SBarry Smith for (i=0; i<n; i++) { 1858f1a2a5eSBarry Smith if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 1864b9ad928SBarry Smith else { 1874b9ad928SBarry Smith x[i] = 1.0; 188cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1894b9ad928SBarry Smith } 1904b9ad928SBarry Smith } 1914b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 1924b9ad928SBarry Smith } 1934b9ad928SBarry Smith if (zeroflag) { 194ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith PetscFunctionReturn(0); 1974b9ad928SBarry Smith } 1984b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1994b9ad928SBarry Smith /* 2004b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 2014b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 2024b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith Input Parameter: 2054b9ad928SBarry Smith . pc - the preconditioner context 2064b9ad928SBarry Smith */ 2074b9ad928SBarry Smith #undef __FUNCT__ 2084b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2096849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2104b9ad928SBarry Smith { 211dfbe8321SBarry Smith PetscErrorCode ierr; 2124b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith PetscFunctionBegin; 21523ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 21652e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 2174b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2184b9ad928SBarry Smith PetscFunctionReturn(0); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2214b9ad928SBarry Smith /* 2224b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2234b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2244b9ad928SBarry Smith right 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_NonSymmetric" 2316849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2324b9ad928SBarry Smith { 233dfbe8321SBarry Smith PetscErrorCode ierr; 2344b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith PetscFunctionBegin; 23723ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 23852e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 2394b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2404b9ad928SBarry Smith PetscFunctionReturn(0); 2414b9ad928SBarry Smith } 2424b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2434b9ad928SBarry Smith /* 2444b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith Input Parameters: 2474b9ad928SBarry Smith . pc - the preconditioner context 2484b9ad928SBarry Smith . x - input vector 2494b9ad928SBarry Smith 2504b9ad928SBarry Smith Output Parameter: 2514b9ad928SBarry Smith . y - output vector 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith Application Interface Routine: PCApply() 2544b9ad928SBarry Smith */ 2554b9ad928SBarry Smith #undef __FUNCT__ 2564b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2576849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2584b9ad928SBarry Smith { 2594b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 260dfbe8321SBarry Smith PetscErrorCode ierr; 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith PetscFunctionBegin; 2634b9ad928SBarry Smith if (!jac->diag) { 2644b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2654b9ad928SBarry Smith } 2662dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 2674b9ad928SBarry Smith PetscFunctionReturn(0); 2684b9ad928SBarry Smith } 2694b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2704b9ad928SBarry Smith /* 2714b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2724b9ad928SBarry Smith symmetric preconditioner to a vector. 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith Input Parameters: 2754b9ad928SBarry Smith . pc - the preconditioner context 2764b9ad928SBarry Smith . x - input vector 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith Output Parameter: 2794b9ad928SBarry Smith . y - output vector 2804b9ad928SBarry Smith 2814b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2824b9ad928SBarry Smith */ 2834b9ad928SBarry Smith #undef __FUNCT__ 2844b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 2856849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2864b9ad928SBarry Smith { 287dfbe8321SBarry Smith PetscErrorCode ierr; 2884b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2894b9ad928SBarry Smith 2904b9ad928SBarry Smith PetscFunctionBegin; 2914b9ad928SBarry Smith if (!jac->diagsqrt) { 2924b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 2934b9ad928SBarry Smith } 2942dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 298a06653b4SBarry Smith #undef __FUNCT__ 299a06653b4SBarry Smith #define __FUNCT__ "PCReset_Jacobi" 300a06653b4SBarry Smith static PetscErrorCode PCReset_Jacobi(PC pc) 301a06653b4SBarry Smith { 302a06653b4SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 303a06653b4SBarry Smith PetscErrorCode ierr; 304a06653b4SBarry Smith 305a06653b4SBarry Smith PetscFunctionBegin; 3066bf464f9SBarry Smith ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 3076bf464f9SBarry Smith ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 308a06653b4SBarry Smith PetscFunctionReturn(0); 309a06653b4SBarry Smith } 310a06653b4SBarry Smith 3114b9ad928SBarry Smith /* 3124b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 3134b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith Input Parameter: 3164b9ad928SBarry Smith . pc - the preconditioner context 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith Application Interface Routine: PCDestroy() 3194b9ad928SBarry Smith */ 3204b9ad928SBarry Smith #undef __FUNCT__ 3214b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3226849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3234b9ad928SBarry Smith { 324dfbe8321SBarry Smith PetscErrorCode ierr; 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith PetscFunctionBegin; 327a06653b4SBarry Smith ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith /* 3304b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3314b9ad928SBarry Smith */ 332c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 3334b9ad928SBarry Smith PetscFunctionReturn(0); 3344b9ad928SBarry Smith } 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith #undef __FUNCT__ 3374b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3386849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 3394b9ad928SBarry Smith { 3404b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 341dfbe8321SBarry Smith PetscErrorCode ierr; 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith PetscFunctionBegin; 3444b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 345acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3460298fd71SBarry Smith &jac->userowmax,NULL);CHKERRQ(ierr); 347acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 3480298fd71SBarry Smith &jac->userowsum,NULL);CHKERRQ(ierr); 349acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 3500298fd71SBarry Smith &jac->useabs,NULL);CHKERRQ(ierr); 3514b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3524b9ad928SBarry Smith PetscFunctionReturn(0); 3534b9ad928SBarry Smith } 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3564b9ad928SBarry Smith /* 3574b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3584b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3594b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith Input Parameter: 3624b9ad928SBarry Smith . pc - the preconditioner context 3634b9ad928SBarry Smith 3644b9ad928SBarry Smith Application Interface Routine: PCCreate() 3654b9ad928SBarry Smith */ 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith /*MC 3685a46d3faSBarry Smith PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Options Database Key: 371cd47f5d9SBarry Smith + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3724b9ad928SBarry Smith rather than the diagonal 373079a340bSHong Zhang . -pc_jacobi_rowsum - use the row sums as the scaling factor, 37486697f06SMatthew Knepley rather than the diagonal 375cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3764b9ad928SBarry Smith 3774b9ad928SBarry Smith Level: beginner 3784b9ad928SBarry Smith 3794b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3804b9ad928SBarry Smith 381b037da10SBarry Smith Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 3824b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 38786697f06SMatthew Knepley PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 3884b9ad928SBarry Smith M*/ 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith #undef __FUNCT__ 3914b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 392*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode 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; 425a06653b4SBarry Smith pc->ops->reset = PCReset_Jacobi; 4264b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4274b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4284b9ad928SBarry Smith pc->ops->view = 0; 4294b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4304b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4314b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 4322fa5cd67SKarl Rupp 43300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 43400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseRowSum_C","PCJacobiSetUseRowSum_Jacobi",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 43500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 4364b9ad928SBarry Smith PetscFunctionReturn(0); 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 445ad4df100SBarry Smith Logically 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 @*/ 4617087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseAbs(PC pc) 462cd47f5d9SBarry Smith { 4634ac538c5SBarry Smith PetscErrorCode ierr; 464cd47f5d9SBarry Smith 465cd47f5d9SBarry Smith PetscFunctionBegin; 4660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4674ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC),(pc));CHKERRQ(ierr); 468cd47f5d9SBarry Smith PetscFunctionReturn(0); 469cd47f5d9SBarry Smith } 470cd47f5d9SBarry Smith 4714b9ad928SBarry Smith #undef __FUNCT__ 4724b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4734b9ad928SBarry Smith /*@ 4744b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4754b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4764b9ad928SBarry Smith the diagonal entry 4774b9ad928SBarry Smith 478ad4df100SBarry Smith Logically Collective on PC 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith Input Parameters: 4814b9ad928SBarry Smith . pc - the preconditioner context 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith 4844b9ad928SBarry Smith Options Database Key: 4854b9ad928SBarry Smith . -pc_jacobi_rowmax 4864b9ad928SBarry Smith 4874b9ad928SBarry Smith Level: intermediate 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith Concepts: Jacobi preconditioner 4904b9ad928SBarry Smith 491cd47f5d9SBarry Smith .seealso: PCJacobiaUseAbs() 4924b9ad928SBarry Smith @*/ 4937087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowMax(PC pc) 4944b9ad928SBarry Smith { 4954ac538c5SBarry Smith PetscErrorCode ierr; 4964b9ad928SBarry Smith 4974b9ad928SBarry Smith PetscFunctionBegin; 4980700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4994ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowMax_C",(PC),(pc));CHKERRQ(ierr); 5004b9ad928SBarry Smith PetscFunctionReturn(0); 5014b9ad928SBarry Smith } 5024b9ad928SBarry Smith 50386697f06SMatthew Knepley #undef __FUNCT__ 50486697f06SMatthew Knepley #define __FUNCT__ "PCJacobiSetUseRowSum" 50586697f06SMatthew Knepley /*@ 50686697f06SMatthew Knepley PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 50786697f06SMatthew Knepley sum of each row as the diagonal preconditioner, instead of 50886697f06SMatthew Knepley the diagonal entry 50986697f06SMatthew Knepley 510ad4df100SBarry Smith Logical Collective on PC 51186697f06SMatthew Knepley 51286697f06SMatthew Knepley Input Parameters: 51386697f06SMatthew Knepley . pc - the preconditioner context 51486697f06SMatthew Knepley 51586697f06SMatthew Knepley 51686697f06SMatthew Knepley Options Database Key: 51786697f06SMatthew Knepley . -pc_jacobi_rowsum 51886697f06SMatthew Knepley 51986697f06SMatthew Knepley Level: intermediate 52086697f06SMatthew Knepley 52186697f06SMatthew Knepley Concepts: Jacobi preconditioner 52286697f06SMatthew Knepley 52386697f06SMatthew Knepley .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 52486697f06SMatthew Knepley @*/ 5257087cfbeSBarry Smith PetscErrorCode PCJacobiSetUseRowSum(PC pc) 52686697f06SMatthew Knepley { 5274ac538c5SBarry Smith PetscErrorCode ierr; 52886697f06SMatthew Knepley 52986697f06SMatthew Knepley PetscFunctionBegin; 5300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5314ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCJacobiSetUseRowSum_C",(PC),(pc));CHKERRQ(ierr); 53286697f06SMatthew Knepley PetscFunctionReturn(0); 53386697f06SMatthew Knepley } 53486697f06SMatthew Knepley 535