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