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