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