1 2 #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "PCFactorSetUpMatSolverPackage" 6 /*@ 7 PCFactorSetUpMatSolverPackage - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may 8 set the options for that particular factorization object. 9 10 Input Parameter: 11 . pc - the preconditioner context 12 13 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. 14 15 .seealso: PCFactorSetMatSolverPackage(), PCFactorGetMatrix() 16 17 Level: intermediate 18 19 @*/ 20 PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 26 ierr = PetscTryMethod(pc,"PCFactorSetUpMatSolverPackage_C",(PC),(pc));CHKERRQ(ierr); 27 PetscFunctionReturn(0); 28 } 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "PCFactorSetZeroPivot" 32 /*@ 33 PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 34 35 Logically Collective on PC 36 37 Input Parameters: 38 + pc - the preconditioner context 39 - zero - all pivots smaller than this will be considered zero 40 41 Options Database Key: 42 . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 43 44 Level: intermediate 45 46 .keywords: PC, set, factorization, direct, fill 47 48 .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount() 49 @*/ 50 PetscErrorCode PCFactorSetZeroPivot(PC pc,PetscReal zero) 51 { 52 PetscErrorCode ierr; 53 54 PetscFunctionBegin; 55 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 56 PetscValidLogicalCollectiveReal(pc,zero,2); 57 ierr = PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));CHKERRQ(ierr); 58 PetscFunctionReturn(0); 59 } 60 61 #undef __FUNCT__ 62 #define __FUNCT__ "PCFactorSetShiftType" 63 /*@ 64 PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 65 numerical factorization, thus the matrix has nonzero pivots 66 67 Logically Collective on PC 68 69 Input Parameters: 70 + pc - the preconditioner context 71 - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS 72 73 Options Database Key: 74 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 75 76 Level: intermediate 77 78 .keywords: PC, set, factorization, 79 80 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount() 81 @*/ 82 PetscErrorCode PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype) 83 { 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88 PetscValidLogicalCollectiveEnum(pc,shifttype,2); 89 ierr = PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "PCFactorSetShiftAmount" 95 /*@ 96 PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 97 numerical factorization, thus the matrix has nonzero pivots 98 99 Logically Collective on PC 100 101 Input Parameters: 102 + pc - the preconditioner context 103 - shiftamount - amount of shift 104 105 Options Database Key: 106 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 107 108 Level: intermediate 109 110 .keywords: PC, set, factorization, 111 112 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType() 113 @*/ 114 PetscErrorCode PCFactorSetShiftAmount(PC pc,PetscReal shiftamount) 115 { 116 PetscErrorCode ierr; 117 118 PetscFunctionBegin; 119 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 120 PetscValidLogicalCollectiveReal(pc,shiftamount,2); 121 ierr = PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "PCFactorSetDropTolerance" 127 /* 128 PCFactorSetDropTolerance - The preconditioner will use an ILU 129 based on a drop tolerance. (Under development) 130 131 Logically Collective on PC 132 133 Input Parameters: 134 + pc - the preconditioner context 135 . dt - the drop tolerance, try from 1.e-10 to .1 136 . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 137 - maxrowcount - the max number of nonzeros allowed in a row, best value 138 depends on the number of nonzeros in row of original matrix 139 140 Options Database Key: 141 . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 142 143 Level: intermediate 144 145 There are NO default values for the 3 parameters, you must set them with reasonable values for your 146 matrix. We don't know how to compute reasonable values. 147 148 .keywords: PC, levels, reordering, factorization, incomplete, ILU 149 */ 150 PetscErrorCode PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 151 { 152 PetscErrorCode ierr; 153 154 PetscFunctionBegin; 155 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 156 PetscValidLogicalCollectiveReal(pc,dtcol,2); 157 PetscValidLogicalCollectiveInt(pc,maxrowcount,3); 158 ierr = PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 #undef __FUNCT__ 163 #define __FUNCT__ "PCFactorSetLevels" 164 /*@ 165 PCFactorSetLevels - Sets the number of levels of fill to use. 166 167 Logically Collective on PC 168 169 Input Parameters: 170 + pc - the preconditioner context 171 - levels - number of levels of fill 172 173 Options Database Key: 174 . -pc_factor_levels <levels> - Sets fill level 175 176 Level: intermediate 177 178 .keywords: PC, levels, fill, factorization, incomplete, ILU 179 @*/ 180 PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels) 181 { 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 186 if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 187 PetscValidLogicalCollectiveInt(pc,levels,2); 188 ierr = PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 194 /*@ 195 PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 196 treated as level 0 fill even if there is no non-zero location. 197 198 Logically Collective on PC 199 200 Input Parameters: 201 + pc - the preconditioner context 202 203 Options Database Key: 204 . -pc_factor_diagonal_fill 205 206 Notes: 207 Does not apply with 0 fill. 208 209 Level: intermediate 210 211 .keywords: PC, levels, fill, factorization, incomplete, ILU 212 @*/ 213 PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 219 ierr = PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC),(pc));CHKERRQ(ierr); 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 225 /*@ 226 PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 227 228 Logically Collective on PC 229 230 Input Parameters: 231 + pc - the preconditioner context 232 - tol - diagonal entries smaller than this in absolute value are considered zero 233 234 Options Database Key: 235 . -pc_factor_nonzeros_along_diagonal 236 237 Level: intermediate 238 239 .keywords: PC, set, factorization, direct, fill 240 241 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 242 @*/ 243 PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 244 { 245 PetscErrorCode ierr; 246 247 PetscFunctionBegin; 248 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 249 PetscValidLogicalCollectiveReal(pc,rtol,2); 250 ierr = PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "PCFactorSetMatSolverPackage" 256 /*@C 257 PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 258 259 Logically Collective on PC 260 261 Input Parameters: 262 + pc - the preconditioner context 263 - stype - for example, superlu, superlu_dist 264 265 Options Database Key: 266 . -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse 267 268 Level: intermediate 269 270 Note: 271 By default this will use the PETSc factorization if it exists 272 273 274 .keywords: PC, set, factorization, direct, fill 275 276 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 277 278 @*/ 279 PetscErrorCode PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 280 { 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 285 ierr = PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "PCFactorGetMatSolverPackage" 291 /*@C 292 PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization 293 294 Not Collective 295 296 Input Parameter: 297 . pc - the preconditioner context 298 299 Output Parameter: 300 . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package) 301 302 Level: intermediate 303 304 305 .keywords: PC, set, factorization, direct, fill 306 307 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 308 309 @*/ 310 PetscErrorCode PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype) 311 { 312 PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*); 313 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 316 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);CHKERRQ(ierr); 317 if (f) { 318 ierr = (*f)(pc,stype);CHKERRQ(ierr); 319 } else { 320 *stype = NULL; 321 } 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "PCFactorSetFill" 327 /*@ 328 PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 329 fill = number nonzeros in factor/number nonzeros in original matrix. 330 331 Not Collective, each process can expect a different amount of fill 332 333 Input Parameters: 334 + pc - the preconditioner context 335 - fill - amount of expected fill 336 337 Options Database Key: 338 . -pc_factor_fill <fill> - Sets fill amount 339 340 Level: intermediate 341 342 Note: 343 For sparse matrix factorizations it is difficult to predict how much 344 fill to expect. By running with the option -info PETSc will print the 345 actual amount of fill used; allowing you to set the value accurately for 346 future runs. Default PETSc uses a value of 5.0 347 348 This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels. 349 350 351 .keywords: PC, set, factorization, direct, fill 352 353 @*/ 354 PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill) 355 { 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 360 if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 361 ierr = PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 #undef __FUNCT__ 366 #define __FUNCT__ "PCFactorSetUseInPlace" 367 /*@ 368 PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 369 For dense matrices, this enables the solution of much larger problems. 370 For sparse matrices the factorization cannot be done truly in-place 371 so this does not save memory during the factorization, but after the matrix 372 is factored, the original unfactored matrix is freed, thus recovering that 373 space. 374 375 Logically Collective on PC 376 377 Input Parameters: 378 . pc - the preconditioner context 379 380 Options Database Key: 381 . -pc_factor_in_place - Activates in-place factorization 382 383 Notes: 384 PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 385 a different matrix is provided for the multiply and the preconditioner in 386 a call to KSPSetOperators(). 387 This is because the Krylov space methods require an application of the 388 matrix multiplication, which is not possible here because the matrix has 389 been factored in-place, replacing the original matrix. 390 391 Level: intermediate 392 393 .keywords: PC, set, factorization, direct, inplace, in-place, LU 394 395 .seealso: PCILUSetUseInPlace() 396 @*/ 397 PetscErrorCode PCFactorSetUseInPlace(PC pc) 398 { 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 403 ierr = PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC),(pc));CHKERRQ(ierr); 404 PetscFunctionReturn(0); 405 } 406 407 #undef __FUNCT__ 408 #define __FUNCT__ "PCFactorSetMatOrderingType" 409 /*@C 410 PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 411 be used in the LU factorization. 412 413 Logically Collective on PC 414 415 Input Parameters: 416 + pc - the preconditioner context 417 - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM 418 419 Options Database Key: 420 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 421 422 Level: intermediate 423 424 Notes: nested dissection is used by default 425 426 For Cholesky and ICC and the SBAIJ format reorderings are not available, 427 since only the upper triangular part of the matrix is stored. You can use the 428 SeqAIJ format in this case to get reorderings. 429 430 @*/ 431 PetscErrorCode PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 432 { 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 437 ierr = PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 #undef __FUNCT__ 442 #define __FUNCT__ "PCFactorSetColumnPivot" 443 /*@ 444 PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 445 For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 446 it is never done. For the MATLAB and SuperLU factorization this is used. 447 448 Logically Collective on PC 449 450 Input Parameters: 451 + pc - the preconditioner context 452 - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 453 454 Options Database Key: 455 . -pc_factor_pivoting <dtcol> 456 457 Level: intermediate 458 459 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 460 @*/ 461 PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol) 462 { 463 PetscErrorCode ierr; 464 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 467 PetscValidLogicalCollectiveReal(pc,dtcol,2); 468 ierr = PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));CHKERRQ(ierr); 469 PetscFunctionReturn(0); 470 } 471 472 #undef __FUNCT__ 473 #define __FUNCT__ "PCFactorSetPivotInBlocks" 474 /*@ 475 PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 476 with BAIJ or SBAIJ matrices 477 478 Logically Collective on PC 479 480 Input Parameters: 481 + pc - the preconditioner context 482 - pivot - PETSC_TRUE or PETSC_FALSE 483 484 Options Database Key: 485 . -pc_factor_pivot_in_blocks <true,false> 486 487 Level: intermediate 488 489 .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot() 490 @*/ 491 PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot) 492 { 493 PetscErrorCode ierr; 494 495 PetscFunctionBegin; 496 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 497 PetscValidLogicalCollectiveBool(pc,pivot,2); 498 ierr = PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));CHKERRQ(ierr); 499 PetscFunctionReturn(0); 500 } 501 502 #undef __FUNCT__ 503 #define __FUNCT__ "PCFactorSetReuseFill" 504 /*@ 505 PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 506 this causes later ones to use the fill ratio computed in the initial factorization. 507 508 Logically Collective on PC 509 510 Input Parameters: 511 + pc - the preconditioner context 512 - flag - PETSC_TRUE to reuse else PETSC_FALSE 513 514 Options Database Key: 515 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 516 517 Level: intermediate 518 519 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 520 521 .seealso: PCFactorSetReuseOrdering() 522 @*/ 523 PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag) 524 { 525 PetscErrorCode ierr; 526 527 PetscFunctionBegin; 528 PetscValidHeaderSpecific(pc,PC_CLASSID,2); 529 PetscValidLogicalCollectiveBool(pc,flag,2); 530 ierr = PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 531 PetscFunctionReturn(0); 532 } 533