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