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