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 Note: 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: `PCCHOLESKY`, `PCLU`, `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: `PCCHOLESKY`, `PCLU`, `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: `PCCHOLESKY`, `PCLU`, `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 or `PETSC_DECIDE` for the default 155 156 Options Database Key: 157 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 158 159 Level: intermediate 160 161 .seealso: `PCCHOLESKY`, `PCLU`, ``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 `PCILU` 173 based on a drop tolerance. 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 Note: 190 There are NO default values for the 3 parameters, you must set them with reasonable values for your 191 matrix. We don't know how to compute reasonable values. 192 193 .seealso: `PCILU` 194 @*/ 195 PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) { 196 PetscFunctionBegin; 197 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 198 PetscValidLogicalCollectiveReal(pc, dtcol, 3); 199 PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 200 PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 201 PetscFunctionReturn(0); 202 } 203 204 /*@ 205 PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 206 207 Not Collective 208 209 Input Parameters: 210 . pc - the preconditioner context 211 212 Output Parameter: 213 . pivot - the tolerance 214 215 Level: intermediate 216 217 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 218 @*/ 219 PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) { 220 PetscFunctionBegin; 221 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 222 PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 223 PetscFunctionReturn(0); 224 } 225 226 /*@ 227 PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 228 229 Not Collective 230 231 Input Parameters: 232 . pc - the preconditioner context 233 234 Output Parameter: 235 . shift - how much to shift the diagonal entry 236 237 Level: intermediate 238 239 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 240 @*/ 241 PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) { 242 PetscFunctionBegin; 243 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 244 PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 245 PetscFunctionReturn(0); 246 } 247 248 /*@ 249 PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 250 251 Not Collective 252 253 Input Parameter: 254 . pc - the preconditioner context 255 256 Output Parameter: 257 . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 258 259 Level: intermediate 260 261 .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 262 @*/ 263 PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) { 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 266 PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 267 PetscFunctionReturn(0); 268 } 269 270 /*@ 271 PCFactorGetLevels - Gets the number of levels of fill to use. 272 273 Logically Collective on pc 274 275 Input Parameters: 276 . pc - the preconditioner context 277 278 Output Parameter: 279 . levels - number of levels of fill 280 281 Level: intermediate 282 283 .seealso: `PCILU`, `PCICC`, `PCFactorSetLevels()` 284 @*/ 285 PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) { 286 PetscFunctionBegin; 287 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 288 PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 289 PetscFunctionReturn(0); 290 } 291 292 /*@ 293 PCFactorSetLevels - Sets the number of levels of fill to use. 294 295 Logically Collective on pc 296 297 Input Parameters: 298 + pc - the preconditioner context 299 - levels - number of levels of fill 300 301 Options Database Key: 302 . -pc_factor_levels <levels> - Sets fill level 303 304 Level: intermediate 305 306 .seealso: `PCILU`, `PCICC`, `PCFactorGetLevels()` 307 @*/ 308 PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) { 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 311 PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 312 PetscValidLogicalCollectiveInt(pc, levels, 2); 313 PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 314 PetscFunctionReturn(0); 315 } 316 317 /*@ 318 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 319 treated as level 0 fill even if there is no non-zero location. 320 321 Logically Collective on pc 322 323 Input Parameters: 324 + pc - the preconditioner context 325 - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 326 327 Options Database Key: 328 . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 329 330 Note: 331 Does not apply with 0 fill. 332 333 Level: intermediate 334 335 .seealso: `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 336 @*/ 337 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) { 338 PetscFunctionBegin; 339 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 340 PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 341 PetscFunctionReturn(0); 342 } 343 344 /*@ 345 PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 346 treated as level 0 fill even if there is no non-zero location. 347 348 Logically Collective on pc 349 350 Input Parameter: 351 . pc - the preconditioner context 352 353 Output Parameter: 354 . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 355 356 Note: 357 Does not apply with 0 fill. 358 359 Level: intermediate 360 361 .seealso: `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 362 @*/ 363 PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) { 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 366 PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 367 PetscFunctionReturn(0); 368 } 369 370 /*@ 371 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 372 373 Logically Collective on pc 374 375 Input Parameters: 376 + pc - the preconditioner context 377 - tol - diagonal entries smaller than this in absolute value are considered zero 378 379 Options Database Key: 380 . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 381 382 Level: intermediate 383 384 .seealso: `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 385 @*/ 386 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) { 387 PetscFunctionBegin; 388 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 389 PetscValidLogicalCollectiveReal(pc, rtol, 2); 390 PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 391 PetscFunctionReturn(0); 392 } 393 394 /*@C 395 PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 396 397 Logically Collective on pc 398 399 Input Parameters: 400 + pc - the preconditioner context 401 - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 402 403 Options Database Key: 404 . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 405 406 Level: intermediate 407 408 Note: 409 By default this will use the PETSc factorization if it exists 410 411 .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 412 `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 413 @*/ 414 PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) { 415 PetscFunctionBegin; 416 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 417 PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 418 PetscFunctionReturn(0); 419 } 420 421 /*@C 422 PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization 423 424 Not Collective 425 426 Input Parameter: 427 . pc - the preconditioner context 428 429 Output Parameter: 430 . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 431 432 Level: intermediate 433 434 .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 435 `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 436 @*/ 437 PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) { 438 PetscErrorCode (*f)(PC, MatSolverType *); 439 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 442 PetscValidPointer(stype, 2); 443 PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 444 if (f) PetscCall((*f)(pc, stype)); 445 else *stype = NULL; 446 PetscFunctionReturn(0); 447 } 448 449 /*@ 450 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 451 fill = number nonzeros in factor/number nonzeros in original matrix. 452 453 Not Collective, each process can expect a different amount of fill 454 455 Input Parameters: 456 + pc - the preconditioner context 457 - fill - amount of expected fill 458 459 Options Database Key: 460 . -pc_factor_fill <fill> - Sets fill amount 461 462 Level: intermediate 463 464 Notes: 465 For sparse matrix factorizations it is difficult to predict how much 466 fill to expect. By running with the option -info PETSc will print the 467 actual amount of fill used; allowing you to set the value accurately for 468 future runs. Default PETSc uses a value of 5.0 469 470 This is ignored for most solver packages 471 472 This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels. 473 474 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()` 475 @*/ 476 PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) { 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 479 PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0"); 480 PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 481 PetscFunctionReturn(0); 482 } 483 484 /*@ 485 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 486 For dense matrices, this enables the solution of much larger problems. 487 For sparse matrices the factorization cannot be done truly in-place 488 so this does not save memory during the factorization, but after the matrix 489 is factored, the original unfactored matrix is freed, thus recovering that 490 space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 491 492 Logically Collective on pc 493 494 Input Parameters: 495 + pc - the preconditioner context 496 - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 497 498 Options Database Key: 499 . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization 500 501 Note: 502 `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when 503 a different matrix is provided for the multiply and the preconditioner in 504 a call to `KSPSetOperators()`. 505 This is because the Krylov space methods require an application of the 506 matrix multiplication, which is not possible here because the matrix has 507 been factored in-place, replacing the original matrix. 508 509 Level: intermediate 510 511 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()` 512 @*/ 513 PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) { 514 PetscFunctionBegin; 515 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 516 PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 517 PetscFunctionReturn(0); 518 } 519 520 /*@ 521 PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 522 523 Logically Collective on pc 524 525 Input Parameter: 526 . pc - the preconditioner context 527 528 Output Parameter: 529 . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 530 531 Level: intermediate 532 533 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()` 534 @*/ 535 PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) { 536 PetscFunctionBegin; 537 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 538 PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 539 PetscFunctionReturn(0); 540 } 541 542 /*@C 543 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 544 be used in the `PCLU`, `PCCHOLESKY`, `PCILU`, or `PCICC` preconditioners 545 546 Logically Collective on pc 547 548 Input Parameters: 549 + pc - the preconditioner context 550 - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM` 551 552 Options Database Key: 553 . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 554 555 Level: intermediate 556 557 Notes: 558 Nested dissection is used by default for some of PETSc's sparse matrix formats 559 560 For `PCCHOLESKY` and `PCICC` and the `MATSBAIJ` format the only reordering available is natural since only the upper half of the matrix is stored 561 and reordering this matrix is very expensive. 562 563 You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering. 564 565 `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others. 566 567 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM` 568 @*/ 569 PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) { 570 PetscFunctionBegin; 571 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 572 PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 573 PetscFunctionReturn(0); 574 } 575 576 /*@ 577 PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 578 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 579 it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used. 580 581 Logically Collective on pc 582 583 Input Parameters: 584 + pc - the preconditioner context 585 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 586 587 Options Database Key: 588 . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 589 590 Level: intermediate 591 592 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 593 @*/ 594 PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) { 595 PetscFunctionBegin; 596 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 597 PetscValidLogicalCollectiveReal(pc, dtcol, 2); 598 PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 599 PetscFunctionReturn(0); 600 } 601 602 /*@ 603 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 604 with `MATBAIJ` or `MATSBAIJ` matrices 605 606 Logically Collective on pc 607 608 Input Parameters: 609 + pc - the preconditioner context 610 - pivot - `PETSC_TRUE` or `PETSC_FALSE` 611 612 Options Database Key: 613 . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ` 614 615 Level: intermediate 616 617 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 618 @*/ 619 PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) { 620 PetscFunctionBegin; 621 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 622 PetscValidLogicalCollectiveBool(pc, pivot, 2); 623 PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 624 PetscFunctionReturn(0); 625 } 626 627 /*@ 628 PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 629 this causes later ones to use the fill ratio computed in the initial factorization. 630 631 Logically Collective on pc 632 633 Input Parameters: 634 + pc - the preconditioner context 635 - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 636 637 Options Database Key: 638 . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 639 640 Level: intermediate 641 642 .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()` 643 @*/ 644 PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) { 645 PetscFunctionBegin; 646 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 647 PetscValidLogicalCollectiveBool(pc, flag, 2); 648 PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 649 PetscFunctionReturn(0); 650 } 651 652 PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) { 653 PC_Factor *fact = (PC_Factor *)pc->data; 654 655 PetscFunctionBegin; 656 PetscCall(MatFactorInfoInitialize(&fact->info)); 657 fact->factortype = ftype; 658 fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 659 fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 660 fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 661 fact->info.pivotinblocks = 1.0; 662 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 663 664 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 665 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 666 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 667 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 668 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 669 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 670 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 671 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 672 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 673 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 674 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 675 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 676 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 677 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 678 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 679 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 680 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 681 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 682 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 683 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 684 PetscFunctionReturn(0); 685 } 686 687 PetscErrorCode PCFactorClearComposedFunctions(PC pc) { 688 PetscFunctionBegin; 689 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 690 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 691 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 692 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 693 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 694 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 695 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 696 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 697 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 698 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 699 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 700 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 701 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 702 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 703 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 704 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 705 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 706 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 707 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 708 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 709 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 710 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 711 PetscFunctionReturn(0); 712 } 713