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