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