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 const char *const PCJacobiTypes[] = {"DIAGONAL","ROWMAX","ROWSUM","PCJacobiType","PC_JACOBI_",0}; 54 55 /* 56 Private context (data structure) for the Jacobi preconditioner. 57 */ 58 typedef struct { 59 Vec diag; /* vector containing the reciprocals of the diagonal elements of the preconditioner matrix */ 60 Vec diagsqrt; /* vector containing the reciprocals of the square roots of 61 the diagonal elements of the preconditioner matrix (used 62 only for symmetric preconditioner application) */ 63 PetscBool userowmax; /* set with PCJacobiSetType() */ 64 PetscBool userowsum; 65 PetscBool useabs; /* use the absolute values of the diagonal entries */ 66 } PC_Jacobi; 67 68 #undef __FUNCT__ 69 #define __FUNCT__ "PCJacobiSetType_Jacobi" 70 static PetscErrorCode PCJacobiSetType_Jacobi(PC pc,PCJacobiType type) 71 { 72 PC_Jacobi *j = (PC_Jacobi*)pc->data; 73 74 PetscFunctionBegin; 75 j->userowmax = PETSC_FALSE; 76 j->userowsum = PETSC_FALSE; 77 if (type == PC_JACOBI_ROWMAX) { 78 j->userowmax = PETSC_TRUE; 79 } else if (type == PC_JACOBI_ROWSUM) { 80 j->userowsum = PETSC_TRUE; 81 } 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "PCJacobiGetType_Jacobi" 87 static PetscErrorCode PCJacobiGetType_Jacobi(PC pc,PCJacobiType *type) 88 { 89 PC_Jacobi *j = (PC_Jacobi*)pc->data; 90 91 PetscFunctionBegin; 92 if (j->userowmax) { 93 *type = PC_JACOBI_ROWMAX; 94 } else if (j->userowsum) { 95 *type = PC_JACOBI_ROWSUM; 96 } else { 97 *type = PC_JACOBI_DIAGONAL; 98 } 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "PCJacobiSetUseAbs_Jacobi" 104 static PetscErrorCode PCJacobiSetUseAbs_Jacobi(PC pc,PetscBool flg) 105 { 106 PC_Jacobi *j = (PC_Jacobi*)pc->data; 107 108 PetscFunctionBegin; 109 j->useabs = flg; 110 PetscFunctionReturn(0); 111 } 112 113 #undef __FUNCT__ 114 #define __FUNCT__ "PCJacobiGetUseAbs_Jacobi" 115 static PetscErrorCode PCJacobiGetUseAbs_Jacobi(PC pc,PetscBool *flg) 116 { 117 PC_Jacobi *j = (PC_Jacobi*)pc->data; 118 119 PetscFunctionBegin; 120 *flg = j->useabs; 121 PetscFunctionReturn(0); 122 } 123 124 /* -------------------------------------------------------------------------- */ 125 /* 126 PCSetUp_Jacobi - Prepares for the use of the Jacobi preconditioner 127 by setting data structures and options. 128 129 Input Parameter: 130 . pc - the preconditioner context 131 132 Application Interface Routine: PCSetUp() 133 134 Notes: 135 The interface routine PCSetUp() is not usually called directly by 136 the user, but instead is called by PCApply() if necessary. 137 */ 138 #undef __FUNCT__ 139 #define __FUNCT__ "PCSetUp_Jacobi" 140 static PetscErrorCode PCSetUp_Jacobi(PC pc) 141 { 142 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 143 Vec diag,diagsqrt; 144 PetscErrorCode ierr; 145 PetscInt n,i; 146 PetscScalar *x; 147 PetscBool zeroflag = PETSC_FALSE; 148 149 PetscFunctionBegin; 150 /* 151 For most preconditioners the code would begin here something like 152 153 if (pc->setupcalled == 0) { allocate space the first time this is ever called 154 ierr = MatCreateVecs(pc->mat,&jac->diag);CHKERRQ(ierr); 155 PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag); 156 } 157 158 But for this preconditioner we want to support use of both the matrix' diagonal 159 elements (for left or right preconditioning) and square root of diagonal elements 160 (for symmetric preconditioning). Hence we do not allocate space here, since we 161 don't know at this point which will be needed (diag and/or diagsqrt) until the user 162 applies the preconditioner, and we don't want to allocate BOTH unless we need 163 them both. Thus, the diag and diagsqrt are allocated in PCSetUp_Jacobi_NonSymmetric() 164 and PCSetUp_Jacobi_Symmetric(), respectively. 165 */ 166 167 /* 168 Here we set up the preconditioner; that is, we copy the diagonal values from 169 the matrix and put them into a format to make them quick to apply as a preconditioner. 170 */ 171 diag = jac->diag; 172 diagsqrt = jac->diagsqrt; 173 174 if (diag) { 175 if (jac->userowmax) { 176 ierr = MatGetRowMaxAbs(pc->pmat,diag,NULL);CHKERRQ(ierr); 177 } else if (jac->userowsum) { 178 ierr = MatGetRowSum(pc->pmat,diag);CHKERRQ(ierr); 179 } else { 180 ierr = MatGetDiagonal(pc->pmat,diag);CHKERRQ(ierr); 181 } 182 ierr = VecReciprocal(diag);CHKERRQ(ierr); 183 ierr = VecGetLocalSize(diag,&n);CHKERRQ(ierr); 184 ierr = VecGetArray(diag,&x);CHKERRQ(ierr); 185 if (jac->useabs) { 186 for (i=0; i<n; i++) x[i] = PetscAbsScalar(x[i]); 187 } 188 for (i=0; i<n; i++) { 189 if (x[i] == 0.0) { 190 x[i] = 1.0; 191 zeroflag = PETSC_TRUE; 192 } 193 } 194 ierr = VecRestoreArray(diag,&x);CHKERRQ(ierr); 195 } 196 if (diagsqrt) { 197 if (jac->userowmax) { 198 ierr = MatGetRowMaxAbs(pc->pmat,diagsqrt,NULL);CHKERRQ(ierr); 199 } else if (jac->userowsum) { 200 ierr = MatGetRowSum(pc->pmat,diagsqrt);CHKERRQ(ierr); 201 } else { 202 ierr = MatGetDiagonal(pc->pmat,diagsqrt);CHKERRQ(ierr); 203 } 204 ierr = VecGetLocalSize(diagsqrt,&n);CHKERRQ(ierr); 205 ierr = VecGetArray(diagsqrt,&x);CHKERRQ(ierr); 206 for (i=0; i<n; i++) { 207 if (x[i] != 0.0) x[i] = 1.0/PetscSqrtReal(PetscAbsScalar(x[i])); 208 else { 209 x[i] = 1.0; 210 zeroflag = PETSC_TRUE; 211 } 212 } 213 ierr = VecRestoreArray(diagsqrt,&x);CHKERRQ(ierr); 214 } 215 if (zeroflag) { 216 ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr); 217 } 218 PetscFunctionReturn(0); 219 } 220 /* -------------------------------------------------------------------------- */ 221 /* 222 PCSetUp_Jacobi_Symmetric - Allocates the vector needed to store the 223 inverse of the square root of the diagonal entries of the matrix. This 224 is used for symmetric application of the Jacobi preconditioner. 225 226 Input Parameter: 227 . pc - the preconditioner context 228 */ 229 #undef __FUNCT__ 230 #define __FUNCT__ "PCSetUp_Jacobi_Symmetric" 231 static PetscErrorCode PCSetUp_Jacobi_Symmetric(PC pc) 232 { 233 PetscErrorCode ierr; 234 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 235 236 PetscFunctionBegin; 237 ierr = MatCreateVecs(pc->pmat,&jac->diagsqrt,0);CHKERRQ(ierr); 238 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diagsqrt);CHKERRQ(ierr); 239 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 /* -------------------------------------------------------------------------- */ 243 /* 244 PCSetUp_Jacobi_NonSymmetric - Allocates the vector needed to store the 245 inverse of the diagonal entries of the matrix. This is used for left of 246 right application of the Jacobi preconditioner. 247 248 Input Parameter: 249 . pc - the preconditioner context 250 */ 251 #undef __FUNCT__ 252 #define __FUNCT__ "PCSetUp_Jacobi_NonSymmetric" 253 static PetscErrorCode PCSetUp_Jacobi_NonSymmetric(PC pc) 254 { 255 PetscErrorCode ierr; 256 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 257 258 PetscFunctionBegin; 259 ierr = MatCreateVecs(pc->pmat,&jac->diag,0);CHKERRQ(ierr); 260 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);CHKERRQ(ierr); 261 ierr = PCSetUp_Jacobi(pc);CHKERRQ(ierr); 262 PetscFunctionReturn(0); 263 } 264 /* -------------------------------------------------------------------------- */ 265 /* 266 PCApply_Jacobi - Applies the Jacobi preconditioner to a vector. 267 268 Input Parameters: 269 . pc - the preconditioner context 270 . x - input vector 271 272 Output Parameter: 273 . y - output vector 274 275 Application Interface Routine: PCApply() 276 */ 277 #undef __FUNCT__ 278 #define __FUNCT__ "PCApply_Jacobi" 279 static PetscErrorCode PCApply_Jacobi(PC pc,Vec x,Vec y) 280 { 281 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 282 PetscErrorCode ierr; 283 284 PetscFunctionBegin; 285 if (!jac->diag) { 286 ierr = PCSetUp_Jacobi_NonSymmetric(pc);CHKERRQ(ierr); 287 } 288 ierr = VecPointwiseMult(y,x,jac->diag);CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291 /* -------------------------------------------------------------------------- */ 292 /* 293 PCApplySymmetricLeftOrRight_Jacobi - Applies the left or right part of a 294 symmetric preconditioner to a vector. 295 296 Input Parameters: 297 . pc - the preconditioner context 298 . x - input vector 299 300 Output Parameter: 301 . y - output vector 302 303 Application Interface Routines: PCApplySymmetricLeft(), PCApplySymmetricRight() 304 */ 305 #undef __FUNCT__ 306 #define __FUNCT__ "PCApplySymmetricLeftOrRight_Jacobi" 307 static PetscErrorCode PCApplySymmetricLeftOrRight_Jacobi(PC pc,Vec x,Vec y) 308 { 309 PetscErrorCode ierr; 310 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 311 312 PetscFunctionBegin; 313 if (!jac->diagsqrt) { 314 ierr = PCSetUp_Jacobi_Symmetric(pc);CHKERRQ(ierr); 315 } 316 VecPointwiseMult(y,x,jac->diagsqrt); 317 PetscFunctionReturn(0); 318 } 319 /* -------------------------------------------------------------------------- */ 320 #undef __FUNCT__ 321 #define __FUNCT__ "PCReset_Jacobi" 322 static PetscErrorCode PCReset_Jacobi(PC pc) 323 { 324 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = VecDestroy(&jac->diag);CHKERRQ(ierr); 329 ierr = VecDestroy(&jac->diagsqrt);CHKERRQ(ierr); 330 PetscFunctionReturn(0); 331 } 332 333 /* 334 PCDestroy_Jacobi - Destroys the private context for the Jacobi preconditioner 335 that was created with PCCreate_Jacobi(). 336 337 Input Parameter: 338 . pc - the preconditioner context 339 340 Application Interface Routine: PCDestroy() 341 */ 342 #undef __FUNCT__ 343 #define __FUNCT__ "PCDestroy_Jacobi" 344 static PetscErrorCode PCDestroy_Jacobi(PC pc) 345 { 346 PetscErrorCode ierr; 347 348 PetscFunctionBegin; 349 ierr = PCReset_Jacobi(pc);CHKERRQ(ierr); 350 351 /* 352 Free the private data structure that was hanging off the PC 353 */ 354 ierr = PetscFree(pc->data);CHKERRQ(ierr); 355 PetscFunctionReturn(0); 356 } 357 358 #undef __FUNCT__ 359 #define __FUNCT__ "PCSetFromOptions_Jacobi" 360 static PetscErrorCode PCSetFromOptions_Jacobi(PC pc) 361 { 362 PC_Jacobi *jac = (PC_Jacobi*)pc->data; 363 PetscErrorCode ierr; 364 PetscBool flg; 365 PCJacobiType deflt,type; 366 367 PetscFunctionBegin; 368 ierr = PCJacobiGetType(pc,&deflt);CHKERRQ(ierr); 369 ierr = PetscOptionsHead("Jacobi options");CHKERRQ(ierr); 370 ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCJacobiSetType",PCJacobiTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 371 if (flg) { 372 ierr = PCJacobiSetType(pc,type);CHKERRQ(ierr); 373 } 374 ierr = PetscOptionsBool("-pc_jacobi_abs","Use absolute values of diagaonal entries","PCJacobiSetUseAbs",jac->useabs,&jac->useabs,NULL);CHKERRQ(ierr); 375 ierr = PetscOptionsTail();CHKERRQ(ierr); 376 PetscFunctionReturn(0); 377 } 378 379 /* -------------------------------------------------------------------------- */ 380 /* 381 PCCreate_Jacobi - Creates a Jacobi preconditioner context, PC_Jacobi, 382 and sets this as the private data within the generic preconditioning 383 context, PC, that was created within PCCreate(). 384 385 Input Parameter: 386 . pc - the preconditioner context 387 388 Application Interface Routine: PCCreate() 389 */ 390 391 /*MC 392 PCJACOBI - Jacobi (i.e. diagonal scaling preconditioning) 393 394 Options Database Key: 395 + -pc_jacobi_type <diagonal,rowmax,rowsum> 396 - -pc_jacobi_abs - use the absolute value of the diagaonl entry 397 398 Level: beginner 399 400 Concepts: Jacobi, diagonal scaling, preconditioners 401 402 Notes: By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric 403 can scale each side of the matrix by the squareroot of the diagonal entries. 404 405 Zero entries along the diagonal are replaced with the value 1.0 406 407 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 408 PCJacobiSetType(), PCJacobiSetUseAbs(), PCJacobiGetUseAbs() 409 M*/ 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "PCCreate_Jacobi" 413 PETSC_EXTERN PetscErrorCode PCCreate_Jacobi(PC pc) 414 { 415 PC_Jacobi *jac; 416 PetscErrorCode ierr; 417 418 PetscFunctionBegin; 419 /* 420 Creates the private data structure for this preconditioner and 421 attach it to the PC object. 422 */ 423 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 424 pc->data = (void*)jac; 425 426 /* 427 Initialize the pointers to vectors to ZERO; these will be used to store 428 diagonal entries of the matrix for fast preconditioner application. 429 */ 430 jac->diag = 0; 431 jac->diagsqrt = 0; 432 jac->userowmax = PETSC_FALSE; 433 jac->userowsum = PETSC_FALSE; 434 jac->useabs = PETSC_FALSE; 435 436 /* 437 Set the pointers for the functions that are provided above. 438 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 439 are called, they will automatically call these functions. Note we 440 choose not to provide a couple of these functions since they are 441 not needed. 442 */ 443 pc->ops->apply = PCApply_Jacobi; 444 pc->ops->applytranspose = PCApply_Jacobi; 445 pc->ops->setup = PCSetUp_Jacobi; 446 pc->ops->reset = PCReset_Jacobi; 447 pc->ops->destroy = PCDestroy_Jacobi; 448 pc->ops->setfromoptions = PCSetFromOptions_Jacobi; 449 pc->ops->view = 0; 450 pc->ops->applyrichardson = 0; 451 pc->ops->applysymmetricleft = PCApplySymmetricLeftOrRight_Jacobi; 452 pc->ops->applysymmetricright = PCApplySymmetricLeftOrRight_Jacobi; 453 454 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetType_C",PCJacobiSetType_Jacobi);CHKERRQ(ierr); 455 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetType_C",PCJacobiGetType_Jacobi);CHKERRQ(ierr); 456 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiSetUseAbs_C",PCJacobiSetUseAbs_Jacobi);CHKERRQ(ierr); 457 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCJacobiGetUseAbs_C",PCJacobiGetUseAbs_Jacobi);CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNCT__ 462 #define __FUNCT__ "PCJacobiSetUseAbs" 463 /*@ 464 PCJacobiSetUseAbs - Causes the Jacobi preconditioner to use the 465 absolute values of the digonal divisors in the preconditioner 466 467 Logically Collective on PC 468 469 Input Parameters: 470 + pc - the preconditioner context 471 - flg - whether to use absolute values or not 472 473 Options Database Key: 474 . -pc_jacobi_abs 475 476 Notes: This takes affect at the next construction of the preconditioner 477 478 Level: intermediate 479 480 Concepts: Jacobi preconditioner 481 482 .seealso: PCJacobiaSetType(), PCJacobiGetUseAbs() 483 484 @*/ 485 PetscErrorCode PCJacobiSetUseAbs(PC pc,PetscBool flg) 486 { 487 PetscErrorCode ierr; 488 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 491 ierr = PetscTryMethod(pc,"PCJacobiSetUseAbs_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 492 PetscFunctionReturn(0); 493 } 494 495 #undef __FUNCT__ 496 #define __FUNCT__ "PCJacobiGetUseAbs" 497 /*@ 498 PCJacobiGetUseAbs - Determines if the Jacobi preconditioner uses the 499 absolute values of the digonal divisors in the preconditioner 500 501 Logically Collective on PC 502 503 Input Parameter: 504 . pc - the preconditioner context 505 506 Output Parameter: 507 . flg - whether to use absolute values or not 508 509 Options Database Key: 510 . -pc_jacobi_abs 511 512 Level: intermediate 513 514 Concepts: Jacobi preconditioner 515 516 .seealso: PCJacobiaSetType(), PCJacobiSetUseAbs(), PCJacobiGetType() 517 518 @*/ 519 PetscErrorCode PCJacobiGetUseAbs(PC pc,PetscBool *flg) 520 { 521 PetscErrorCode ierr; 522 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 525 ierr = PetscUseMethod(pc,"PCJacobiGetUseAbs_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr); 526 PetscFunctionReturn(0); 527 } 528 529 #undef __FUNCT__ 530 #define __FUNCT__ "PCJacobiSetType" 531 /*@ 532 PCJacobiSetType - Causes the Jacobi preconditioner to use either the diagonal, the maximum entry in each row, 533 of the sum of rows entries for the diagonal preconditioner 534 535 Logically Collective on PC 536 537 Input Parameters: 538 + pc - the preconditioner context 539 - type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 540 541 Options Database Key: 542 . -pc_jacobi_type <diagonal,rowmax,rowsum> 543 544 Level: intermediate 545 546 Concepts: Jacobi preconditioner 547 548 .seealso: PCJacobiaUseAbs(), PCJacobiGetType() 549 @*/ 550 PetscErrorCode PCJacobiSetType(PC pc,PCJacobiType type) 551 { 552 PetscErrorCode ierr; 553 554 PetscFunctionBegin; 555 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 556 ierr = PetscTryMethod(pc,"PCJacobiSetType_C",(PC,PCJacobiType),(pc,type));CHKERRQ(ierr); 557 PetscFunctionReturn(0); 558 } 559 560 #undef __FUNCT__ 561 #define __FUNCT__ "PCJacobiGetType" 562 /*@ 563 PCJacobiGetType - Gets how the diagonal matrix is produced for the preconditioner 564 565 Not Collective on PC 566 567 Input Parameter: 568 . pc - the preconditioner context 569 570 Output Parameter: 571 . type - PC_JACOBI_DIAGONAL, PC_JACOBI_ROWMAX, PC_JACOBI_ROWSUM 572 573 Level: intermediate 574 575 Concepts: Jacobi preconditioner 576 577 .seealso: PCJacobiaUseAbs(), PCJacobiSetType() 578 @*/ 579 PetscErrorCode PCJacobiGetType(PC pc,PCJacobiType *type) 580 { 581 PetscErrorCode ierr; 582 583 PetscFunctionBegin; 584 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 585 ierr = PetscUseMethod(pc,"PCJacobiGetType_C",(PC,PCJacobiType*),(pc,type));CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588