1dba47a55SKris Buschelman #define PETSCKSP_DLL 25e8efad8SHong Zhang 3e1222f9fSBarry Smith #include "src/ksp/pc/impls/factor/factor.h" /*I "petscpc.h" I*/ 45e8efad8SHong Zhang 5ee45ca4aSHong Zhang #undef __FUNCT__ 6ee45ca4aSHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot" 7ee45ca4aSHong Zhang /*@ 8ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 9ee45ca4aSHong Zhang 10ee45ca4aSHong Zhang Collective on PC 11ee45ca4aSHong Zhang 12ee45ca4aSHong Zhang Input Parameters: 13afaefe49SHong Zhang + pc - the preconditioner context 14afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 15ee45ca4aSHong Zhang 16ee45ca4aSHong Zhang Options Database Key: 17ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 18ee45ca4aSHong Zhang 19ee45ca4aSHong Zhang Level: intermediate 20ee45ca4aSHong Zhang 21ee45ca4aSHong Zhang .keywords: PC, set, factorization, direct, fill 22ee45ca4aSHong Zhang 23ee45ca4aSHong Zhang .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd() 24ee45ca4aSHong Zhang @*/ 25dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero) 26ee45ca4aSHong Zhang { 27afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 28afaefe49SHong Zhang 29ee45ca4aSHong Zhang PetscFunctionBegin; 30afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 31afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 32afaefe49SHong Zhang if (f) { 33afaefe49SHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 34afaefe49SHong Zhang } 35ee45ca4aSHong Zhang PetscFunctionReturn(0); 36ee45ca4aSHong Zhang } 37ee45ca4aSHong Zhang 385e8efad8SHong Zhang #undef __FUNCT__ 395e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero" 405e8efad8SHong Zhang /*@ 415e8efad8SHong Zhang PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during 425e8efad8SHong Zhang numerical factorization, thus the matrix has nonzero pivots 435e8efad8SHong Zhang 445e8efad8SHong Zhang Collective on PC 455e8efad8SHong Zhang 465e8efad8SHong Zhang Input Parameters: 47afaefe49SHong Zhang + pc - the preconditioner context 48afaefe49SHong Zhang - shift - amount of shift 495e8efad8SHong Zhang 505e8efad8SHong Zhang Options Database Key: 519f95998fSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 525e8efad8SHong Zhang 535e8efad8SHong Zhang Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero 545e8efad8SHong Zhang pivot, then the shift is doubled until this is alleviated. 555e8efad8SHong Zhang 565e8efad8SHong Zhang Level: intermediate 575e8efad8SHong Zhang 585e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill 595e8efad8SHong Zhang 60ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd() 615e8efad8SHong Zhang @*/ 62dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift) 635e8efad8SHong Zhang { 64afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 65afaefe49SHong Zhang 665e8efad8SHong Zhang PetscFunctionBegin; 67afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 68afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr); 69afaefe49SHong Zhang if (f) { 70afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 715e8efad8SHong Zhang } 725e8efad8SHong Zhang PetscFunctionReturn(0); 735e8efad8SHong Zhang } 740a29876aSHong Zhang 750a29876aSHong Zhang #undef __FUNCT__ 760a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd" 770a29876aSHong Zhang /*@ 780a29876aSHong Zhang PCFactorSetShiftPd - specify whether to use Manteuffel shifting. 790a29876aSHong Zhang If a matrix factorisation breaks down because of nonpositive pivots, 800a29876aSHong Zhang adding sufficient identity to the diagonal will remedy this. 810a29876aSHong Zhang Setting this causes a bisection method to find the minimum shift that 820a29876aSHong Zhang will lead to a well-defined matrix factor. 830a29876aSHong Zhang 840a29876aSHong Zhang Collective on PC 850a29876aSHong Zhang 860a29876aSHong Zhang Input parameters: 87afaefe49SHong Zhang + pc - the preconditioner context 88afaefe49SHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 890a29876aSHong Zhang 900a29876aSHong Zhang Options Database Key: 919f95998fSHong Zhang . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 929f95998fSHong Zhang is optional with PETSC_TRUE being the default 930a29876aSHong Zhang 940a29876aSHong Zhang Level: intermediate 950a29876aSHong Zhang 960a29876aSHong Zhang .keywords: PC, indefinite, factorization 970a29876aSHong Zhang 98ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero() 990a29876aSHong Zhang @*/ 100dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift) 1010a29876aSHong Zhang { 102afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 103afaefe49SHong Zhang 1040a29876aSHong Zhang PetscFunctionBegin; 105afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 106afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr); 107afaefe49SHong Zhang if (f) { 108afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 109afaefe49SHong Zhang } 1100a29876aSHong Zhang PetscFunctionReturn(0); 1110a29876aSHong Zhang } 112*85317021SBarry Smith 113*85317021SBarry Smith #undef __FUNCT__ 114*85317021SBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance" 115*85317021SBarry Smith /*@ 116*85317021SBarry Smith PCFactorSetUseDropTolerance - The preconditioner will use an ILU 117*85317021SBarry Smith based on a drop tolerance. 118*85317021SBarry Smith 119*85317021SBarry Smith Collective on PC 120*85317021SBarry Smith 121*85317021SBarry Smith Input Parameters: 122*85317021SBarry Smith + pc - the preconditioner context 123*85317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 124*85317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 125*85317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 126*85317021SBarry Smith depends on the number of nonzeros in row of original matrix 127*85317021SBarry Smith 128*85317021SBarry Smith Options Database Key: 129*85317021SBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 130*85317021SBarry Smith 131*85317021SBarry Smith Level: intermediate 132*85317021SBarry Smith 133*85317021SBarry Smith Notes: 134*85317021SBarry Smith This uses the iludt() code of Saad's SPARSKIT package 135*85317021SBarry Smith 136*85317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 137*85317021SBarry Smith matrix. We don't know how to compute reasonable values. 138*85317021SBarry Smith 139*85317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU 140*85317021SBarry Smith @*/ 141*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 142*85317021SBarry Smith { 143*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 144*85317021SBarry Smith 145*85317021SBarry Smith PetscFunctionBegin; 146*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 147*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 148*85317021SBarry Smith if (f) { 149*85317021SBarry Smith ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 150*85317021SBarry Smith } 151*85317021SBarry Smith PetscFunctionReturn(0); 152*85317021SBarry Smith } 153*85317021SBarry Smith 154*85317021SBarry Smith #undef __FUNCT__ 155*85317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels" 156*85317021SBarry Smith /*@ 157*85317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 158*85317021SBarry Smith 159*85317021SBarry Smith Collective on PC 160*85317021SBarry Smith 161*85317021SBarry Smith Input Parameters: 162*85317021SBarry Smith + pc - the preconditioner context 163*85317021SBarry Smith - levels - number of levels of fill 164*85317021SBarry Smith 165*85317021SBarry Smith Options Database Key: 166*85317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 167*85317021SBarry Smith 168*85317021SBarry Smith Level: intermediate 169*85317021SBarry Smith 170*85317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU 171*85317021SBarry Smith @*/ 172*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 173*85317021SBarry Smith { 174*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 175*85317021SBarry Smith 176*85317021SBarry Smith PetscFunctionBegin; 177*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 178*85317021SBarry Smith if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 179*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 180*85317021SBarry Smith if (f) { 181*85317021SBarry Smith ierr = (*f)(pc,levels);CHKERRQ(ierr); 182*85317021SBarry Smith } 183*85317021SBarry Smith PetscFunctionReturn(0); 184*85317021SBarry Smith } 185*85317021SBarry Smith 186*85317021SBarry Smith #undef __FUNCT__ 187*85317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 188*85317021SBarry Smith /*@ 189*85317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 190*85317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 191*85317021SBarry Smith 192*85317021SBarry Smith Collective on PC 193*85317021SBarry Smith 194*85317021SBarry Smith Input Parameters: 195*85317021SBarry Smith + pc - the preconditioner context 196*85317021SBarry Smith 197*85317021SBarry Smith Options Database Key: 198*85317021SBarry Smith . -pc_factor_diagonal_fill 199*85317021SBarry Smith 200*85317021SBarry Smith Notes: 201*85317021SBarry Smith Does not apply with 0 fill. 202*85317021SBarry Smith 203*85317021SBarry Smith Level: intermediate 204*85317021SBarry Smith 205*85317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU 206*85317021SBarry Smith @*/ 207*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 208*85317021SBarry Smith { 209*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC); 210*85317021SBarry Smith 211*85317021SBarry Smith PetscFunctionBegin; 212*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 213*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 214*85317021SBarry Smith if (f) { 215*85317021SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 216*85317021SBarry Smith } 217*85317021SBarry Smith PetscFunctionReturn(0); 218*85317021SBarry Smith } 219*85317021SBarry Smith 220*85317021SBarry Smith #undef __FUNCT__ 221*85317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks" 222*85317021SBarry Smith /*@ 223*85317021SBarry Smith PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot. 224*85317021SBarry Smith 225*85317021SBarry Smith Collective on PC 226*85317021SBarry Smith 227*85317021SBarry Smith Input Parameters: 228*85317021SBarry Smith + pc - the preconditioner context 229*85317021SBarry Smith - shift - amount to shift or PETSC_DECIDE 230*85317021SBarry Smith 231*85317021SBarry Smith Options Database Key: 232*85317021SBarry Smith . -pc_factor_shift_in_blocks <shift> 233*85317021SBarry Smith 234*85317021SBarry Smith Level: intermediate 235*85317021SBarry Smith 236*85317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks() 237*85317021SBarry Smith @*/ 238*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift) 239*85317021SBarry Smith { 240*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 241*85317021SBarry Smith 242*85317021SBarry Smith PetscFunctionBegin; 243*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 244*85317021SBarry Smith if (f) { 245*85317021SBarry Smith ierr = (*f)(pc,shift);CHKERRQ(ierr); 246*85317021SBarry Smith } 247*85317021SBarry Smith PetscFunctionReturn(0); 248*85317021SBarry Smith } 249*85317021SBarry Smith 250*85317021SBarry Smith #undef __FUNCT__ 251*85317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 252*85317021SBarry Smith /*@ 253*85317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 254*85317021SBarry Smith 255*85317021SBarry Smith Collective on PC 256*85317021SBarry Smith 257*85317021SBarry Smith Input Parameters: 258*85317021SBarry Smith + pc - the preconditioner context 259*85317021SBarry Smith - tol - diagonal entries smaller than this in absolute value are considered zero 260*85317021SBarry Smith 261*85317021SBarry Smith Options Database Key: 262*85317021SBarry Smith . -pc_factor_nonzeros_along_diagonal 263*85317021SBarry Smith 264*85317021SBarry Smith Level: intermediate 265*85317021SBarry Smith 266*85317021SBarry Smith .keywords: PC, set, factorization, direct, fill 267*85317021SBarry Smith 268*85317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 269*85317021SBarry Smith @*/ 270*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 271*85317021SBarry Smith { 272*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 273*85317021SBarry Smith 274*85317021SBarry Smith PetscFunctionBegin; 275*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 276*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 277*85317021SBarry Smith if (f) { 278*85317021SBarry Smith ierr = (*f)(pc,rtol);CHKERRQ(ierr); 279*85317021SBarry Smith } 280*85317021SBarry Smith PetscFunctionReturn(0); 281*85317021SBarry Smith } 282*85317021SBarry Smith 283*85317021SBarry Smith #undef __FUNCT__ 284*85317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage" 285*85317021SBarry Smith /*@ 286*85317021SBarry Smith PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 287*85317021SBarry Smith 288*85317021SBarry Smith Collective on PC 289*85317021SBarry Smith 290*85317021SBarry Smith Input Parameters: 291*85317021SBarry Smith + pc - the preconditioner context 292*85317021SBarry Smith - stype - for example, spooles, superlu, superlu_dist 293*85317021SBarry Smith 294*85317021SBarry Smith Options Database Key: 295*85317021SBarry Smith . -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps 296*85317021SBarry Smith 297*85317021SBarry Smith Level: intermediate 298*85317021SBarry Smith 299*85317021SBarry Smith Note: 300*85317021SBarry Smith By default this will use the PETSc factorization if it exists 301*85317021SBarry Smith 302*85317021SBarry Smith Use PCFactorGetMatrix(pc,&mat); MatFactorGetPackage(mat,&package); to see what package is being used 303*85317021SBarry Smith 304*85317021SBarry Smith 305*85317021SBarry Smith .keywords: PC, set, factorization, direct, fill 306*85317021SBarry Smith 307*85317021SBarry Smith .seealso: MatGetFactor(), MatSolverPackage 308*85317021SBarry Smith 309*85317021SBarry Smith @*/ 310*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 311*85317021SBarry Smith { 312*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,const MatSolverPackage); 313*85317021SBarry Smith 314*85317021SBarry Smith PetscFunctionBegin; 315*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 316*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr); 317*85317021SBarry Smith if (f) { 318*85317021SBarry Smith ierr = (*f)(pc,stype);CHKERRQ(ierr); 319*85317021SBarry Smith } 320*85317021SBarry Smith PetscFunctionReturn(0); 321*85317021SBarry Smith } 322*85317021SBarry Smith 323*85317021SBarry Smith #undef __FUNCT__ 324*85317021SBarry Smith #define __FUNCT__ "PCFactorSetFill" 325*85317021SBarry Smith /*@ 326*85317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 327*85317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 328*85317021SBarry Smith 329*85317021SBarry Smith Collective on PC 330*85317021SBarry Smith 331*85317021SBarry Smith Input Parameters: 332*85317021SBarry Smith + pc - the preconditioner context 333*85317021SBarry Smith - fill - amount of expected fill 334*85317021SBarry Smith 335*85317021SBarry Smith Options Database Key: 336*85317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 337*85317021SBarry Smith 338*85317021SBarry Smith Level: intermediate 339*85317021SBarry Smith 340*85317021SBarry Smith Note: 341*85317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 342*85317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 343*85317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 344*85317021SBarry Smith future runs. Default PETSc uses a value of 5.0 345*85317021SBarry Smith 346*85317021SBarry Smith .keywords: PC, set, factorization, direct, fill 347*85317021SBarry Smith 348*85317021SBarry Smith @*/ 349*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 350*85317021SBarry Smith { 351*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 352*85317021SBarry Smith 353*85317021SBarry Smith PetscFunctionBegin; 354*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 355*85317021SBarry Smith if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 356*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 357*85317021SBarry Smith if (f) { 358*85317021SBarry Smith ierr = (*f)(pc,fill);CHKERRQ(ierr); 359*85317021SBarry Smith } 360*85317021SBarry Smith PetscFunctionReturn(0); 361*85317021SBarry Smith } 362*85317021SBarry Smith 363*85317021SBarry Smith #undef __FUNCT__ 364*85317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace" 365*85317021SBarry Smith /*@ 366*85317021SBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 367*85317021SBarry Smith For dense matrices, this enables the solution of much larger problems. 368*85317021SBarry Smith For sparse matrices the factorization cannot be done truly in-place 369*85317021SBarry Smith so this does not save memory during the factorization, but after the matrix 370*85317021SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 371*85317021SBarry Smith space. 372*85317021SBarry Smith 373*85317021SBarry Smith Collective on PC 374*85317021SBarry Smith 375*85317021SBarry Smith Input Parameters: 376*85317021SBarry Smith . pc - the preconditioner context 377*85317021SBarry Smith 378*85317021SBarry Smith Options Database Key: 379*85317021SBarry Smith . -pc_factor_in_place - Activates in-place factorization 380*85317021SBarry Smith 381*85317021SBarry Smith Notes: 382*85317021SBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 383*85317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 384*85317021SBarry Smith a call to KSPSetOperators(). 385*85317021SBarry Smith This is because the Krylov space methods require an application of the 386*85317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 387*85317021SBarry Smith been factored in-place, replacing the original matrix. 388*85317021SBarry Smith 389*85317021SBarry Smith Level: intermediate 390*85317021SBarry Smith 391*85317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU 392*85317021SBarry Smith 393*85317021SBarry Smith .seealso: PCILUSetUseInPlace() 394*85317021SBarry Smith @*/ 395*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 396*85317021SBarry Smith { 397*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC); 398*85317021SBarry Smith 399*85317021SBarry Smith PetscFunctionBegin; 400*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 401*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 402*85317021SBarry Smith if (f) { 403*85317021SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 404*85317021SBarry Smith } 405*85317021SBarry Smith PetscFunctionReturn(0); 406*85317021SBarry Smith } 407*85317021SBarry Smith 408*85317021SBarry Smith #undef __FUNCT__ 409*85317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType" 410*85317021SBarry Smith /*@C 411*85317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 412*85317021SBarry Smith be used in the LU factorization. 413*85317021SBarry Smith 414*85317021SBarry Smith Collective on PC 415*85317021SBarry Smith 416*85317021SBarry Smith Input Parameters: 417*85317021SBarry Smith + pc - the preconditioner context 418*85317021SBarry Smith - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 419*85317021SBarry Smith 420*85317021SBarry Smith Options Database Key: 421*85317021SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 422*85317021SBarry Smith 423*85317021SBarry Smith Level: intermediate 424*85317021SBarry Smith 425*85317021SBarry Smith Notes: nested dissection is used by default 426*85317021SBarry Smith 427*85317021SBarry Smith For Cholesky and ICC and the SBAIJ format reorderings are not available, 428*85317021SBarry Smith since only the upper triangular part of the matrix is stored. You can use the 429*85317021SBarry Smith SeqAIJ format in this case to get reorderings. 430*85317021SBarry Smith 431*85317021SBarry Smith @*/ 432*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering) 433*85317021SBarry Smith { 434*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,const MatOrderingType); 435*85317021SBarry Smith 436*85317021SBarry Smith PetscFunctionBegin; 437*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 438*85317021SBarry Smith if (f) { 439*85317021SBarry Smith ierr = (*f)(pc,ordering);CHKERRQ(ierr); 440*85317021SBarry Smith } 441*85317021SBarry Smith PetscFunctionReturn(0); 442*85317021SBarry Smith } 443*85317021SBarry Smith 444*85317021SBarry Smith #undef __FUNCT__ 445*85317021SBarry Smith #define __FUNCT__ "PCFactorSetPivoting" 446*85317021SBarry Smith /*@ 447*85317021SBarry Smith PCFactorSetPivoting - Determines when pivoting is done during LU. 448*85317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 449*85317021SBarry Smith it is never done. For the Matlab and SuperLU factorization this is used. 450*85317021SBarry Smith 451*85317021SBarry Smith Collective on PC 452*85317021SBarry Smith 453*85317021SBarry Smith Input Parameters: 454*85317021SBarry Smith + pc - the preconditioner context 455*85317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 456*85317021SBarry Smith 457*85317021SBarry Smith Options Database Key: 458*85317021SBarry Smith . -pc_factor_pivoting <dtcol> 459*85317021SBarry Smith 460*85317021SBarry Smith Level: intermediate 461*85317021SBarry Smith 462*85317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 463*85317021SBarry Smith @*/ 464*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 465*85317021SBarry Smith { 466*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 467*85317021SBarry Smith 468*85317021SBarry Smith PetscFunctionBegin; 469*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 470*85317021SBarry Smith if (f) { 471*85317021SBarry Smith ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 472*85317021SBarry Smith } 473*85317021SBarry Smith PetscFunctionReturn(0); 474*85317021SBarry Smith } 475*85317021SBarry Smith 476*85317021SBarry Smith #undef __FUNCT__ 477*85317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks" 478*85317021SBarry Smith /*@ 479*85317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 480*85317021SBarry Smith with BAIJ or SBAIJ matrices 481*85317021SBarry Smith 482*85317021SBarry Smith Collective on PC 483*85317021SBarry Smith 484*85317021SBarry Smith Input Parameters: 485*85317021SBarry Smith + pc - the preconditioner context 486*85317021SBarry Smith - pivot - PETSC_TRUE or PETSC_FALSE 487*85317021SBarry Smith 488*85317021SBarry Smith Options Database Key: 489*85317021SBarry Smith . -pc_factor_pivot_in_blocks <true,false> 490*85317021SBarry Smith 491*85317021SBarry Smith Level: intermediate 492*85317021SBarry Smith 493*85317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 494*85317021SBarry Smith @*/ 495*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 496*85317021SBarry Smith { 497*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 498*85317021SBarry Smith 499*85317021SBarry Smith PetscFunctionBegin; 500*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 501*85317021SBarry Smith if (f) { 502*85317021SBarry Smith ierr = (*f)(pc,pivot);CHKERRQ(ierr); 503*85317021SBarry Smith } 504*85317021SBarry Smith PetscFunctionReturn(0); 505*85317021SBarry Smith } 506*85317021SBarry Smith 507*85317021SBarry Smith #undef __FUNCT__ 508*85317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill" 509*85317021SBarry Smith /*@ 510*85317021SBarry Smith PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 511*85317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 512*85317021SBarry Smith 513*85317021SBarry Smith Collective on PC 514*85317021SBarry Smith 515*85317021SBarry Smith Input Parameters: 516*85317021SBarry Smith + pc - the preconditioner context 517*85317021SBarry Smith - flag - PETSC_TRUE to reuse else PETSC_FALSE 518*85317021SBarry Smith 519*85317021SBarry Smith Options Database Key: 520*85317021SBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 521*85317021SBarry Smith 522*85317021SBarry Smith Level: intermediate 523*85317021SBarry Smith 524*85317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 525*85317021SBarry Smith 526*85317021SBarry Smith .seealso: PCFactorSetReuseOrdering() 527*85317021SBarry Smith @*/ 528*85317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag) 529*85317021SBarry Smith { 530*85317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 531*85317021SBarry Smith 532*85317021SBarry Smith PetscFunctionBegin; 533*85317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,2); 534*85317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 535*85317021SBarry Smith if (f) { 536*85317021SBarry Smith ierr = (*f)(pc,flag);CHKERRQ(ierr); 537*85317021SBarry Smith } 538*85317021SBarry Smith PetscFunctionReturn(0); 539*85317021SBarry Smith } 540