14b9ad928SBarry Smith 24b9ad928SBarry Smith /* -------------------------------------------------------------------- 34b9ad928SBarry Smith 44b9ad928SBarry Smith This file implements a Jacobi preconditioner for matrices that use 54b9ad928SBarry Smith the Mat interface (various matrix formats). Actually, the only 64b9ad928SBarry Smith matrix operation used here is MatGetDiagonal(), which extracts 74b9ad928SBarry Smith diagonal elements of the preconditioning matrix. 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 514b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 524b9ad928SBarry Smith 534b9ad928SBarry Smith /* 544b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 554b9ad928SBarry Smith */ 564b9ad928SBarry Smith typedef struct { 574b9ad928SBarry Smith Vec diag; /* vector containing the reciprocals of the diagonal elements 584b9ad928SBarry Smith of the preconditioner matrix */ 594b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 604b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 614b9ad928SBarry Smith only for symmetric preconditioner application) */ 624b9ad928SBarry Smith PetscTruth userowmax; 634b9ad928SBarry Smith } PC_Jacobi; 644b9ad928SBarry Smith 654b9ad928SBarry Smith EXTERN_C_BEGIN 664b9ad928SBarry Smith #undef __FUNCT__ 674b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 68dfbe8321SBarry Smith 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 EXTERN_C_END 784b9ad928SBarry Smith 794b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 804b9ad928SBarry Smith /* 814b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 824b9ad928SBarry Smith by setting data structures and options. 834b9ad928SBarry Smith 844b9ad928SBarry Smith Input Parameter: 854b9ad928SBarry Smith . pc - the preconditioner context 864b9ad928SBarry Smith 874b9ad928SBarry Smith Application Interface Routine: PCSetUp() 884b9ad928SBarry Smith 894b9ad928SBarry Smith Notes: 904b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 914b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 924b9ad928SBarry Smith */ 934b9ad928SBarry Smith #undef __FUNCT__ 944b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 95*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi(PC pc) 964b9ad928SBarry Smith { 974b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 984b9ad928SBarry Smith Vec diag,diagsqrt; 99dfbe8321SBarry Smith PetscErrorCode ierr; 100dfbe8321SBarry Smith int n,i,zeroflag=0; 1014b9ad928SBarry Smith PetscScalar *x; 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith PetscFunctionBegin; 1044b9ad928SBarry Smith 1054b9ad928SBarry Smith /* 1064b9ad928SBarry Smith For most preconditioners the code would begin here something like 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 10923ce1328SBarry Smith ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 1104b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 1114b9ad928SBarry Smith } 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 1144b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 1154b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 1164b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 1174b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 1184b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 1194b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 1204b9ad928SBarry Smith */ 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith /* 1234b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 1244b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 1254b9ad928SBarry Smith */ 1264b9ad928SBarry Smith diag = jac->diag; 1274b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith if (diag) { 1304b9ad928SBarry Smith if (jac->userowmax) { 1314b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diag);CHKERRQ(ierr); 1324b9ad928SBarry Smith } else { 1334b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 1344b9ad928SBarry Smith } 1354b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 1364b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 1374b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 1384b9ad928SBarry Smith for (i=0; i<n; i++) { 1394b9ad928SBarry Smith if (x[i] == 0.0) { 1404b9ad928SBarry Smith x[i] = 1.0; 1414b9ad928SBarry Smith zeroflag = 1; 1424b9ad928SBarry Smith } 1434b9ad928SBarry Smith } 1444b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 1454b9ad928SBarry Smith } 1464b9ad928SBarry Smith if (diagsqrt) { 1474b9ad928SBarry Smith if (jac->userowmax) { 1484b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diagsqrt);CHKERRQ(ierr); 1494b9ad928SBarry Smith } else { 1504b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 1514b9ad928SBarry Smith } 1524b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 1534b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 1544b9ad928SBarry Smith for (i=0; i<n; i++) { 1554b9ad928SBarry Smith if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 1564b9ad928SBarry Smith else { 1574b9ad928SBarry Smith x[i] = 1.0; 1584b9ad928SBarry Smith zeroflag = 1; 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith } 1614b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 1624b9ad928SBarry Smith } 1634b9ad928SBarry Smith if (zeroflag) { 1644b9ad928SBarry Smith PetscLogInfo(pc,"PCSetUp_Jacobi:Zero detected in diagonal of matrix, using 1 at those locations\n"); 1654b9ad928SBarry Smith } 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith PetscFunctionReturn(0); 1684b9ad928SBarry Smith } 1694b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1704b9ad928SBarry Smith /* 1714b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 1724b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 1734b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith Input Parameter: 1764b9ad928SBarry Smith . pc - the preconditioner context 1774b9ad928SBarry Smith */ 1784b9ad928SBarry Smith #undef __FUNCT__ 1794b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 180*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 1814b9ad928SBarry Smith { 182dfbe8321SBarry Smith PetscErrorCode ierr; 1834b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith PetscFunctionBegin; 1864b9ad928SBarry Smith 18723ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 1884b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diagsqrt); 1894b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 1904b9ad928SBarry Smith PetscFunctionReturn(0); 1914b9ad928SBarry Smith } 1924b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1934b9ad928SBarry Smith /* 1944b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 1954b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 1964b9ad928SBarry Smith right application of the Jacobi preconditioner. 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith Input Parameter: 1994b9ad928SBarry Smith . pc - the preconditioner context 2004b9ad928SBarry Smith */ 2014b9ad928SBarry Smith #undef __FUNCT__ 2024b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 203*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 2044b9ad928SBarry Smith { 205dfbe8321SBarry Smith PetscErrorCode ierr; 2064b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith PetscFunctionBegin; 2094b9ad928SBarry Smith 21023ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 2114b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 2124b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 2134b9ad928SBarry Smith PetscFunctionReturn(0); 2144b9ad928SBarry Smith } 2154b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2164b9ad928SBarry Smith /* 2174b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 2184b9ad928SBarry Smith 2194b9ad928SBarry Smith Input Parameters: 2204b9ad928SBarry Smith . pc - the preconditioner context 2214b9ad928SBarry Smith . x - input vector 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith Output Parameter: 2244b9ad928SBarry Smith . y - output vector 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Application Interface Routine: PCApply() 2274b9ad928SBarry Smith */ 2284b9ad928SBarry Smith #undef __FUNCT__ 2294b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 230*6849ba73SBarry Smith static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 2314b9ad928SBarry Smith { 2324b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 233dfbe8321SBarry Smith PetscErrorCode ierr; 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith PetscFunctionBegin; 2364b9ad928SBarry Smith if (!jac->diag) { 2374b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 2384b9ad928SBarry Smith } 2394b9ad928SBarry Smith ierr = VecPointwiseMult(x,jac->diag,y);CHKERRQ(ierr); 2404b9ad928SBarry Smith PetscFunctionReturn(0); 2414b9ad928SBarry Smith } 2424b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2434b9ad928SBarry Smith /* 2444b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 2454b9ad928SBarry Smith symmetric preconditioner to a vector. 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith Input Parameters: 2484b9ad928SBarry Smith . pc - the preconditioner context 2494b9ad928SBarry Smith . x - input vector 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith Output Parameter: 2524b9ad928SBarry Smith . y - output vector 2534b9ad928SBarry Smith 2544b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 2554b9ad928SBarry Smith */ 2564b9ad928SBarry Smith #undef __FUNCT__ 2574b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 258*6849ba73SBarry Smith static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 2594b9ad928SBarry Smith { 260dfbe8321SBarry Smith PetscErrorCode ierr; 2614b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 2624b9ad928SBarry Smith 2634b9ad928SBarry Smith PetscFunctionBegin; 2644b9ad928SBarry Smith if (!jac->diagsqrt) { 2654b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 2664b9ad928SBarry Smith } 2674b9ad928SBarry Smith VecPointwiseMult(x,jac->diagsqrt,y); 2684b9ad928SBarry Smith PetscFunctionReturn(0); 2694b9ad928SBarry Smith } 2704b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2714b9ad928SBarry Smith /* 2724b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 2734b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 2744b9ad928SBarry Smith 2754b9ad928SBarry Smith Input Parameter: 2764b9ad928SBarry Smith . pc - the preconditioner context 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith Application Interface Routine: PCDestroy() 2794b9ad928SBarry Smith */ 2804b9ad928SBarry Smith #undef __FUNCT__ 2814b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 282*6849ba73SBarry Smith static PetscErrorCode PCDestroy_Jacobi(PC pc) 2834b9ad928SBarry Smith { 2844b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 285dfbe8321SBarry Smith PetscErrorCode ierr; 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith PetscFunctionBegin; 2884b9ad928SBarry Smith if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 2894b9ad928SBarry Smith if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 2904b9ad928SBarry Smith 2914b9ad928SBarry Smith /* 2924b9ad928SBarry Smith Free the private data structure that was hanging off the PC 2934b9ad928SBarry Smith */ 2944b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 2984b9ad928SBarry Smith #undef __FUNCT__ 2994b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 300*6849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_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 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 3074b9ad928SBarry Smith ierr = PetscOptionsLogical("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 3084b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 3094b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3104b9ad928SBarry Smith PetscFunctionReturn(0); 3114b9ad928SBarry Smith } 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 3144b9ad928SBarry Smith /* 3154b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 3164b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 3174b9ad928SBarry Smith context, PC, that was created within PCCreate(). 3184b9ad928SBarry Smith 3194b9ad928SBarry Smith Input Parameter: 3204b9ad928SBarry Smith . pc - the preconditioner context 3214b9ad928SBarry Smith 3224b9ad928SBarry Smith Application Interface Routine: PCCreate() 3234b9ad928SBarry Smith */ 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith /*MC 3264b9ad928SBarry Smith PCJacobi - Jacobi (i.e. diagonal scaling preconditioning) 3274b9ad928SBarry Smith 3284b9ad928SBarry Smith Options Database Key: 3294b9ad928SBarry Smith . -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 3304b9ad928SBarry Smith rather than the diagonal 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Level: beginner 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 3374b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3424b9ad928SBarry Smith PCJacobiSetUseRowMax(), 3434b9ad928SBarry Smith M*/ 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith EXTERN_C_BEGIN 3464b9ad928SBarry Smith #undef __FUNCT__ 3474b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 348dfbe8321SBarry Smith PetscErrorCode PCCreate_Jacobi(PC pc) 3494b9ad928SBarry Smith { 3504b9ad928SBarry Smith PC_Jacobi *jac; 351dfbe8321SBarry Smith PetscErrorCode ierr; 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith PetscFunctionBegin; 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith /* 3564b9ad928SBarry Smith Creates the private data structure for this preconditioner and 3574b9ad928SBarry Smith attach it to the PC object. 3584b9ad928SBarry Smith */ 3594b9ad928SBarry Smith ierr = PetscNew(PC_Jacobi,&jac);CHKERRQ(ierr); 3604b9ad928SBarry Smith pc->data = (void*)jac; 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith /* 3634b9ad928SBarry Smith Logs the memory usage; this is not needed but allows PETSc to 3644b9ad928SBarry Smith monitor how much memory is being used for various purposes. 3654b9ad928SBarry Smith */ 3664b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Jacobi)); 3674b9ad928SBarry Smith 3684b9ad928SBarry Smith /* 3694b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 3704b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 3714b9ad928SBarry Smith */ 3724b9ad928SBarry Smith jac->diag = 0; 3734b9ad928SBarry Smith jac->diagsqrt = 0; 3744b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith /* 3774b9ad928SBarry Smith Set the pointers for the functions that are provided above. 3784b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 3794b9ad928SBarry Smith are called, they will automatically call these functions. Note we 3804b9ad928SBarry Smith choose not to provide a couple of these functions since they are 3814b9ad928SBarry Smith not needed. 3824b9ad928SBarry Smith */ 3834b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 3844b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 3854b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 3864b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 3874b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 3884b9ad928SBarry Smith pc->ops->view = 0; 3894b9ad928SBarry Smith pc->ops->applyrichardson = 0; 3904b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 3914b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 3924b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi", 3934b9ad928SBarry Smith PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 3944b9ad928SBarry Smith PetscFunctionReturn(0); 3954b9ad928SBarry Smith } 3964b9ad928SBarry Smith EXTERN_C_END 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith #undef __FUNCT__ 3994b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 4004b9ad928SBarry Smith /*@ 4014b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 4024b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 4034b9ad928SBarry Smith the diagonal entry 4044b9ad928SBarry Smith 4054b9ad928SBarry Smith Collective on PC 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith Input Parameters: 4084b9ad928SBarry Smith . pc - the preconditioner context 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith 4114b9ad928SBarry Smith Options Database Key: 4124b9ad928SBarry Smith . -pc_jacobi_rowmax 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith Level: intermediate 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith Concepts: Jacobi preconditioner 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith @*/ 419dfbe8321SBarry Smith PetscErrorCode PCJacobiSetUseRowMax(PC pc) 4204b9ad928SBarry Smith { 421dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith PetscFunctionBegin; 4244482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4254b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 4264b9ad928SBarry Smith if (f) { 4274b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 4284b9ad928SBarry Smith } 4294b9ad928SBarry Smith PetscFunctionReturn(0); 4304b9ad928SBarry Smith } 4314b9ad928SBarry Smith 432