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