1 2 /* -------------------------------------------------------------------- 3 4 This file implements a Jacobi preconditioner in PETSc as part of PC. 5 You can use this as a starting point for implementing your own 6 preconditioner that is not provided with PETSc. (You might also consider 7 just using PCSHELL) 8 9 The following basic routines are required for each preconditioner. 10 PCCreate_XXX() - Creates a preconditioner context 11 PCSetFromOptions_XXX() - Sets runtime options 12 PCApply_XXX() - Applies the preconditioner 13 PCDestroy_XXX() - Destroys the preconditioner context 14 where the suffix "_XXX" denotes a particular implementation, in 15 this case we use _Jacobi (e.g., PCCreate_Jacobi, PCApply_Jacobi). 16 These routines are actually called via the common user interface 17 routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(), 18 so the application code interface remains identical for all 19 preconditioners. 20 21 Another key routine is: 22 PCSetUp_XXX() - Prepares for the use of a preconditioner 23 by setting data structures and options. The interface routine PCSetUp() 24 is not usually called directly by the user, but instead is called by 25 PCApply() if necessary. 26 27 Additional basic routines are: 28 PCView_XXX() - Prints details of runtime options that 29 have actually been used. 30 These are called by application codes via the interface routines 31 PCView(). 32 33 The various types of solvers (preconditioners, Krylov subspace methods, 34 nonlinear solvers, timesteppers) are all organized similarly, so the 35 above description applies to these categories also. One exception is 36 that the analogues of PCApply() for these components are KSPSolve(), 37 SNESSolve(), and TSSolve(). 38 39 Additional optional functionality unique to preconditioners is left and 40 right symmetric preconditioner application via PCApplySymmetricLeft() 41 and PCApplySymmetricRight(). The Jacobi implementation is 42 PCApplySymmetricLeftOrRight_Jacobi(). 43 44 -------------------------------------------------------------------- */ 45 46 /* 47 Include files needed for the Jacobi preconditioner: 48 pcimpl.h - private include file intended for use by all preconditioners 49 */ 50 51 #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 52 53 /* 54 Private context (data structure) for the Jacobi preconditioner. 55 */ 56 typedef struct { 57 Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 58 Vec diagsqrt; /* vector containing the reciprocals of the square roots of 59 the diagonal elements of the preconditioner matrix (used 60 only for symmetric preconditioner application) */ 61 PetscBool userowmax; 62 PetscBool userowsum; 63 PetscBool useabs; /* use the absolute values of the diagonal entries */ 64 } PC_Jacobi; 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCJacobiSetUseRowMax_Jacobi" 68 static PetscErrorCode PCJacobiSetUseRowMax_Jacobi(PC pc) 69 { 70 PC_Jacobi *j; 71 72 PetscFunctionBegin; 73 j = (PC_Jacobi*)pc->data; 74 j->userowmax = PETSC_TRUE; 75 PetscFunctionReturn(0); 76 } 77 78 #undef __FUNCT__ 79 #define __FUNCT__ "PCJacobiSetUseRowSum_Jacobi" 80 static PetscErrorCode PCJacobiSetUseRowSum_Jacobi(PC pc) 81 { 82 PC_Jacobi *j; 83 84 PetscFunctionBegin; 85 j = (PC_Jacobi*)pc->data; 86 j->userowsum = PETSC_TRUE; 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 92 static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc) 93 { 94 PC_Jacobi *j; 95 96 PetscFunctionBegin; 97 j = (PC_Jacobi*)pc->data; 98 j->useabs = PETSC_TRUE; 99 PetscFunctionReturn(0); 100 } 101 102 /* -------------------------------------------------------------------------- */ 103 /* 104 PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 105 by setting data structures and options. 106 107 Input Parameter: 108 . pc - the preconditioner context 109 110 Application Interface Routine: PCSetUp() 111 112 Notes: 113 The interface routine PCSetUp() is not usually called directly by 114 the user, but instead is called by PCApply() if necessary. 115 */ 116 #undef __FUNCT__ 117 #define __FUNCT__ "PCSetUp_Jacobi" 118 static PetscErrorCode PCSetUp_Jacobi(PC pc) 119 { 120 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 121 Vec diag,diagsqrt; 122 PetscErrorCode ierr; 123 PetscInt n,i; 124 PetscScalar *x; 125 PetscBool zeroflag = PETSC_FALSE; 126 127 PetscFunctionBegin; 128 /* 129 For most preconditioners the code would begin here something like 130 131 if (pc->setupcalled == 0) { allocate space the first time this is ever called 132 ierr = MatGetVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 133 PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag); 134 } 135 136 But for this preconditioner we want to support use of both the matrix' diagonal 137 elements (for left or right preconditioning) and square root of diagonal elements 138 (for symmetric preconditioning). Hence we do not allocate space here, since we 139 don't know at this point which will be needed (diag and/or diagsqrt) until the user 140 applies the preconditioner, and we don't want to allocate BOTH unless we need 141 them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 142 and PCSetUp_Jacobi_Symmetric(), respectively. 143 */ 144 145 /* 146 Here we set up the preconditioner; that is, we copy the diagonal values from 147 the matrix and put them into a format to make them quick to apply as a preconditioner. 148 */ 149 diag = jac->diag; 150 diagsqrt = jac->diagsqrt; 151 152 if (diag) { 153 if (jac->userowmax) { 154 ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr); 155 } else if (jac->userowsum) { 156 ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 157 } else { 158 ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 159 } 160 ierr = VecReciprocal(diag);CHKERRQ(ierr); 161 ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 162 ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 163 if (jac->useabs) { 164 for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 165 } 166 for (i=0; i<n; i++) { 167 if (x[i] == 0.0) { 168 x[i] = 1.0; 169 zeroflag = PETSC_TRUE; 170 } 171 } 172 ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 173 } 174 if (diagsqrt) { 175 if (jac->userowmax) { 176 ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr); 177 } else if (jac->userowsum) { 178 ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 179 } else { 180 ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 181 } 182 ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 183 ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 184 for (i=0; i<n; i++) { 185 if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 186 else { 187 x[i] = 1.0; 188 zeroflag = PETSC_TRUE; 189 } 190 } 191 ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 192 } 193 if (zeroflag) { 194 ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 195 } 196 PetscFunctionReturn(0); 197 } 198 /* -------------------------------------------------------------------------- */ 199 /* 200 PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 201 inverse of the square root of the diagonal entries of the matrix. This 202 is used for symmetric application of the Jacobi preconditioner. 203 204 Input Parameter: 205 . pc - the preconditioner context 206 */ 207 #undef __FUNCT__ 208 #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 209 static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 210 { 211 PetscErrorCode ierr; 212 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 213 214 PetscFunctionBegin; 215 ierr = MatGetVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 216 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt);CHKERRQ(ierr); 217 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 218 PetscFunctionReturn(0); 219 } 220 /* -------------------------------------------------------------------------- */ 221 /* 222 PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 223 inverse of the diagonal entries of the matrix. This is used for left of 224 right application of the Jacobi preconditioner. 225 226 Input Parameter: 227 . pc - the preconditioner context 228 */ 229 #undef __FUNCT__ 230 #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 231 static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 232 { 233 PetscErrorCode ierr; 234 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 235 236 PetscFunctionBegin; 237 ierr = MatGetVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 238 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);CHKERRQ(ierr); 239 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 /* -------------------------------------------------------------------------- */ 243 /* 244 PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 245 246 Input Parameters: 247 . pc - the preconditioner context 248 . x - input vector 249 250 Output Parameter: 251 . y - output vector 252 253 Application Interface Routine: PCApply() 254 */ 255 #undef __FUNCT__ 256 #define __FUNCT__ "PCApply_Jacobi" 257 static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 258 { 259 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 260 PetscErrorCode ierr; 261 262 PetscFunctionBegin; 263 if (!jac->diag) { 264 ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 265 } 266 ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 267 PetscFunctionReturn(0); 268 } 269 /* -------------------------------------------------------------------------- */ 270 /* 271 PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 272 symmetric preconditioner to a vector. 273 274 Input Parameters: 275 . pc - the preconditioner context 276 . x - input vector 277 278 Output Parameter: 279 . y - output vector 280 281 Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 282 */ 283 #undef __FUNCT__ 284 #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 285 static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 286 { 287 PetscErrorCode ierr; 288 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 289 290 PetscFunctionBegin; 291 if (!jac->diagsqrt) { 292 ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 293 } 294 VecPointwiseMult(y,x,jac->diagsqrt); 295 PetscFunctionReturn(0); 296 } 297 /* -------------------------------------------------------------------------- */ 298 #undef __FUNCT__ 299 #define __FUNCT__ "PCReset_Jacobi" 300 static PetscErrorCode PCReset_Jacobi(PC pc) 301 { 302 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 303 PetscErrorCode ierr; 304 305 PetscFunctionBegin; 306 ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 307 ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 308 PetscFunctionReturn(0); 309 } 310 311 /* 312 PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 313 that was created with PCCreate_Jacobi(). 314 315 Input Parameter: 316 . pc - the preconditioner context 317 318 Application Interface Routine: PCDestroy() 319 */ 320 #undef __FUNCT__ 321 #define __FUNCT__ "PCDestroy_Jacobi" 322 static PetscErrorCode PCDestroy_Jacobi(PC pc) 323 { 324 PetscErrorCode ierr; 325 326 PetscFunctionBegin; 327 ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 328 329 /* 330 Free the private data structure that was hanging off the PC 331 */ 332 ierr = PetscFree(pc->data);CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "PCSetFromOptions_Jacobi" 338 static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 339 { 340 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 341 PetscErrorCode ierr; 342 343 PetscFunctionBegin; 344 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 345 ierr = PetscOptionsBool("-pc_jacobi_rowmax","Use row maximums for diagonal","PCJacobiSetUseRowMax",jac->userowmax, 346 &jac->userowmax,NULL);CHKERRQ(ierr); 347 ierr = PetscOptionsBool("-pc_jacobi_rowsum","Use row sums for diagonal","PCJacobiSetUseRowSum",jac->userowsum, 348 &jac->userowsum,NULL);CHKERRQ(ierr); 349 ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs, 350 &jac->useabs,NULL);CHKERRQ(ierr); 351 ierr = PetscOptionsTail();CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 /* -------------------------------------------------------------------------- */ 356 /* 357 PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 358 and sets this as the private data within the generic preconditioning 359 context, PC, that was created within PCCreate(). 360 361 Input Parameter: 362 . pc - the preconditioner context 363 364 Application Interface Routine: PCCreate() 365 */ 366 367 /*MC 368 PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 369 370 Options Database Key: 371 + -pc_jacobi_rowmax - use the maximum absolute value in each row as the scaling factor, 372 rather than the diagonal 373 . -pc_jacobi_rowsum - use the row sums as the scaling factor, 374 rather than the diagonal 375 - -pc_jacobi_abs - use the absolute value of the diagaonl entry 376 377 Level: beginner 378 379 Concepts: Jacobi, diagonal scaling, preconditioners 380 381 Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 382 can scale each side of the matrix by the squareroot of the diagonal entries. 383 384 Zero entries along the diagonal are replaced with the value 1.0 385 386 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 387 PCJacobiSetUseRowMax(), PCJacobiSetUseRowSum(), PCJacobiSetUseAbs() 388 M*/ 389 390 #undef __FUNCT__ 391 #define __FUNCT__ "PCCreate_Jacobi" 392 PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) 393 { 394 PC_Jacobi *jac; 395 PetscErrorCode ierr; 396 397 PetscFunctionBegin; 398 /* 399 Creates the private data structure for this preconditioner and 400 attach it to the PC object. 401 */ 402 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 403 pc->data = (void*)jac; 404 405 /* 406 Initialize the pointers to vectors to ZERO; these will be used to store 407 diagonal entries of the matrix for fast preconditioner application. 408 */ 409 jac->diag = 0; 410 jac->diagsqrt = 0; 411 jac->userowmax = PETSC_FALSE; 412 jac->userowsum = PETSC_FALSE; 413 jac->useabs = PETSC_FALSE; 414 415 /* 416 Set the pointers for the functions that are provided above. 417 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 418 are called, they will automatically call these functions. Note we 419 choose not to provide a couple of these functions since they are 420 not needed. 421 */ 422 pc->ops->apply = PCApply_Jacobi; 423 pc->ops->applytranspose = PCApply_Jacobi; 424 pc->ops->setup = PCSetUp_Jacobi; 425 pc->ops->reset = PCReset_Jacobi; 426 pc->ops->destroy = PCDestroy_Jacobi; 427 pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 428 pc->ops->view = 0; 429 pc->ops->applyrichardson = 0; 430 pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 431 pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 432 433 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseRowMax_C",PCJacobiSetUseRowMax_Jacobi);CHKERRQ(ierr); 434 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseRowSum_C",PCJacobiSetUseRowSum_Jacobi);CHKERRQ(ierr); 435 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "PCJacobiSetUseAbs" 441 /*@ 442 PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 443 absolute value of the diagonal to for the preconditioner 444 445 Logically Collective on PC 446 447 Input Parameters: 448 . pc - the preconditioner context 449 450 451 Options Database Key: 452 . -pc_jacobi_abs 453 454 Level: intermediate 455 456 Concepts: Jacobi preconditioner 457 458 .seealso: PCJacobiaUseRowMax(), PCJacobiaUseRowSum() 459 460 @*/ 461 PetscErrorCode PCJacobiSetUseAbs(PC pc) 462 { 463 PetscErrorCode ierr; 464 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 467 ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC),(pc));CHKERRQ(ierr); 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "PCJacobiSetUseRowMax" 473 /*@ 474 PCJacobiSetUseRowMax - Causes the Jacobi preconditioner to use the 475 maximum entry in each row as the diagonal preconditioner, instead of 476 the diagonal entry 477 478 Logically Collective on PC 479 480 Input Parameters: 481 . pc - the preconditioner context 482 483 484 Options Database Key: 485 . -pc_jacobi_rowmax 486 487 Level: intermediate 488 489 Concepts: Jacobi preconditioner 490 491 .seealso: PCJacobiaUseAbs() 492 @*/ 493 PetscErrorCode PCJacobiSetUseRowMax(PC pc) 494 { 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 499 ierr = PetscTryMethod(pc,"PCJacobiSetUseRowMax_C",(PC),(pc));CHKERRQ(ierr); 500 PetscFunctionReturn(0); 501 } 502 503 #undef __FUNCT__ 504 #define __FUNCT__ "PCJacobiSetUseRowSum" 505 /*@ 506 PCJacobiSetUseRowSum - Causes the Jacobi preconditioner to use the 507 sum of each row as the diagonal preconditioner, instead of 508 the diagonal entry 509 510 Logical Collective on PC 511 512 Input Parameters: 513 . pc - the preconditioner context 514 515 516 Options Database Key: 517 . -pc_jacobi_rowsum 518 519 Level: intermediate 520 521 Concepts: Jacobi preconditioner 522 523 .seealso: PCJacobiaUseAbs(), PCJacobiaUseRowSum() 524 @*/ 525 PetscErrorCode PCJacobiSetUseRowSum(PC pc) 526 { 527 PetscErrorCode ierr; 528 529 PetscFunctionBegin; 530 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 531 ierr = PetscTryMethod(pc,"PCJacobiSetUseRowSum_C",(PC),(pc));CHKERRQ(ierr); 532 PetscFunctionReturn(0); 533 } 534 535