1*4b9ad928SBarry Smith /*$Id: jacobi.c,v 1.75 2001/08/07 03:03:32 balay Exp $*/ 2*4b9ad928SBarry Smith 3*4b9ad928SBarry Smith /* -------------------------------------------------------------------- 4*4b9ad928SBarry Smith 5*4b9ad928SBarry Smith This file implements a Jacobi preconditioner for matrices that use 6*4b9ad928SBarry Smith the Mat interface (various matrix formats). Actually, the only 7*4b9ad928SBarry Smith matrix operation used here is MatGetDiagonal(), which extracts 8*4b9ad928SBarry Smith diagonal elements of the preconditioning matrix. 9*4b9ad928SBarry Smith 10*4b9ad928SBarry Smith The following basic routines are required for each preconditioner. 11*4b9ad928SBarry Smith PCCreate_XXX() - Creates a preconditioner context 12*4b9ad928SBarry Smith PCSetFromOptions_XXX() - Sets runtime options 13*4b9ad928SBarry Smith PCApply_XXX() - Applies the preconditioner 14*4b9ad928SBarry Smith PCDestroy_XXX() - Destroys the preconditioner context 15*4b9ad928SBarry Smith where the suffix "_XXX" denotes a particular implementation, in 16*4b9ad928SBarry Smith this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi). 17*4b9ad928SBarry Smith These routines are actually called via the common user interface 18*4b9ad928SBarry Smith routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 19*4b9ad928SBarry Smith so the application code interface remains identical for all 20*4b9ad928SBarry Smith preconditioners. 21*4b9ad928SBarry Smith 22*4b9ad928SBarry Smith Another key routine is: 23*4b9ad928SBarry Smith PCSetUp_XXX() - Prepares for the use of a preconditioner 24*4b9ad928SBarry Smith by setting data structures and options. The interface routine PCSetUp() 25*4b9ad928SBarry Smith is not usually called directly by the user, but instead is called by 26*4b9ad928SBarry Smith PCApply() if necessary. 27*4b9ad928SBarry Smith 28*4b9ad928SBarry Smith Additional basic routines are: 29*4b9ad928SBarry Smith PCView_XXX() - Prints details of runtime options that 30*4b9ad928SBarry Smith have actually been used. 31*4b9ad928SBarry Smith These are called by application codes via the interface routines 32*4b9ad928SBarry Smith PCView(). 33*4b9ad928SBarry Smith 34*4b9ad928SBarry Smith The various types of solvers (preconditioners, Krylov subspace methods, 35*4b9ad928SBarry Smith nonlinear solvers, timesteppers) are all organized similarly, so the 36*4b9ad928SBarry Smith above description applies to these categories also. One exception is 37*4b9ad928SBarry Smith that the analogues of PCApply() for these components are KSPSolve(), 38*4b9ad928SBarry Smith SNESSolve(), and TSSolve(). 39*4b9ad928SBarry Smith 40*4b9ad928SBarry Smith Additional optional functionality unique to preconditioners is left and 41*4b9ad928SBarry Smith right symmetric preconditioner application via PCApplySymmetricLeft() 42*4b9ad928SBarry Smith and PCApplySymmetricRight(). The Jacobi implementation is 43*4b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi(). 44*4b9ad928SBarry Smith 45*4b9ad928SBarry Smith -------------------------------------------------------------------- */ 46*4b9ad928SBarry Smith 47*4b9ad928SBarry Smith /* 48*4b9ad928SBarry Smith Include files needed for the Jacobi preconditioner: 49*4b9ad928SBarry Smith pcimpl.h - private include file intended for use by all preconditioners 50*4b9ad928SBarry Smith */ 51*4b9ad928SBarry Smith 52*4b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 53*4b9ad928SBarry Smith 54*4b9ad928SBarry Smith /* 55*4b9ad928SBarry Smith Private context (data structure) for the Jacobi preconditioner. 56*4b9ad928SBarry Smith */ 57*4b9ad928SBarry Smith typedef struct { 58*4b9ad928SBarry Smith Vec diag; /* vector containing the reciprocals of the diagonal elements 59*4b9ad928SBarry Smith of the preconditioner matrix */ 60*4b9ad928SBarry Smith Vec diagsqrt; /* vector containing the reciprocals of the square roots of 61*4b9ad928SBarry Smith the diagonal elements of the preconditioner matrix (used 62*4b9ad928SBarry Smith only for symmetric preconditioner application) */ 63*4b9ad928SBarry Smith PetscTruth userowmax; 64*4b9ad928SBarry Smith } PC_Jacobi; 65*4b9ad928SBarry Smith 66*4b9ad928SBarry Smith EXTERN_C_BEGIN 67*4b9ad928SBarry Smith #undef __FUNCT__ 68*4b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 69*4b9ad928SBarry Smith int PCJacobiSetUseRowMax_Jacobi(PC pc) 70*4b9ad928SBarry Smith { 71*4b9ad928SBarry Smith PC_Jacobi *j; 72*4b9ad928SBarry Smith 73*4b9ad928SBarry Smith PetscFunctionBegin; 74*4b9ad928SBarry Smith j = (PC_Jacobi*)pc->data; 75*4b9ad928SBarry Smith j->userowmax = PETSC_TRUE; 76*4b9ad928SBarry Smith PetscFunctionReturn(0); 77*4b9ad928SBarry Smith } 78*4b9ad928SBarry Smith EXTERN_C_END 79*4b9ad928SBarry Smith 80*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 81*4b9ad928SBarry Smith /* 82*4b9ad928SBarry Smith PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 83*4b9ad928SBarry Smith by setting data structures and options. 84*4b9ad928SBarry Smith 85*4b9ad928SBarry Smith Input Parameter: 86*4b9ad928SBarry Smith . pc - the preconditioner context 87*4b9ad928SBarry Smith 88*4b9ad928SBarry Smith Application Interface Routine: PCSetUp() 89*4b9ad928SBarry Smith 90*4b9ad928SBarry Smith Notes: 91*4b9ad928SBarry Smith The interface routine PCSetUp() is not usually called directly by 92*4b9ad928SBarry Smith the user, but instead is called by PCApply() if necessary. 93*4b9ad928SBarry Smith */ 94*4b9ad928SBarry Smith #undef __FUNCT__ 95*4b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi" 96*4b9ad928SBarry Smith static int PCSetUp_Jacobi(PC pc) 97*4b9ad928SBarry Smith { 98*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 99*4b9ad928SBarry Smith Vec diag,diagsqrt; 100*4b9ad928SBarry Smith int ierr,n,i,zeroflag=0; 101*4b9ad928SBarry Smith PetscScalar *x; 102*4b9ad928SBarry Smith 103*4b9ad928SBarry Smith PetscFunctionBegin; 104*4b9ad928SBarry Smith 105*4b9ad928SBarry Smith /* 106*4b9ad928SBarry Smith For most preconditioners the code would begin here something like 107*4b9ad928SBarry Smith 108*4b9ad928SBarry Smith if (pc->setupcalled == 0) { allocate space the first time this is ever called 109*4b9ad928SBarry Smith ierr = VecDuplicate(pc->vec,&jac->diag);CHKERRQ(ierr); 110*4b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 111*4b9ad928SBarry Smith } 112*4b9ad928SBarry Smith 113*4b9ad928SBarry Smith But for this preconditioner we want to support use of both the matrix' diagonal 114*4b9ad928SBarry Smith elements (for left or right preconditioning) and square root of diagonal elements 115*4b9ad928SBarry Smith (for symmetric preconditioning). Hence we do not allocate space here, since we 116*4b9ad928SBarry Smith don't know at this point which will be needed (diag and/or diagsqrt) until the user 117*4b9ad928SBarry Smith applies the preconditioner, and we don't want to allocate BOTH unless we need 118*4b9ad928SBarry Smith them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 119*4b9ad928SBarry Smith and PCSetUp_Jacobi_Symmetric(), respectively. 120*4b9ad928SBarry Smith */ 121*4b9ad928SBarry Smith 122*4b9ad928SBarry Smith /* 123*4b9ad928SBarry Smith Here we set up the preconditioner; that is, we copy the diagonal values from 124*4b9ad928SBarry Smith the matrix and put them into a format to make them quick to apply as a preconditioner. 125*4b9ad928SBarry Smith */ 126*4b9ad928SBarry Smith diag = jac->diag; 127*4b9ad928SBarry Smith diagsqrt = jac->diagsqrt; 128*4b9ad928SBarry Smith 129*4b9ad928SBarry Smith if (diag) { 130*4b9ad928SBarry Smith if (jac->userowmax) { 131*4b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diag);CHKERRQ(ierr); 132*4b9ad928SBarry Smith } else { 133*4b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 134*4b9ad928SBarry Smith } 135*4b9ad928SBarry Smith ierr = VecReciprocal(diag);CHKERRQ(ierr); 136*4b9ad928SBarry Smith ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 137*4b9ad928SBarry Smith ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 138*4b9ad928SBarry Smith for (i=0; i<n; i++) { 139*4b9ad928SBarry Smith if (x[i] == 0.0) { 140*4b9ad928SBarry Smith x[i] = 1.0; 141*4b9ad928SBarry Smith zeroflag = 1; 142*4b9ad928SBarry Smith } 143*4b9ad928SBarry Smith } 144*4b9ad928SBarry Smith ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 145*4b9ad928SBarry Smith } 146*4b9ad928SBarry Smith if (diagsqrt) { 147*4b9ad928SBarry Smith if (jac->userowmax) { 148*4b9ad928SBarry Smith ierr = MatGetRowMax(pc->pmat,diagsqrt);CHKERRQ(ierr); 149*4b9ad928SBarry Smith } else { 150*4b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 151*4b9ad928SBarry Smith } 152*4b9ad928SBarry Smith ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 153*4b9ad928SBarry Smith ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 154*4b9ad928SBarry Smith for (i=0; i<n; i++) { 155*4b9ad928SBarry Smith if (x[i] != 0.0) x[i] = 1.0/sqrt(PetscAbsScalar(x[i])); 156*4b9ad928SBarry Smith else { 157*4b9ad928SBarry Smith x[i] = 1.0; 158*4b9ad928SBarry Smith zeroflag = 1; 159*4b9ad928SBarry Smith } 160*4b9ad928SBarry Smith } 161*4b9ad928SBarry Smith ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 162*4b9ad928SBarry Smith } 163*4b9ad928SBarry Smith if (zeroflag) { 164*4b9ad928SBarry Smith PetscLogInfo(pc,"PCSetUp_Jacobi:Zero detected in diagonal of matrix, using 1 at those locations\n"); 165*4b9ad928SBarry Smith } 166*4b9ad928SBarry Smith 167*4b9ad928SBarry Smith PetscFunctionReturn(0); 168*4b9ad928SBarry Smith } 169*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 170*4b9ad928SBarry Smith /* 171*4b9ad928SBarry Smith PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 172*4b9ad928SBarry Smith inverse of the square root of the diagonal entries of the matrix. This 173*4b9ad928SBarry Smith is used for symmetric application of the Jacobi preconditioner. 174*4b9ad928SBarry Smith 175*4b9ad928SBarry Smith Input Parameter: 176*4b9ad928SBarry Smith . pc - the preconditioner context 177*4b9ad928SBarry Smith */ 178*4b9ad928SBarry Smith #undef __FUNCT__ 179*4b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 180*4b9ad928SBarry Smith static int PCSetUp_Jacobi_Symmetric(PC pc) 181*4b9ad928SBarry Smith { 182*4b9ad928SBarry Smith int ierr; 183*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 184*4b9ad928SBarry Smith 185*4b9ad928SBarry Smith PetscFunctionBegin; 186*4b9ad928SBarry Smith 187*4b9ad928SBarry Smith ierr = VecDuplicate(pc->vec,&jac->diagsqrt);CHKERRQ(ierr); 188*4b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diagsqrt); 189*4b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 190*4b9ad928SBarry Smith PetscFunctionReturn(0); 191*4b9ad928SBarry Smith } 192*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 193*4b9ad928SBarry Smith /* 194*4b9ad928SBarry Smith PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 195*4b9ad928SBarry Smith inverse of the diagonal entries of the matrix. This is used for left of 196*4b9ad928SBarry Smith right application of the Jacobi preconditioner. 197*4b9ad928SBarry Smith 198*4b9ad928SBarry Smith Input Parameter: 199*4b9ad928SBarry Smith . pc - the preconditioner context 200*4b9ad928SBarry Smith */ 201*4b9ad928SBarry Smith #undef __FUNCT__ 202*4b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 203*4b9ad928SBarry Smith static int PCSetUp_Jacobi_NonSymmetric(PC pc) 204*4b9ad928SBarry Smith { 205*4b9ad928SBarry Smith int ierr; 206*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 207*4b9ad928SBarry Smith 208*4b9ad928SBarry Smith PetscFunctionBegin; 209*4b9ad928SBarry Smith 210*4b9ad928SBarry Smith ierr = VecDuplicate(pc->vec,&jac->diag);CHKERRQ(ierr); 211*4b9ad928SBarry Smith PetscLogObjectParent(pc,jac->diag); 212*4b9ad928SBarry Smith ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 213*4b9ad928SBarry Smith PetscFunctionReturn(0); 214*4b9ad928SBarry Smith } 215*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 216*4b9ad928SBarry Smith /* 217*4b9ad928SBarry Smith PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 218*4b9ad928SBarry Smith 219*4b9ad928SBarry Smith Input Parameters: 220*4b9ad928SBarry Smith . pc - the preconditioner context 221*4b9ad928SBarry Smith . x - input vector 222*4b9ad928SBarry Smith 223*4b9ad928SBarry Smith Output Parameter: 224*4b9ad928SBarry Smith . y - output vector 225*4b9ad928SBarry Smith 226*4b9ad928SBarry Smith Application Interface Routine: PCApply() 227*4b9ad928SBarry Smith */ 228*4b9ad928SBarry Smith #undef __FUNCT__ 229*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_Jacobi" 230*4b9ad928SBarry Smith static int PCApply_Jacobi(PC pc,Vec x,Vec y) 231*4b9ad928SBarry Smith { 232*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 233*4b9ad928SBarry Smith int ierr; 234*4b9ad928SBarry Smith 235*4b9ad928SBarry Smith PetscFunctionBegin; 236*4b9ad928SBarry Smith if (!jac->diag) { 237*4b9ad928SBarry Smith ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 238*4b9ad928SBarry Smith } 239*4b9ad928SBarry Smith ierr = VecPointwiseMult(x,jac->diag,y);CHKERRQ(ierr); 240*4b9ad928SBarry Smith PetscFunctionReturn(0); 241*4b9ad928SBarry Smith } 242*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 243*4b9ad928SBarry Smith /* 244*4b9ad928SBarry Smith PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 245*4b9ad928SBarry Smith symmetric preconditioner to a vector. 246*4b9ad928SBarry Smith 247*4b9ad928SBarry Smith Input Parameters: 248*4b9ad928SBarry Smith . pc - the preconditioner context 249*4b9ad928SBarry Smith . x - input vector 250*4b9ad928SBarry Smith 251*4b9ad928SBarry Smith Output Parameter: 252*4b9ad928SBarry Smith . y - output vector 253*4b9ad928SBarry Smith 254*4b9ad928SBarry Smith Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 255*4b9ad928SBarry Smith */ 256*4b9ad928SBarry Smith #undef __FUNCT__ 257*4b9ad928SBarry Smith #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 258*4b9ad928SBarry Smith static int PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 259*4b9ad928SBarry Smith { 260*4b9ad928SBarry Smith int ierr; 261*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 262*4b9ad928SBarry Smith 263*4b9ad928SBarry Smith PetscFunctionBegin; 264*4b9ad928SBarry Smith if (!jac->diagsqrt) { 265*4b9ad928SBarry Smith ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 266*4b9ad928SBarry Smith } 267*4b9ad928SBarry Smith VecPointwiseMult(x,jac->diagsqrt,y); 268*4b9ad928SBarry Smith PetscFunctionReturn(0); 269*4b9ad928SBarry Smith } 270*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 271*4b9ad928SBarry Smith /* 272*4b9ad928SBarry Smith PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 273*4b9ad928SBarry Smith that was created with PCCreate_Jacobi(). 274*4b9ad928SBarry Smith 275*4b9ad928SBarry Smith Input Parameter: 276*4b9ad928SBarry Smith . pc - the preconditioner context 277*4b9ad928SBarry Smith 278*4b9ad928SBarry Smith Application Interface Routine: PCDestroy() 279*4b9ad928SBarry Smith */ 280*4b9ad928SBarry Smith #undef __FUNCT__ 281*4b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Jacobi" 282*4b9ad928SBarry Smith static int PCDestroy_Jacobi(PC pc) 283*4b9ad928SBarry Smith { 284*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 285*4b9ad928SBarry Smith int ierr; 286*4b9ad928SBarry Smith 287*4b9ad928SBarry Smith PetscFunctionBegin; 288*4b9ad928SBarry Smith if (jac->diag) {ierr = VecDestroy(jac->diag);CHKERRQ(ierr);} 289*4b9ad928SBarry Smith if (jac->diagsqrt) {ierr = VecDestroy(jac->diagsqrt);CHKERRQ(ierr);} 290*4b9ad928SBarry Smith 291*4b9ad928SBarry Smith /* 292*4b9ad928SBarry Smith Free the private data structure that was hanging off the PC 293*4b9ad928SBarry Smith */ 294*4b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 295*4b9ad928SBarry Smith PetscFunctionReturn(0); 296*4b9ad928SBarry Smith } 297*4b9ad928SBarry Smith 298*4b9ad928SBarry Smith #undef __FUNCT__ 299*4b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Jacobi" 300*4b9ad928SBarry Smith static int PCSetFromOptions_Jacobi(PC pc) 301*4b9ad928SBarry Smith { 302*4b9ad928SBarry Smith PC_Jacobi *jac = (PC_Jacobi*)pc->data; 303*4b9ad928SBarry Smith int ierr; 304*4b9ad928SBarry Smith 305*4b9ad928SBarry Smith PetscFunctionBegin; 306*4b9ad928SBarry Smith ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 307*4b9ad928SBarry Smith ierr = PetscOptionsLogical("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 308*4b9ad928SBarry Smith &jac->userowmax,PETSC_NULL);CHKERRQ(ierr); 309*4b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 310*4b9ad928SBarry Smith PetscFunctionReturn(0); 311*4b9ad928SBarry Smith } 312*4b9ad928SBarry Smith 313*4b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 314*4b9ad928SBarry Smith /* 315*4b9ad928SBarry Smith PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 316*4b9ad928SBarry Smith and sets this as the private data within the generic preconditioning 317*4b9ad928SBarry Smith context, PC, that was created within PCCreate(). 318*4b9ad928SBarry Smith 319*4b9ad928SBarry Smith Input Parameter: 320*4b9ad928SBarry Smith . pc - the preconditioner context 321*4b9ad928SBarry Smith 322*4b9ad928SBarry Smith Application Interface Routine: PCCreate() 323*4b9ad928SBarry Smith */ 324*4b9ad928SBarry Smith 325*4b9ad928SBarry Smith /*MC 326*4b9ad928SBarry Smith PCJacobi - Jacobi (i.e. diagonal scaling preconditioning) 327*4b9ad928SBarry Smith 328*4b9ad928SBarry Smith Options Database Key: 329*4b9ad928SBarry Smith . -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 330*4b9ad928SBarry Smith rather than the diagonal 331*4b9ad928SBarry Smith 332*4b9ad928SBarry Smith Level: beginner 333*4b9ad928SBarry Smith 334*4b9ad928SBarry Smith Concepts: Jacobi, diagonal scaling, preconditioners 335*4b9ad928SBarry Smith 336*4b9ad928SBarry Smith Notes: By using KSPSetPreconditionerSide(ksp,PC_SYMMETRIC) or -ksp_symmetric_pc you 337*4b9ad928SBarry Smith can scale each side of the matrix by the squareroot of the diagonal entries. 338*4b9ad928SBarry Smith 339*4b9ad928SBarry Smith Zero entries along the diagonal are replaced with the value 1.0 340*4b9ad928SBarry Smith 341*4b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 342*4b9ad928SBarry Smith PCJacobiSetUseRowMax(), 343*4b9ad928SBarry Smith M*/ 344*4b9ad928SBarry Smith 345*4b9ad928SBarry Smith EXTERN_C_BEGIN 346*4b9ad928SBarry Smith #undef __FUNCT__ 347*4b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Jacobi" 348*4b9ad928SBarry Smith int PCCreate_Jacobi(PC pc) 349*4b9ad928SBarry Smith { 350*4b9ad928SBarry Smith PC_Jacobi *jac; 351*4b9ad928SBarry Smith int ierr; 352*4b9ad928SBarry Smith 353*4b9ad928SBarry Smith PetscFunctionBegin; 354*4b9ad928SBarry Smith 355*4b9ad928SBarry Smith /* 356*4b9ad928SBarry Smith Creates the private data structure for this preconditioner and 357*4b9ad928SBarry Smith attach it to the PC object. 358*4b9ad928SBarry Smith */ 359*4b9ad928SBarry Smith ierr = PetscNew(PC_Jacobi,&jac);CHKERRQ(ierr); 360*4b9ad928SBarry Smith pc->data = (void*)jac; 361*4b9ad928SBarry Smith 362*4b9ad928SBarry Smith /* 363*4b9ad928SBarry Smith Logs the memory usage; this is not needed but allows PETSc to 364*4b9ad928SBarry Smith monitor how much memory is being used for various purposes. 365*4b9ad928SBarry Smith */ 366*4b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Jacobi)); 367*4b9ad928SBarry Smith 368*4b9ad928SBarry Smith /* 369*4b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 370*4b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 371*4b9ad928SBarry Smith */ 372*4b9ad928SBarry Smith jac->diag = 0; 373*4b9ad928SBarry Smith jac->diagsqrt = 0; 374*4b9ad928SBarry Smith jac->userowmax = PETSC_FALSE; 375*4b9ad928SBarry Smith 376*4b9ad928SBarry Smith /* 377*4b9ad928SBarry Smith Set the pointers for the functions that are provided above. 378*4b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 379*4b9ad928SBarry Smith are called, they will automatically call these functions. Note we 380*4b9ad928SBarry Smith choose not to provide a couple of these functions since they are 381*4b9ad928SBarry Smith not needed. 382*4b9ad928SBarry Smith */ 383*4b9ad928SBarry Smith pc->ops->apply = PCApply_Jacobi; 384*4b9ad928SBarry Smith pc->ops->applytranspose = PCApply_Jacobi; 385*4b9ad928SBarry Smith pc->ops->setup = PCSetUp_Jacobi; 386*4b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Jacobi; 387*4b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 388*4b9ad928SBarry Smith pc->ops->view = 0; 389*4b9ad928SBarry Smith pc->ops->applyrichardson = 0; 390*4b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 391*4b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 392*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCJacobiSetUseRowMax_C","PCJacobiSetUseRowMax_Jacobi", 393*4b9ad928SBarry Smith PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 394*4b9ad928SBarry Smith PetscFunctionReturn(0); 395*4b9ad928SBarry Smith } 396*4b9ad928SBarry Smith EXTERN_C_END 397*4b9ad928SBarry Smith 398*4b9ad928SBarry Smith #undef __FUNCT__ 399*4b9ad928SBarry Smith #define __FUNCT__ "PCJacobiSetUseRowMax" 400*4b9ad928SBarry Smith /*@ 401*4b9ad928SBarry Smith PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 402*4b9ad928SBarry Smith maximum entry in each row as the diagonal preconditioner, instead of 403*4b9ad928SBarry Smith the diagonal entry 404*4b9ad928SBarry Smith 405*4b9ad928SBarry Smith Collective on PC 406*4b9ad928SBarry Smith 407*4b9ad928SBarry Smith Input Parameters: 408*4b9ad928SBarry Smith . pc - the preconditioner context 409*4b9ad928SBarry Smith 410*4b9ad928SBarry Smith 411*4b9ad928SBarry Smith Options Database Key: 412*4b9ad928SBarry Smith . -pc_jacobi_rowmax 413*4b9ad928SBarry Smith 414*4b9ad928SBarry Smith Level: intermediate 415*4b9ad928SBarry Smith 416*4b9ad928SBarry Smith Concepts: Jacobi preconditioner 417*4b9ad928SBarry Smith 418*4b9ad928SBarry Smith @*/ 419*4b9ad928SBarry Smith int PCJacobiSetUseRowMax(PC pc) 420*4b9ad928SBarry Smith { 421*4b9ad928SBarry Smith int ierr,(*f)(PC); 422*4b9ad928SBarry Smith 423*4b9ad928SBarry Smith PetscFunctionBegin; 424*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 425*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCJacobiSetRowMax_C",(void (**)(void))&f);CHKERRQ(ierr); 426*4b9ad928SBarry Smith if (f) { 427*4b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 428*4b9ad928SBarry Smith } 429*4b9ad928SBarry Smith PetscFunctionReturn(0); 430*4b9ad928SBarry Smith } 431*4b9ad928SBarry Smith 432