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