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