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