1 2 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 3 #include <petsc/private/matimpl.h> 4 5 /* 6 If an ordering is not yet set and the matrix is available determine a default ordering 7 */ 8 PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc) { 9 PetscBool foundmtype, flg; 10 const char *prefix; 11 12 PetscFunctionBegin; 13 if (pc->pmat) { 14 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 15 PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 16 PC_Factor *fact = (PC_Factor *)pc->data; 17 PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL)); 18 if (foundmtype) { 19 if (!fact->fact) { 20 PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); 21 } else if (!fact->fact->assembled) { 22 PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg)); 23 if (!flg) { 24 Mat B; 25 PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B)); 26 PetscCall(MatHeaderReplace(fact->fact, &B)); 27 } 28 } 29 if (!fact->ordering) { 30 PetscBool canuseordering; 31 MatOrderingType otype; 32 33 PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering)); 34 if (canuseordering) { 35 PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype)); 36 } else otype = MATORDERINGEXTERNAL; 37 PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering)); 38 } 39 } 40 } 41 PetscFunctionReturn(0); 42 } 43 44 static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag) { 45 PC_Factor *lu = (PC_Factor *)pc->data; 46 47 PetscFunctionBegin; 48 lu->reuseordering = flag; 49 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag) { 53 PC_Factor *lu = (PC_Factor *)pc->data; 54 55 PetscFunctionBegin; 56 lu->reusefill = flag; 57 PetscFunctionReturn(0); 58 } 59 60 static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg) { 61 PC_Factor *dir = (PC_Factor *)pc->data; 62 63 PetscFunctionBegin; 64 dir->inplace = flg; 65 PetscFunctionReturn(0); 66 } 67 68 static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg) { 69 PC_Factor *dir = (PC_Factor *)pc->data; 70 71 PetscFunctionBegin; 72 *flg = dir->inplace; 73 PetscFunctionReturn(0); 74 } 75 76 /*@ 77 PCFactorSetUpMatSolverType - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may 78 set the options for that particular factorization object. 79 80 Input Parameter: 81 . pc - the preconditioner context 82 83 Notes: 84 After you have called this function (which has to be after the KSPSetOperators() or PCSetOperators()) you can call PCFactorGetMatrix() and then set factor options on that matrix. 85 86 Level: intermediate 87 88 .seealso: `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 89 @*/ 90 PetscErrorCode PCFactorSetUpMatSolverType(PC pc) { 91 PetscFunctionBegin; 92 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 93 PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 94 PetscFunctionReturn(0); 95 } 96 97 /*@ 98 PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 99 100 Logically Collective on PC 101 102 Input Parameters: 103 + pc - the preconditioner context 104 - zero - all pivots smaller than this will be considered zero 105 106 Options Database Key: 107 . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 108 109 Level: intermediate 110 111 .seealso: `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 112 @*/ 113 PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) { 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 116 PetscValidLogicalCollectiveReal(pc, zero, 2); 117 PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 118 PetscFunctionReturn(0); 119 } 120 121 /*@ 122 PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 123 numerical factorization, thus the matrix has nonzero pivots 124 125 Logically Collective on PC 126 127 Input Parameters: 128 + pc - the preconditioner context 129 - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS 130 131 Options Database Key: 132 . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 133 134 Level: intermediate 135 136 .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 137 @*/ 138 PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) { 139 PetscFunctionBegin; 140 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 141 PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 142 PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 143 PetscFunctionReturn(0); 144 } 145 146 /*@ 147 PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 148 numerical factorization, thus the matrix has nonzero pivots 149 150 Logically Collective on PC 151 152 Input Parameters: 153 + pc - the preconditioner context 154 - shiftamount - amount of shift 155 156 Options Database Key: 157 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 158 159 Level: intermediate 160 161 .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 162 @*/ 163 PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) { 164 PetscFunctionBegin; 165 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 166 PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 167 PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 168 PetscFunctionReturn(0); 169 } 170 171 /*@ 172 PCFactorSetDropTolerance - The preconditioner will use an ILU 173 based on a drop tolerance. (Under development) 174 175 Logically Collective on PC 176 177 Input Parameters: 178 + pc - the preconditioner context 179 . dt - the drop tolerance, try from 1.e-10 to .1 180 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 181 - maxrowcount - the max number of nonzeros allowed in a row, best value 182 depends on the number of nonzeros in row of original matrix 183 184 Options Database Key: 185 . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 186 187 Level: intermediate 188 189 There are NO default values for the 3 parameters, you must set them with reasonable values for your 190 matrix. We don't know how to compute reasonable values. 191 192 @*/ 193 PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) { 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 196 PetscValidLogicalCollectiveReal(pc, dtcol, 3); 197 PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 198 PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 199 PetscFunctionReturn(0); 200 } 201 202 /*@ 203 PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 204 205 Not Collective 206 207 Input Parameters: 208 . pc - the preconditioner context 209 210 Output Parameter: 211 . pivot - the tolerance 212 213 Level: intermediate 214 215 .seealso: `PCFactorSetZeroPivot()` 216 @*/ 217 PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) { 218 PetscFunctionBegin; 219 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 220 PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 221 PetscFunctionReturn(0); 222 } 223 224 /*@ 225 PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 226 227 Not Collective 228 229 Input Parameters: 230 . pc - the preconditioner context 231 232 Output Parameter: 233 . shift - how much to shift the diagonal entry 234 235 Level: intermediate 236 237 .seealso: `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 238 @*/ 239 PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) { 240 PetscFunctionBegin; 241 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 242 PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 243 PetscFunctionReturn(0); 244 } 245 246 /*@ 247 PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 248 249 Not Collective 250 251 Input Parameters: 252 . pc - the preconditioner context 253 254 Output Parameter: 255 . type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS 256 257 Level: intermediate 258 259 .seealso: `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 260 @*/ 261 PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) { 262 PetscFunctionBegin; 263 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 264 PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 265 PetscFunctionReturn(0); 266 } 267 268 /*@ 269 PCFactorGetLevels - Gets the number of levels of fill to use. 270 271 Logically Collective on PC 272 273 Input Parameters: 274 . pc - the preconditioner context 275 276 Output Parameter: 277 . levels - number of levels of fill 278 279 Level: intermediate 280 281 @*/ 282 PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) { 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 285 PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 286 PetscFunctionReturn(0); 287 } 288 289 /*@ 290 PCFactorSetLevels - Sets the number of levels of fill to use. 291 292 Logically Collective on PC 293 294 Input Parameters: 295 + pc - the preconditioner context 296 - levels - number of levels of fill 297 298 Options Database Key: 299 . -pc_factor_levels <levels> - Sets fill level 300 301 Level: intermediate 302 303 @*/ 304 PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) { 305 PetscFunctionBegin; 306 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 307 PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 308 PetscValidLogicalCollectiveInt(pc, levels, 2); 309 PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 310 PetscFunctionReturn(0); 311 } 312 313 /*@ 314 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 315 treated as level 0 fill even if there is no non-zero location. 316 317 Logically Collective on PC 318 319 Input Parameters: 320 + pc - the preconditioner context 321 - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 322 323 Options Database Key: 324 . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 325 326 Notes: 327 Does not apply with 0 fill. 328 329 Level: intermediate 330 331 .seealso: `PCFactorGetAllowDiagonalFill()` 332 @*/ 333 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) { 334 PetscFunctionBegin; 335 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 336 PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 337 PetscFunctionReturn(0); 338 } 339 340 /*@ 341 PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 342 treated as level 0 fill even if there is no non-zero location. 343 344 Logically Collective on PC 345 346 Input Parameter: 347 . pc - the preconditioner context 348 349 Output Parameter: 350 . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 351 352 Notes: 353 Does not apply with 0 fill. 354 355 Level: intermediate 356 357 .seealso: `PCFactorSetAllowDiagonalFill()` 358 @*/ 359 PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) { 360 PetscFunctionBegin; 361 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 362 PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 363 PetscFunctionReturn(0); 364 } 365 366 /*@ 367 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 368 369 Logically Collective on PC 370 371 Input Parameters: 372 + pc - the preconditioner context 373 - tol - diagonal entries smaller than this in absolute value are considered zero 374 375 Options Database Key: 376 . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 377 378 Level: intermediate 379 380 .seealso: `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 381 @*/ 382 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) { 383 PetscFunctionBegin; 384 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 385 PetscValidLogicalCollectiveReal(pc, rtol, 2); 386 PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 387 PetscFunctionReturn(0); 388 } 389 390 /*@C 391 PCFactorSetMatSolverType - sets the software that is used to perform the factorization 392 393 Logically Collective on PC 394 395 Input Parameters: 396 + pc - the preconditioner context 397 - stype - for example, superlu, superlu_dist 398 399 Options Database Key: 400 . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 401 402 Level: intermediate 403 404 Note: 405 By default this will use the PETSc factorization if it exists 406 407 .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()` 408 @*/ 409 PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) { 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 412 PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 413 PetscFunctionReturn(0); 414 } 415 416 /*@C 417 PCFactorGetMatSolverType - gets the software that is used to perform the factorization 418 419 Not Collective 420 421 Input Parameter: 422 . pc - the preconditioner context 423 424 Output Parameter: 425 . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package) 426 427 Level: intermediate 428 429 .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()` 430 @*/ 431 PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) { 432 PetscErrorCode (*f)(PC, MatSolverType *); 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 436 PetscValidPointer(stype, 2); 437 PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 438 if (f) PetscCall((*f)(pc, stype)); 439 else *stype = NULL; 440 PetscFunctionReturn(0); 441 } 442 443 /*@ 444 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 445 fill = number nonzeros in factor/number nonzeros in original matrix. 446 447 Not Collective, each process can expect a different amount of fill 448 449 Input Parameters: 450 + pc - the preconditioner context 451 - fill - amount of expected fill 452 453 Options Database Key: 454 . -pc_factor_fill <fill> - Sets fill amount 455 456 Level: intermediate 457 458 Note: 459 For sparse matrix factorizations it is difficult to predict how much 460 fill to expect. By running with the option -info PETSc will print the 461 actual amount of fill used; allowing you to set the value accurately for 462 future runs. Default PETSc uses a value of 5.0 463 464 This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels. 465 466 @*/ 467 PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) { 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 470 PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0"); 471 PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 472 PetscFunctionReturn(0); 473 } 474 475 /*@ 476 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 477 For dense matrices, this enables the solution of much larger problems. 478 For sparse matrices the factorization cannot be done truly in-place 479 so this does not save memory during the factorization, but after the matrix 480 is factored, the original unfactored matrix is freed, thus recovering that 481 space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 482 483 Logically Collective on PC 484 485 Input Parameters: 486 + pc - the preconditioner context 487 - flg - PETSC_TRUE to enable, PETSC_FALSE to disable 488 489 Options Database Key: 490 . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization 491 492 Notes: 493 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 494 a different matrix is provided for the multiply and the preconditioner in 495 a call to KSPSetOperators(). 496 This is because the Krylov space methods require an application of the 497 matrix multiplication, which is not possible here because the matrix has 498 been factored in-place, replacing the original matrix. 499 500 Level: intermediate 501 502 .seealso: `PCFactorGetUseInPlace()` 503 @*/ 504 PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) { 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 507 PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 508 PetscFunctionReturn(0); 509 } 510 511 /*@ 512 PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 513 514 Logically Collective on PC 515 516 Input Parameter: 517 . pc - the preconditioner context 518 519 Output Parameter: 520 . flg - PETSC_TRUE to enable, PETSC_FALSE to disable 521 522 Level: intermediate 523 524 .seealso: `PCFactorSetUseInPlace()` 525 @*/ 526 PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) { 527 PetscFunctionBegin; 528 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 529 PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 530 PetscFunctionReturn(0); 531 } 532 533 /*@C 534 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 535 be used in the LU, ILU, Cholesky, and ICC factorizations. 536 537 Logically Collective on PC 538 539 Input Parameters: 540 + pc - the preconditioner context 541 - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM 542 543 Options Database Key: 544 . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 545 546 Level: intermediate 547 548 Notes: 549 Nested dissection is used by default for some of PETSc's sparse matrix formats 550 551 For Cholesky and ICC and the SBAIJ format the only reordering available is natural since only the upper half of the matrix is stored 552 and reordering this matrix is very expensive. 553 554 You can use a SeqAIJ matrix with Cholesky and ICC and use any ordering. 555 556 MATORDERINGEXTERNAL means PETSc will not compute an ordering and the package will use its own ordering, usable with MATSOLVERCHOLMOD, MATSOLVERUMFPACK, and others. 557 558 .seealso: `MatOrderingType` 559 560 @*/ 561 PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) { 562 PetscFunctionBegin; 563 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 564 PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 565 PetscFunctionReturn(0); 566 } 567 568 /*@ 569 PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 570 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 571 it is never done. For the MATLAB and SuperLU factorization this is used. 572 573 Logically Collective on PC 574 575 Input Parameters: 576 + pc - the preconditioner context 577 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 578 579 Options Database Key: 580 . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 581 582 Level: intermediate 583 584 .seealso: `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 585 @*/ 586 PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) { 587 PetscFunctionBegin; 588 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 589 PetscValidLogicalCollectiveReal(pc, dtcol, 2); 590 PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 591 PetscFunctionReturn(0); 592 } 593 594 /*@ 595 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 596 with BAIJ or SBAIJ matrices 597 598 Logically Collective on PC 599 600 Input Parameters: 601 + pc - the preconditioner context 602 - pivot - PETSC_TRUE or PETSC_FALSE 603 604 Options Database Key: 605 . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for BAIJ and SBAIJ 606 607 Level: intermediate 608 609 .seealso: `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 610 @*/ 611 PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) { 612 PetscFunctionBegin; 613 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 614 PetscValidLogicalCollectiveBool(pc, pivot, 2); 615 PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 616 PetscFunctionReturn(0); 617 } 618 619 /*@ 620 PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 621 this causes later ones to use the fill ratio computed in the initial factorization. 622 623 Logically Collective on PC 624 625 Input Parameters: 626 + pc - the preconditioner context 627 - flag - PETSC_TRUE to reuse else PETSC_FALSE 628 629 Options Database Key: 630 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 631 632 Level: intermediate 633 634 .seealso: `PCFactorSetReuseOrdering()` 635 @*/ 636 PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) { 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 639 PetscValidLogicalCollectiveBool(pc, flag, 2); 640 PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 641 PetscFunctionReturn(0); 642 } 643 644 PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) { 645 PC_Factor *fact = (PC_Factor *)pc->data; 646 647 PetscFunctionBegin; 648 PetscCall(MatFactorInfoInitialize(&fact->info)); 649 fact->factortype = ftype; 650 fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 651 fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 652 fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 653 fact->info.pivotinblocks = 1.0; 654 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 655 656 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 657 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 658 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 659 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 660 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 661 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 662 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 663 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 664 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 665 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 666 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 667 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 668 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 669 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 670 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 671 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 672 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 673 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 674 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 675 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 676 PetscFunctionReturn(0); 677 } 678 679 PetscErrorCode PCFactorClearComposedFunctions(PC pc) { 680 PetscFunctionBegin; 681 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 682 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 683 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 684 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 685 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 686 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 687 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 688 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 689 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 690 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 691 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 692 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 693 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 694 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 695 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 696 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 697 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 698 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 699 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 700 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 701 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 702 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 703 PetscFunctionReturn(0); 704 } 705