1dba47a55SKris Buschelman #define PETSCKSP_DLL 24b9ad928SBarry Smith 34b9ad928SBarry Smith /* -------------------------------------------------------------------- 44b9ad928SBarry Smith 54b9ad928SBarry Smith This file implements a Jacobi preconditioner for matrices that use 64b9ad928SBarry Smith the Mat interface (various matrix formats). Actually, the only 74b9ad928SBarry Smith matrix operation used here is MatGetDiagonal(), which extracts 84b9ad928SBarry Smith diagonal elements of the preconditioning matrix. 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; 64*cd47f5d9SBarry Smith PetscTruth 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" 70dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT 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 81*cd47f5d9SBarry Smith EXTERN_C_BEGIN 82*cd47f5d9SBarry Smith #undef __FUNCT__ 83*cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 84*cd47f5d9SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs_Jacobi(PC pc) 85*cd47f5d9SBarry Smith { 86*cd47f5d9SBarry Smith PC_Jacobi *j; 87*cd47f5d9SBarry Smith 88*cd47f5d9SBarry Smith PetscFunctionBegin; 89*cd47f5d9SBarry Smith j = (PC_Jacobi*)pc->data; 90*cd47f5d9SBarry Smith j->useabs = PETSC_TRUE; 91*cd47f5d9SBarry Smith PetscFunctionReturn(0); 92*cd47f5d9SBarry Smith } 93*cd47f5d9SBarry Smith EXTERN_C_END 94*cd47f5d9SBarry Smith 954b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 964b9ad928SBarry Smith /* 974b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 984b9ad928SBarry Smith by setting data structures and options. 994b9ad928SBarry Smith 1004b9ad928SBarry Smith Input Parameter: 1014b9ad928SBarry Smith . pc - the preconditioner context 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith Application Interface Routine: PCSetUp() 1044b9ad928SBarry Smith 1054b9ad928SBarry Smith Notes: 1064b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 1074b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 1084b9ad928SBarry Smith */ 1094b9ad928SBarry Smith #undef __FUNCT__ 1104b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 1116849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 1124b9ad928SBarry Smith { 1134b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1144b9ad928SBarry Smith Vec diag,diagsqrt; 115dfbe8321SBarry Smith PetscErrorCode ierr; 116*cd47f5d9SBarry Smith PetscInt n,i; 1174b9ad928SBarry Smith PetscScalar *x; 118*cd47f5d9SBarry Smith PetscTruth zeroflag = PETSC_FALSE; 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith PetscFunctionBegin; 1214b9ad928SBarry Smith /* 1224b9ad928SBarry Smith For most preconditioners the code would begin here something like 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 12523ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1264b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1274b9ad928SBarry Smith } 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1304b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1314b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1324b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1334b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1344b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1354b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1364b9ad928SBarry Smith */ 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith /* 1394b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1404b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1414b9ad928SBarry Smith */ 1424b9ad928SBarry Smith diag = jac->diag; 1434b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith if (diag) { 1464b9ad928SBarry Smith if (jac->userowmax) { 1474b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diag);CHKERRQ(ierr); 1484b9ad928SBarry Smith } else { 1494b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1504b9ad928SBarry Smith } 1514b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1524b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1534b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 154*cd47f5d9SBarry Smith if (jac->useabs) { 155*cd47f5d9SBarry Smith for (i=0; i<n; i++) { 156*cd47f5d9SBarry Smith x[i] = PetscAbsScalar(x[i]); 157*cd47f5d9SBarry Smith } 158*cd47f5d9SBarry Smith } 1594b9ad928SBarry Smith for (i=0; i<n; i++) { 1604b9ad928SBarry Smith if (x[i] == 0.0) { 1614b9ad928SBarry Smith x[i] = 1.0; 162*cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1634b9ad928SBarry Smith } 1644b9ad928SBarry Smith } 1654b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1664b9ad928SBarry Smith } 1674b9ad928SBarry Smith if (diagsqrt) { 1684b9ad928SBarry Smith if (jac->userowmax) { 1694b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diagsqrt);CHKERRQ(ierr); 1704b9ad928SBarry Smith } else { 1714b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1724b9ad928SBarry Smith } 1734b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1744b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1754b9ad928SBarry Smith for (i=0; i<n; i++) { 1764b9ad928SBarry Smith if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 1774b9ad928SBarry Smith else { 1784b9ad928SBarry Smith x[i] = 1.0; 179*cd47f5d9SBarry Smith zeroflag = PETSC_TRUE; 1804b9ad928SBarry Smith } 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith if (zeroflag) { 185ae15b995SBarry Smith ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 1864b9ad928SBarry Smith } 1874b9ad928SBarry Smith PetscFunctionReturn(0); 1884b9ad928SBarry Smith } 1894b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1904b9ad928SBarry Smith /* 1914b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 1924b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 1934b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith Input Parameter: 1964b9ad928SBarry Smith . pc - the preconditioner context 1974b9ad928SBarry Smith */ 1984b9ad928SBarry Smith #undef __FUNCT__ 1994b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 2006849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 2014b9ad928SBarry Smith { 202dfbe8321SBarry Smith PetscErrorCode ierr; 2034b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2044b9ad928SBarry Smith 2054b9ad928SBarry Smith PetscFunctionBegin; 20623ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 20752e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diagsqrt);CHKERRQ(ierr); 2084b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2094b9ad928SBarry Smith PetscFunctionReturn(0); 2104b9ad928SBarry Smith } 2114b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2124b9ad928SBarry Smith /* 2134b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 2144b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 2154b9ad928SBarry Smith right application of the Jacobi preconditioner. 2164b9ad928SBarry Smith 2174b9ad928SBarry Smith Input Parameter: 2184b9ad928SBarry Smith . pc - the preconditioner context 2194b9ad928SBarry Smith */ 2204b9ad928SBarry Smith #undef __FUNCT__ 2214b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 2226849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2234b9ad928SBarry Smith { 224dfbe8321SBarry Smith PetscErrorCode ierr; 2254b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2264b9ad928SBarry Smith 2274b9ad928SBarry Smith PetscFunctionBegin; 22823ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 22952e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,jac->diag);CHKERRQ(ierr); 2304b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2314b9ad928SBarry Smith PetscFunctionReturn(0); 2324b9ad928SBarry Smith } 2334b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2344b9ad928SBarry Smith /* 2354b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith Input Parameters: 2384b9ad928SBarry Smith . pc - the preconditioner context 2394b9ad928SBarry Smith . x - input vector 2404b9ad928SBarry Smith 2414b9ad928SBarry Smith Output Parameter: 2424b9ad928SBarry Smith . y - output vector 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith Application Interface Routine: PCApply() 2454b9ad928SBarry Smith */ 2464b9ad928SBarry Smith #undef __FUNCT__ 2474b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 2486849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2494b9ad928SBarry Smith { 2504b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 251dfbe8321SBarry Smith PetscErrorCode ierr; 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith PetscFunctionBegin; 2544b9ad928SBarry Smith if (!jac->diag) { 2554b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2564b9ad928SBarry Smith } 2572dcb1b2aSMatthew Knepley ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 2584b9ad928SBarry Smith PetscFunctionReturn(0); 2594b9ad928SBarry Smith } 2604b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2614b9ad928SBarry Smith /* 2624b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2634b9ad928SBarry Smith symmetric preconditioner to a vector. 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith Input Parameters: 2664b9ad928SBarry Smith . pc - the preconditioner context 2674b9ad928SBarry Smith . x - input vector 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith Output Parameter: 2704b9ad928SBarry Smith . y - output vector 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2734b9ad928SBarry Smith */ 2744b9ad928SBarry Smith #undef __FUNCT__ 2754b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 2766849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2774b9ad928SBarry Smith { 278dfbe8321SBarry Smith PetscErrorCode ierr; 2794b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2804b9ad928SBarry Smith 2814b9ad928SBarry Smith PetscFunctionBegin; 2824b9ad928SBarry Smith if (!jac->diagsqrt) { 2834b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 2844b9ad928SBarry Smith } 2852dcb1b2aSMatthew Knepley VecPointwiseMult(y,x,jac->diagsqrt); 2864b9ad928SBarry Smith PetscFunctionReturn(0); 2874b9ad928SBarry Smith } 2884b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2894b9ad928SBarry Smith /* 2904b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 2914b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 2924b9ad928SBarry Smith 2934b9ad928SBarry Smith Input Parameter: 2944b9ad928SBarry Smith . pc - the preconditioner context 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith Application Interface Routine: PCDestroy() 2974b9ad928SBarry Smith */ 2984b9ad928SBarry Smith #undef __FUNCT__ 2994b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 3006849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 3014b9ad928SBarry Smith { 3024b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 303dfbe8321SBarry Smith PetscErrorCode ierr; 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith PetscFunctionBegin; 3064b9ad928SBarry Smith if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 3074b9ad928SBarry Smith if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 3084b9ad928SBarry Smith 3094b9ad928SBarry Smith /* 3104b9ad928SBarry Smith Free the private data structure that was hanging off the PC 3114b9ad928SBarry Smith */ 3124b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 3134b9ad928SBarry Smith PetscFunctionReturn(0); 3144b9ad928SBarry Smith } 3154b9ad928SBarry Smith 3164b9ad928SBarry Smith #undef __FUNCT__ 3174b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 3186849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_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 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 3254dc4c822SBarry Smith ierr = PetscOptionsTruth("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3264b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 327*cd47f5d9SBarry Smith ierr = PetscOptionsTruth("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 328*cd47f5d9SBarry Smith &jac->useabs,PETSC_NULL);CHKERRQ(ierr); 3294b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3304b9ad928SBarry Smith PetscFunctionReturn(0); 3314b9ad928SBarry Smith } 3324b9ad928SBarry Smith 3334b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3344b9ad928SBarry Smith /* 3354b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3364b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3374b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Input Parameter: 3404b9ad928SBarry Smith . pc - the preconditioner context 3414b9ad928SBarry Smith 3424b9ad928SBarry Smith Application Interface Routine: PCCreate() 3434b9ad928SBarry Smith */ 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith /*MC 3464b9ad928SBarry Smith PCJacobi - Jacobi (i.e. diagonal scaling preconditioning) 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith Options Database Key: 349*cd47f5d9SBarry Smith + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3504b9ad928SBarry Smith rather than the diagonal 351*cd47f5d9SBarry Smith - -pc_jacobi_abs - use the absolute value of the diagaonl entry 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith Level: beginner 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3564b9ad928SBarry Smith 3574b9ad928SBarry Smith Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 3584b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3594b9ad928SBarry Smith 3604b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 363*cd47f5d9SBarry Smith PCJacobiSetUseRowMax(), PCJacobiSetUseAbs() 3644b9ad928SBarry Smith M*/ 3654b9ad928SBarry Smith 3664b9ad928SBarry Smith EXTERN_C_BEGIN 3674b9ad928SBarry Smith #undef __FUNCT__ 3684b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 369dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Jacobi(PC pc) 3704b9ad928SBarry Smith { 3714b9ad928SBarry Smith PC_Jacobi *jac; 372dfbe8321SBarry Smith PetscErrorCode ierr; 3734b9ad928SBarry Smith 3744b9ad928SBarry Smith PetscFunctionBegin; 3754b9ad928SBarry Smith /* 3764b9ad928SBarry Smith Creates the private data structure for this preconditioner and 3774b9ad928SBarry Smith attach it to the PC object. 3784b9ad928SBarry Smith */ 3794b9ad928SBarry Smith ierr = PetscNew(PC_Jacobi,&jac);CHKERRQ(ierr); 3804b9ad928SBarry Smith pc->data = (void*)jac; 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith /* 3834b9ad928SBarry Smith Logs the memory usage; this is not needed but allows PETSc to 3844b9ad928SBarry Smith monitor how much memory is being used for various purposes. 3854b9ad928SBarry Smith */ 38652e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_Jacobi));CHKERRQ(ierr); 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith /* 3894b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 3904b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 3914b9ad928SBarry Smith */ 3924b9ad928SBarry Smith jac->diag = 0; 3934b9ad928SBarry Smith jac->diagsqrt = 0; 3944b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 395*cd47f5d9SBarry Smith jac->useabs = PETSC_FALSE; 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith /* 3984b9ad928SBarry Smith Set the pointers for the functions that are provided above. 3994b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 4004b9ad928SBarry Smith are called, they will automatically call these functions. Note we 4014b9ad928SBarry Smith choose not to provide a couple of these functions since they are 4024b9ad928SBarry Smith not needed. 4034b9ad928SBarry Smith */ 4044b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 4054b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 4064b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 4074b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 4084b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 4094b9ad928SBarry Smith pc->ops->view = 0; 4104b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4114b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 4124b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 413*cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 414*cd47f5d9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseAbs_C","PCJacobiSetUseAbs_Jacobi",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 4154b9ad928SBarry Smith PetscFunctionReturn(0); 4164b9ad928SBarry Smith } 4174b9ad928SBarry Smith EXTERN_C_END 4184b9ad928SBarry Smith 419*cd47f5d9SBarry Smith 420*cd47f5d9SBarry Smith #undef __FUNCT__ 421*cd47f5d9SBarry Smith #define __FUNCT__ "PCJacobiSetUseAbs" 422*cd47f5d9SBarry Smith /*@ 423*cd47f5d9SBarry Smith PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 424*cd47f5d9SBarry Smith absolute value of the diagonal to for the preconditioner 425*cd47f5d9SBarry Smith 426*cd47f5d9SBarry Smith Collective on PC 427*cd47f5d9SBarry Smith 428*cd47f5d9SBarry Smith Input Parameters: 429*cd47f5d9SBarry Smith . pc - the preconditioner context 430*cd47f5d9SBarry Smith 431*cd47f5d9SBarry Smith 432*cd47f5d9SBarry Smith Options Database Key: 433*cd47f5d9SBarry Smith . -pc_jacobi_abs 434*cd47f5d9SBarry Smith 435*cd47f5d9SBarry Smith Level: intermediate 436*cd47f5d9SBarry Smith 437*cd47f5d9SBarry Smith Concepts: Jacobi preconditioner 438*cd47f5d9SBarry Smith 439*cd47f5d9SBarry Smith .seealso: PCJacobiaUseRowMax() 440*cd47f5d9SBarry Smith 441*cd47f5d9SBarry Smith @*/ 442*cd47f5d9SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC pc) 443*cd47f5d9SBarry Smith { 444*cd47f5d9SBarry Smith PetscErrorCode ierr,(*f)(PC); 445*cd47f5d9SBarry Smith 446*cd47f5d9SBarry Smith PetscFunctionBegin; 447*cd47f5d9SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 448*cd47f5d9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",(void (**)(void))&f);CHKERRQ(ierr); 449*cd47f5d9SBarry Smith if (f) { 450*cd47f5d9SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 451*cd47f5d9SBarry Smith } 452*cd47f5d9SBarry Smith PetscFunctionReturn(0); 453*cd47f5d9SBarry Smith } 454*cd47f5d9SBarry Smith 4554b9ad928SBarry Smith #undef __FUNCT__ 4564b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4574b9ad928SBarry Smith /*@ 4584b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4594b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4604b9ad928SBarry Smith the diagonal entry 4614b9ad928SBarry Smith 4624b9ad928SBarry Smith Collective on PC 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith Input Parameters: 4654b9ad928SBarry Smith . pc - the preconditioner context 4664b9ad928SBarry Smith 4674b9ad928SBarry Smith 4684b9ad928SBarry Smith Options Database Key: 4694b9ad928SBarry Smith . -pc_jacobi_rowmax 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith Level: intermediate 4724b9ad928SBarry Smith 4734b9ad928SBarry Smith Concepts: Jacobi preconditioner 4744b9ad928SBarry Smith 475*cd47f5d9SBarry Smith .seealso: PCJacobiaUseAbs() 4764b9ad928SBarry Smith @*/ 477dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC pc) 4784b9ad928SBarry Smith { 479dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 4804b9ad928SBarry Smith 4814b9ad928SBarry Smith PetscFunctionBegin; 4824482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 483e34fafa9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetUseRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 4844b9ad928SBarry Smith if (f) { 4854b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 4864b9ad928SBarry Smith } 4874b9ad928SBarry Smith PetscFunctionReturn(0); 4884b9ad928SBarry Smith } 4894b9ad928SBarry Smith 490