1dba47a55SKris Buschelman #define PETSCKSP_DLL 25e8efad8SHong Zhang 37c4f633dSBarry 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 23d6e5152cSHong Zhang .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks() 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__ 39d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType" 40*915743fcSHong Zhang /*@ 41*915743fcSHong Zhang PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 42*915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 43*915743fcSHong Zhang 44*915743fcSHong Zhang Collective on PC 45*915743fcSHong Zhang 46*915743fcSHong Zhang Input Parameters: 47*915743fcSHong Zhang + pc - the preconditioner context 48*915743fcSHong Zhang - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS 49*915743fcSHong Zhang 50*915743fcSHong Zhang Options Database Key: 51*915743fcSHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 52*915743fcSHong Zhang 53*915743fcSHong Zhang Level: intermediate 54*915743fcSHong Zhang 55*915743fcSHong Zhang .keywords: PC, set, factorization, 56*915743fcSHong Zhang 57*915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount() 58*915743fcSHong Zhang @*/ 59d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype) 60d90ac83dSHong Zhang { 61d90ac83dSHong Zhang PetscErrorCode ierr,(*f)(PC,MatFactorShiftType); 62d90ac83dSHong Zhang 63d90ac83dSHong Zhang PetscFunctionBegin; 64d90ac83dSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 65d90ac83dSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftType_C",(void (**)(void))&f);CHKERRQ(ierr); 66d90ac83dSHong Zhang if (f) { 67d90ac83dSHong Zhang ierr = (*f)(pc,shifttype);CHKERRQ(ierr); 68d90ac83dSHong Zhang } 69d90ac83dSHong Zhang PetscFunctionReturn(0); 70d90ac83dSHong Zhang } 71d90ac83dSHong Zhang 72d90ac83dSHong Zhang #undef __FUNCT__ 73d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount" 74*915743fcSHong Zhang /*@ 75*915743fcSHong Zhang PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 76*915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 77*915743fcSHong Zhang 78*915743fcSHong Zhang Collective on PC 79*915743fcSHong Zhang 80*915743fcSHong Zhang Input Parameters: 81*915743fcSHong Zhang + pc - the preconditioner context 82*915743fcSHong Zhang - shiftamount - amount of shift 83*915743fcSHong Zhang 84*915743fcSHong Zhang Options Database Key: 85*915743fcSHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 86*915743fcSHong Zhang 87*915743fcSHong Zhang Level: intermediate 88*915743fcSHong Zhang 89*915743fcSHong Zhang .keywords: PC, set, factorization, 90*915743fcSHong Zhang 91*915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType() 92*915743fcSHong Zhang @*/ 93d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC pc,PetscReal shiftamount) 94d90ac83dSHong Zhang { 95d90ac83dSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 96d90ac83dSHong Zhang 97d90ac83dSHong Zhang PetscFunctionBegin; 98d90ac83dSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 99d90ac83dSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",(void (**)(void))&f);CHKERRQ(ierr); 100d90ac83dSHong Zhang if (f) { 101d90ac83dSHong Zhang ierr = (*f)(pc,shiftamount);CHKERRQ(ierr); 102d90ac83dSHong Zhang } 103d90ac83dSHong Zhang PetscFunctionReturn(0); 104d90ac83dSHong Zhang } 105d90ac83dSHong Zhang 106d90ac83dSHong Zhang #undef __FUNCT__ 1075e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero" 1085e8efad8SHong Zhang /*@ 1095e8efad8SHong Zhang PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during 1105e8efad8SHong Zhang numerical factorization, thus the matrix has nonzero pivots 1115e8efad8SHong Zhang 1125e8efad8SHong Zhang Collective on PC 1135e8efad8SHong Zhang 1145e8efad8SHong Zhang Input Parameters: 115afaefe49SHong Zhang + pc - the preconditioner context 116afaefe49SHong Zhang - shift - amount of shift 1175e8efad8SHong Zhang 1185e8efad8SHong Zhang Options Database Key: 1199f95998fSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 1205e8efad8SHong Zhang 1215e8efad8SHong Zhang Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero 1225e8efad8SHong Zhang pivot, then the shift is doubled until this is alleviated. 1235e8efad8SHong Zhang 1245e8efad8SHong Zhang Level: intermediate 1255e8efad8SHong Zhang 1265e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill 1275e8efad8SHong Zhang 128d6e5152cSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks() 1295e8efad8SHong Zhang @*/ 130dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift) 1315e8efad8SHong Zhang { 132afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 133afaefe49SHong Zhang 1345e8efad8SHong Zhang PetscFunctionBegin; 135afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 136afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr); 137afaefe49SHong Zhang if (f) { 138afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 1395e8efad8SHong Zhang } 1405e8efad8SHong Zhang PetscFunctionReturn(0); 1415e8efad8SHong Zhang } 1420a29876aSHong Zhang 1430a29876aSHong Zhang #undef __FUNCT__ 1440a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd" 1450a29876aSHong Zhang /*@ 1460a29876aSHong Zhang PCFactorSetShiftPd - specify whether to use Manteuffel shifting. 1470a29876aSHong Zhang If a matrix factorisation breaks down because of nonpositive pivots, 1480a29876aSHong Zhang adding sufficient identity to the diagonal will remedy this. 1490a29876aSHong Zhang Setting this causes a bisection method to find the minimum shift that 1500a29876aSHong Zhang will lead to a well-defined matrix factor. 1510a29876aSHong Zhang 1520a29876aSHong Zhang Collective on PC 1530a29876aSHong Zhang 1540a29876aSHong Zhang Input parameters: 155afaefe49SHong Zhang + pc - the preconditioner context 156afaefe49SHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 1570a29876aSHong Zhang 1580a29876aSHong Zhang Options Database Key: 15953ae6812SBarry Smith . -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value 16053ae6812SBarry Smith is optional with false being the default 1610a29876aSHong Zhang 1620a29876aSHong Zhang Level: intermediate 1630a29876aSHong Zhang 1640a29876aSHong Zhang .keywords: PC, indefinite, factorization 1650a29876aSHong Zhang 166d6e5152cSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftInBlocks() 1670a29876aSHong Zhang @*/ 168dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift) 1690a29876aSHong Zhang { 170afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 171afaefe49SHong Zhang 1720a29876aSHong Zhang PetscFunctionBegin; 173afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 174afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr); 175afaefe49SHong Zhang if (f) { 176afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 177afaefe49SHong Zhang } 1780a29876aSHong Zhang PetscFunctionReturn(0); 1790a29876aSHong Zhang } 18085317021SBarry Smith 18185317021SBarry Smith #undef __FUNCT__ 182b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance" 18385317021SBarry Smith /*@ 184b7c853c4SBarry Smith PCFactorSetDropTolerance - The preconditioner will use an ILU 18585317021SBarry Smith based on a drop tolerance. 18685317021SBarry Smith 18785317021SBarry Smith Collective on PC 18885317021SBarry Smith 18985317021SBarry Smith Input Parameters: 19085317021SBarry Smith + pc - the preconditioner context 19185317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 19285317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 19385317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 19485317021SBarry Smith depends on the number of nonzeros in row of original matrix 19585317021SBarry Smith 19685317021SBarry Smith Options Database Key: 197b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 19885317021SBarry Smith 19985317021SBarry Smith Level: intermediate 20085317021SBarry Smith 20185317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 20285317021SBarry Smith matrix. We don't know how to compute reasonable values. 20385317021SBarry Smith 20485317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU 20585317021SBarry Smith @*/ 206b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 20785317021SBarry Smith { 20885317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 20985317021SBarry Smith 21085317021SBarry Smith PetscFunctionBegin; 21185317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 212b7c853c4SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 21385317021SBarry Smith if (f) { 21485317021SBarry Smith ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 21585317021SBarry Smith } 21685317021SBarry Smith PetscFunctionReturn(0); 21785317021SBarry Smith } 21885317021SBarry Smith 21985317021SBarry Smith #undef __FUNCT__ 22085317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels" 22185317021SBarry Smith /*@ 22285317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 22385317021SBarry Smith 22485317021SBarry Smith Collective on PC 22585317021SBarry Smith 22685317021SBarry Smith Input Parameters: 22785317021SBarry Smith + pc - the preconditioner context 22885317021SBarry Smith - levels - number of levels of fill 22985317021SBarry Smith 23085317021SBarry Smith Options Database Key: 23185317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 23285317021SBarry Smith 23385317021SBarry Smith Level: intermediate 23485317021SBarry Smith 23585317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU 23685317021SBarry Smith @*/ 23785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 23885317021SBarry Smith { 23985317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 24085317021SBarry Smith 24185317021SBarry Smith PetscFunctionBegin; 24285317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 24385317021SBarry Smith if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 24485317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 24585317021SBarry Smith if (f) { 24685317021SBarry Smith ierr = (*f)(pc,levels);CHKERRQ(ierr); 24785317021SBarry Smith } 24885317021SBarry Smith PetscFunctionReturn(0); 24985317021SBarry Smith } 25085317021SBarry Smith 25185317021SBarry Smith #undef __FUNCT__ 25285317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 25385317021SBarry Smith /*@ 25485317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 25585317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 25685317021SBarry Smith 25785317021SBarry Smith Collective on PC 25885317021SBarry Smith 25985317021SBarry Smith Input Parameters: 26085317021SBarry Smith + pc - the preconditioner context 26185317021SBarry Smith 26285317021SBarry Smith Options Database Key: 26385317021SBarry Smith . -pc_factor_diagonal_fill 26485317021SBarry Smith 26585317021SBarry Smith Notes: 26685317021SBarry Smith Does not apply with 0 fill. 26785317021SBarry Smith 26885317021SBarry Smith Level: intermediate 26985317021SBarry Smith 27085317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU 27185317021SBarry Smith @*/ 27285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 27385317021SBarry Smith { 27485317021SBarry Smith PetscErrorCode ierr,(*f)(PC); 27585317021SBarry Smith 27685317021SBarry Smith PetscFunctionBegin; 27785317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 27885317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 27985317021SBarry Smith if (f) { 28085317021SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 28185317021SBarry Smith } 28285317021SBarry Smith PetscFunctionReturn(0); 28385317021SBarry Smith } 28485317021SBarry Smith 28585317021SBarry Smith #undef __FUNCT__ 28685317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks" 28785317021SBarry Smith /*@ 288d6e5152cSHong Zhang PCFactorSetShiftInBlocks - Shifts the diagonal of any block if matrix factor on that block detects a zero pivot. 28985317021SBarry Smith 29085317021SBarry Smith Collective on PC 29185317021SBarry Smith 29285317021SBarry Smith Input Parameters: 29385317021SBarry Smith + pc - the preconditioner context 29485317021SBarry Smith - shift - amount to shift or PETSC_DECIDE 29585317021SBarry Smith 29685317021SBarry Smith Options Database Key: 29785317021SBarry Smith . -pc_factor_shift_in_blocks <shift> 29885317021SBarry Smith 29985317021SBarry Smith Level: intermediate 30085317021SBarry Smith 3018ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot(), PCFactorSetShiftInBlocks() 30285317021SBarry Smith @*/ 30385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift) 30485317021SBarry Smith { 30585317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 30685317021SBarry Smith 30785317021SBarry Smith PetscFunctionBegin; 30885317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 30985317021SBarry Smith if (f) { 31085317021SBarry Smith ierr = (*f)(pc,shift);CHKERRQ(ierr); 31185317021SBarry Smith } 31285317021SBarry Smith PetscFunctionReturn(0); 31385317021SBarry Smith } 31485317021SBarry Smith 31585317021SBarry Smith #undef __FUNCT__ 31685317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 31785317021SBarry Smith /*@ 31885317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 31985317021SBarry Smith 32085317021SBarry Smith Collective on PC 32185317021SBarry Smith 32285317021SBarry Smith Input Parameters: 32385317021SBarry Smith + pc - the preconditioner context 32485317021SBarry Smith - tol - diagonal entries smaller than this in absolute value are considered zero 32585317021SBarry Smith 32685317021SBarry Smith Options Database Key: 32785317021SBarry Smith . -pc_factor_nonzeros_along_diagonal 32885317021SBarry Smith 32985317021SBarry Smith Level: intermediate 33085317021SBarry Smith 33185317021SBarry Smith .keywords: PC, set, factorization, direct, fill 33285317021SBarry Smith 33385317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 33485317021SBarry Smith @*/ 33585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 33685317021SBarry Smith { 33785317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 33885317021SBarry Smith 33985317021SBarry Smith PetscFunctionBegin; 34085317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 34185317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 34285317021SBarry Smith if (f) { 34385317021SBarry Smith ierr = (*f)(pc,rtol);CHKERRQ(ierr); 34485317021SBarry Smith } 34585317021SBarry Smith PetscFunctionReturn(0); 34685317021SBarry Smith } 34785317021SBarry Smith 34885317021SBarry Smith #undef __FUNCT__ 34985317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage" 350bf6011e8SBarry Smith /*@C 35185317021SBarry Smith PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 35285317021SBarry Smith 35385317021SBarry Smith Collective on PC 35485317021SBarry Smith 35585317021SBarry Smith Input Parameters: 35685317021SBarry Smith + pc - the preconditioner context 35785317021SBarry Smith - stype - for example, spooles, superlu, superlu_dist 35885317021SBarry Smith 35985317021SBarry Smith Options Database Key: 36085317021SBarry Smith . -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps 36185317021SBarry Smith 36285317021SBarry Smith Level: intermediate 36385317021SBarry Smith 36485317021SBarry Smith Note: 36585317021SBarry Smith By default this will use the PETSc factorization if it exists 36685317021SBarry Smith 36785317021SBarry Smith 36885317021SBarry Smith .keywords: PC, set, factorization, direct, fill 36985317021SBarry Smith 3707112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 37185317021SBarry Smith 37285317021SBarry Smith @*/ 37385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 37485317021SBarry Smith { 37585317021SBarry Smith PetscErrorCode ierr,(*f)(PC,const MatSolverPackage); 37685317021SBarry Smith 37785317021SBarry Smith PetscFunctionBegin; 37885317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 37985317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr); 38085317021SBarry Smith if (f) { 38185317021SBarry Smith ierr = (*f)(pc,stype);CHKERRQ(ierr); 38285317021SBarry Smith } 38385317021SBarry Smith PetscFunctionReturn(0); 38485317021SBarry Smith } 38585317021SBarry Smith 38685317021SBarry Smith #undef __FUNCT__ 3877112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage" 388bf6011e8SBarry Smith /*@C 3897112b564SBarry Smith PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization 3907112b564SBarry Smith 3917112b564SBarry Smith Collective on PC 3927112b564SBarry Smith 3937112b564SBarry Smith Input Parameter: 3947112b564SBarry Smith . pc - the preconditioner context 3957112b564SBarry Smith 3967112b564SBarry Smith Output Parameter: 3977112b564SBarry Smith . stype - for example, spooles, superlu, superlu_dist 3987112b564SBarry Smith 3997112b564SBarry Smith Level: intermediate 4007112b564SBarry Smith 4017112b564SBarry Smith 4027112b564SBarry Smith .keywords: PC, set, factorization, direct, fill 4037112b564SBarry Smith 4047112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage() 4057112b564SBarry Smith 4067112b564SBarry Smith @*/ 4077112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype) 4087112b564SBarry Smith { 4097112b564SBarry Smith PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*); 4107112b564SBarry Smith 4117112b564SBarry Smith PetscFunctionBegin; 4127112b564SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4137112b564SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr); 4147112b564SBarry Smith if (f) { 4157112b564SBarry Smith ierr = (*f)(pc,stype);CHKERRQ(ierr); 4167112b564SBarry Smith } 4177112b564SBarry Smith PetscFunctionReturn(0); 4187112b564SBarry Smith } 4197112b564SBarry Smith 4207112b564SBarry Smith #undef __FUNCT__ 42185317021SBarry Smith #define __FUNCT__ "PCFactorSetFill" 42285317021SBarry Smith /*@ 42385317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 42485317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 42585317021SBarry Smith 42685317021SBarry Smith Collective on PC 42785317021SBarry Smith 42885317021SBarry Smith Input Parameters: 42985317021SBarry Smith + pc - the preconditioner context 43085317021SBarry Smith - fill - amount of expected fill 43185317021SBarry Smith 43285317021SBarry Smith Options Database Key: 43385317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 43485317021SBarry Smith 43585317021SBarry Smith Level: intermediate 43685317021SBarry Smith 43785317021SBarry Smith Note: 43885317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 43985317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 44085317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 44185317021SBarry Smith future runs. Default PETSc uses a value of 5.0 44285317021SBarry Smith 44385317021SBarry Smith .keywords: PC, set, factorization, direct, fill 44485317021SBarry Smith 44585317021SBarry Smith @*/ 44685317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 44785317021SBarry Smith { 44885317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 44985317021SBarry Smith 45085317021SBarry Smith PetscFunctionBegin; 45185317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 45285317021SBarry Smith if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 45385317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 45485317021SBarry Smith if (f) { 45585317021SBarry Smith ierr = (*f)(pc,fill);CHKERRQ(ierr); 45685317021SBarry Smith } 45785317021SBarry Smith PetscFunctionReturn(0); 45885317021SBarry Smith } 45985317021SBarry Smith 46085317021SBarry Smith #undef __FUNCT__ 46185317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace" 46285317021SBarry Smith /*@ 46385317021SBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 46485317021SBarry Smith For dense matrices, this enables the solution of much larger problems. 46585317021SBarry Smith For sparse matrices the factorization cannot be done truly in-place 46685317021SBarry Smith so this does not save memory during the factorization, but after the matrix 46785317021SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 46885317021SBarry Smith space. 46985317021SBarry Smith 47085317021SBarry Smith Collective on PC 47185317021SBarry Smith 47285317021SBarry Smith Input Parameters: 47385317021SBarry Smith . pc - the preconditioner context 47485317021SBarry Smith 47585317021SBarry Smith Options Database Key: 47685317021SBarry Smith . -pc_factor_in_place - Activates in-place factorization 47785317021SBarry Smith 47885317021SBarry Smith Notes: 47985317021SBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 48085317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 48185317021SBarry Smith a call to KSPSetOperators(). 48285317021SBarry Smith This is because the Krylov space methods require an application of the 48385317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 48485317021SBarry Smith been factored in-place, replacing the original matrix. 48585317021SBarry Smith 48685317021SBarry Smith Level: intermediate 48785317021SBarry Smith 48885317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU 48985317021SBarry Smith 49085317021SBarry Smith .seealso: PCILUSetUseInPlace() 49185317021SBarry Smith @*/ 49285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 49385317021SBarry Smith { 49485317021SBarry Smith PetscErrorCode ierr,(*f)(PC); 49585317021SBarry Smith 49685317021SBarry Smith PetscFunctionBegin; 49785317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 49885317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 49985317021SBarry Smith if (f) { 50085317021SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 50185317021SBarry Smith } 50285317021SBarry Smith PetscFunctionReturn(0); 50385317021SBarry Smith } 50485317021SBarry Smith 50585317021SBarry Smith #undef __FUNCT__ 50685317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType" 50785317021SBarry Smith /*@C 50885317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 50985317021SBarry Smith be used in the LU factorization. 51085317021SBarry Smith 51185317021SBarry Smith Collective on PC 51285317021SBarry Smith 51385317021SBarry Smith Input Parameters: 51485317021SBarry Smith + pc - the preconditioner context 51585317021SBarry Smith - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 51685317021SBarry Smith 51785317021SBarry Smith Options Database Key: 51885317021SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 51985317021SBarry Smith 52085317021SBarry Smith Level: intermediate 52185317021SBarry Smith 52285317021SBarry Smith Notes: nested dissection is used by default 52385317021SBarry Smith 52485317021SBarry Smith For Cholesky and ICC and the SBAIJ format reorderings are not available, 52585317021SBarry Smith since only the upper triangular part of the matrix is stored. You can use the 52685317021SBarry Smith SeqAIJ format in this case to get reorderings. 52785317021SBarry Smith 52885317021SBarry Smith @*/ 52985317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering) 53085317021SBarry Smith { 53185317021SBarry Smith PetscErrorCode ierr,(*f)(PC,const MatOrderingType); 53285317021SBarry Smith 53385317021SBarry Smith PetscFunctionBegin; 53485317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 53585317021SBarry Smith if (f) { 53685317021SBarry Smith ierr = (*f)(pc,ordering);CHKERRQ(ierr); 53785317021SBarry Smith } 53885317021SBarry Smith PetscFunctionReturn(0); 53985317021SBarry Smith } 54085317021SBarry Smith 54185317021SBarry Smith #undef __FUNCT__ 5428ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot" 54385317021SBarry Smith /*@ 5448ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 54585317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 54685317021SBarry Smith it is never done. For the Matlab and SuperLU factorization this is used. 54785317021SBarry Smith 54885317021SBarry Smith Collective on PC 54985317021SBarry Smith 55085317021SBarry Smith Input Parameters: 55185317021SBarry Smith + pc - the preconditioner context 55285317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 55385317021SBarry Smith 55485317021SBarry Smith Options Database Key: 55585317021SBarry Smith . -pc_factor_pivoting <dtcol> 55685317021SBarry Smith 55785317021SBarry Smith Level: intermediate 55885317021SBarry Smith 55985317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 56085317021SBarry Smith @*/ 5618ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC pc,PetscReal dtcol) 56285317021SBarry Smith { 56385317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 56485317021SBarry Smith 56585317021SBarry Smith PetscFunctionBegin; 5668ff23777SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 56785317021SBarry Smith if (f) { 56885317021SBarry Smith ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 56985317021SBarry Smith } 57085317021SBarry Smith PetscFunctionReturn(0); 57185317021SBarry Smith } 57285317021SBarry Smith 57385317021SBarry Smith #undef __FUNCT__ 57485317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks" 57585317021SBarry Smith /*@ 57685317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 57785317021SBarry Smith with BAIJ or SBAIJ matrices 57885317021SBarry Smith 57985317021SBarry Smith Collective on PC 58085317021SBarry Smith 58185317021SBarry Smith Input Parameters: 58285317021SBarry Smith + pc - the preconditioner context 58385317021SBarry Smith - pivot - PETSC_TRUE or PETSC_FALSE 58485317021SBarry Smith 58585317021SBarry Smith Options Database Key: 58685317021SBarry Smith . -pc_factor_pivot_in_blocks <true,false> 58785317021SBarry Smith 58885317021SBarry Smith Level: intermediate 58985317021SBarry Smith 5908ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot() 59185317021SBarry Smith @*/ 59285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 59385317021SBarry Smith { 59485317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 59585317021SBarry Smith 59685317021SBarry Smith PetscFunctionBegin; 59785317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 59885317021SBarry Smith if (f) { 59985317021SBarry Smith ierr = (*f)(pc,pivot);CHKERRQ(ierr); 60085317021SBarry Smith } 60185317021SBarry Smith PetscFunctionReturn(0); 60285317021SBarry Smith } 60385317021SBarry Smith 60485317021SBarry Smith #undef __FUNCT__ 60585317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill" 60685317021SBarry Smith /*@ 60785317021SBarry Smith PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 60885317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 60985317021SBarry Smith 61085317021SBarry Smith Collective on PC 61185317021SBarry Smith 61285317021SBarry Smith Input Parameters: 61385317021SBarry Smith + pc - the preconditioner context 61485317021SBarry Smith - flag - PETSC_TRUE to reuse else PETSC_FALSE 61585317021SBarry Smith 61685317021SBarry Smith Options Database Key: 61785317021SBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 61885317021SBarry Smith 61985317021SBarry Smith Level: intermediate 62085317021SBarry Smith 62185317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 62285317021SBarry Smith 62385317021SBarry Smith .seealso: PCFactorSetReuseOrdering() 62485317021SBarry Smith @*/ 62585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag) 62685317021SBarry Smith { 62785317021SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 62885317021SBarry Smith 62985317021SBarry Smith PetscFunctionBegin; 63085317021SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,2); 63185317021SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 63285317021SBarry Smith if (f) { 63385317021SBarry Smith ierr = (*f)(pc,flag);CHKERRQ(ierr); 63485317021SBarry Smith } 63585317021SBarry Smith PetscFunctionReturn(0); 63685317021SBarry Smith } 637