1 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 2 #include <petsc/private/matimpl.h> 3 4 /* 5 If an ordering is not yet set and the matrix is available determine a default ordering 6 */ 7 PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc) 8 { 9 PetscBool foundmtype, flg, destroy = PETSC_FALSE; 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 /* factored matrix is not present at this point, we want to create it during PCSetUp. 21 Since this can be called from setfromoptions, we destroy it when we are done with it */ 22 PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); 23 destroy = PETSC_TRUE; 24 } 25 if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS); 26 if (!fact->fact->assembled) { 27 PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg)); 28 if (!flg) { 29 Mat B; 30 PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B)); 31 PetscCall(MatHeaderReplace(fact->fact, &B)); 32 } 33 } 34 if (!fact->ordering) { 35 PetscBool canuseordering; 36 MatOrderingType otype; 37 38 PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering)); 39 if (canuseordering) PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype)); 40 else otype = MATORDERINGEXTERNAL; 41 PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering)); 42 } 43 if (destroy) PetscCall(MatDestroy(&fact->fact)); 44 } 45 } 46 PetscFunctionReturn(PETSC_SUCCESS); 47 } 48 49 static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag) 50 { 51 PC_Factor *lu = (PC_Factor *)pc->data; 52 53 PetscFunctionBegin; 54 lu->reuseordering = flag; 55 PetscFunctionReturn(PETSC_SUCCESS); 56 } 57 58 static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag) 59 { 60 PC_Factor *lu = (PC_Factor *)pc->data; 61 62 PetscFunctionBegin; 63 lu->reusefill = flag; 64 PetscFunctionReturn(PETSC_SUCCESS); 65 } 66 67 static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg) 68 { 69 PC_Factor *dir = (PC_Factor *)pc->data; 70 71 PetscFunctionBegin; 72 dir->inplace = flg; 73 PetscFunctionReturn(PETSC_SUCCESS); 74 } 75 76 static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg) 77 { 78 PC_Factor *dir = (PC_Factor *)pc->data; 79 80 PetscFunctionBegin; 81 *flg = dir->inplace; 82 PetscFunctionReturn(PETSC_SUCCESS); 83 } 84 85 /*@ 86 PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may 87 set the options for that particular factorization object. 88 89 Input Parameter: 90 . pc - the preconditioner context 91 92 Note: 93 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. 94 This function raises an error if the requested combination of solver package and matrix type is not supported. 95 96 Level: intermediate 97 98 .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 99 @*/ 100 PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 101 { 102 PetscFunctionBegin; 103 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 104 PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 105 PetscFunctionReturn(PETSC_SUCCESS); 106 } 107 108 /*@ 109 PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 110 111 Logically Collective 112 113 Input Parameters: 114 + pc - the preconditioner context 115 - zero - all pivots smaller than this will be considered zero 116 117 Options Database Key: 118 . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 119 120 Level: intermediate 121 122 .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 123 @*/ 124 PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) 125 { 126 PetscFunctionBegin; 127 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 128 PetscValidLogicalCollectiveReal(pc, zero, 2); 129 PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 130 PetscFunctionReturn(PETSC_SUCCESS); 131 } 132 133 /*@ 134 PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 135 numerical factorization, thus the matrix has nonzero pivots 136 137 Logically Collective 138 139 Input Parameters: 140 + pc - the preconditioner context 141 - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS` 142 143 Options Database Key: 144 . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 145 146 Level: intermediate 147 148 .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 149 @*/ 150 PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) 151 { 152 PetscFunctionBegin; 153 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 154 PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 155 PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 156 PetscFunctionReturn(PETSC_SUCCESS); 157 } 158 159 /*@ 160 PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 161 numerical factorization, thus the matrix has nonzero pivots 162 163 Logically Collective 164 165 Input Parameters: 166 + pc - the preconditioner context 167 - shiftamount - amount of shift or `PETSC_DECIDE` for the default 168 169 Options Database Key: 170 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 171 172 Level: intermediate 173 174 .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 175 @*/ 176 PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) 177 { 178 PetscFunctionBegin; 179 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 180 PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 181 PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 182 PetscFunctionReturn(PETSC_SUCCESS); 183 } 184 185 /*@ 186 PCFactorSetDropTolerance - The preconditioner will use an `PCILU` 187 based on a drop tolerance. 188 189 Logically Collective 190 191 Input Parameters: 192 + pc - the preconditioner context 193 . dt - the drop tolerance, try from 1.e-10 to .1 194 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 195 - maxrowcount - the max number of nonzeros allowed in a row, best value 196 depends on the number of nonzeros in row of original matrix 197 198 Options Database Key: 199 . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 200 201 Level: intermediate 202 203 Note: 204 There are NO default values for the 3 parameters, you must set them with reasonable values for your 205 matrix. We don't know how to compute reasonable values. 206 207 .seealso: [](ch_ksp), `PCILU` 208 @*/ 209 PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) 210 { 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 213 PetscValidLogicalCollectiveReal(pc, dtcol, 3); 214 PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 215 PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 216 PetscFunctionReturn(PETSC_SUCCESS); 217 } 218 219 /*@ 220 PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 221 222 Not Collective 223 224 Input Parameter: 225 . pc - the preconditioner context 226 227 Output Parameter: 228 . pivot - the tolerance 229 230 Level: intermediate 231 232 .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 233 @*/ 234 PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) 235 { 236 PetscFunctionBegin; 237 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 238 PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 239 PetscFunctionReturn(PETSC_SUCCESS); 240 } 241 242 /*@ 243 PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 244 245 Not Collective 246 247 Input Parameter: 248 . pc - the preconditioner context 249 250 Output Parameter: 251 . shift - how much to shift the diagonal entry 252 253 Level: intermediate 254 255 .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 256 @*/ 257 PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) 258 { 259 PetscFunctionBegin; 260 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 261 PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 262 PetscFunctionReturn(PETSC_SUCCESS); 263 } 264 265 /*@ 266 PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 267 268 Not Collective 269 270 Input Parameter: 271 . pc - the preconditioner context 272 273 Output Parameter: 274 . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 275 276 Level: intermediate 277 278 .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 279 @*/ 280 PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) 281 { 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 284 PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 285 PetscFunctionReturn(PETSC_SUCCESS); 286 } 287 288 /*@ 289 PCFactorGetLevels - Gets the number of levels of fill to use. 290 291 Logically Collective 292 293 Input Parameter: 294 . pc - the preconditioner context 295 296 Output Parameter: 297 . levels - number of levels of fill 298 299 Level: intermediate 300 301 .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetLevels()` 302 @*/ 303 PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) 304 { 305 PetscFunctionBegin; 306 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 307 PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 308 PetscFunctionReturn(PETSC_SUCCESS); 309 } 310 311 /*@ 312 PCFactorSetLevels - Sets the number of levels of fill to use. 313 314 Logically Collective 315 316 Input Parameters: 317 + pc - the preconditioner context 318 - levels - number of levels of fill 319 320 Options Database Key: 321 . -pc_factor_levels <levels> - Sets fill level 322 323 Level: intermediate 324 325 .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetLevels()` 326 @*/ 327 PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) 328 { 329 PetscFunctionBegin; 330 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 331 PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 332 PetscValidLogicalCollectiveInt(pc, levels, 2); 333 PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 334 PetscFunctionReturn(PETSC_SUCCESS); 335 } 336 337 /*@ 338 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 339 treated as level 0 fill even if there is no non-zero location. 340 341 Logically Collective 342 343 Input Parameters: 344 + pc - the preconditioner context 345 - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 346 347 Options Database Key: 348 . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 349 350 Note: 351 Does not apply with 0 fill. 352 353 Level: intermediate 354 355 .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 356 @*/ 357 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) 358 { 359 PetscFunctionBegin; 360 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 361 PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 362 PetscFunctionReturn(PETSC_SUCCESS); 363 } 364 365 /*@ 366 PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 367 treated as level 0 fill even if there is no non-zero location. 368 369 Logically Collective 370 371 Input Parameter: 372 . pc - the preconditioner context 373 374 Output Parameter: 375 . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 376 377 Note: 378 Does not apply with 0 fill. 379 380 Level: intermediate 381 382 .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 383 @*/ 384 PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) 385 { 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 388 PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 389 PetscFunctionReturn(PETSC_SUCCESS); 390 } 391 392 /*@ 393 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 394 395 Logically Collective 396 397 Input Parameters: 398 + pc - the preconditioner context 399 - rtol - diagonal entries smaller than this in absolute value are considered zero 400 401 Options Database Key: 402 . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 403 404 Level: intermediate 405 406 .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftAmount()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 407 @*/ 408 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) 409 { 410 PetscFunctionBegin; 411 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 412 PetscValidLogicalCollectiveReal(pc, rtol, 2); 413 PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 414 PetscFunctionReturn(PETSC_SUCCESS); 415 } 416 417 /*@ 418 PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 419 420 Logically Collective 421 422 Input Parameters: 423 + pc - the preconditioner context 424 - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 425 426 Options Database Key: 427 . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 428 429 Level: intermediate 430 431 Note: 432 The default type is set by searching for available types based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`. 433 Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver. 434 For example if one configuration had --download-mumps while a different one had --download-superlu_dist. 435 436 .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverTypeRegister()`, 437 `MatInitializePackage()`, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `MatSolverTypeGet()` 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 /*@ 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: [](ch_ksp), `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: [](ch_ksp), `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 than 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: [](ch_ksp), `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: [](ch_ksp), `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 /*@ 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: [](ch_ksp), `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: [](ch_ksp), `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: [](ch_ksp), `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: [](ch_ksp), `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