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