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