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