xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision d90ac83d909a2cb66ba66e74550740849e04bb53)
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__
39*d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType"
40*d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
41*d90ac83dSHong Zhang {
42*d90ac83dSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatFactorShiftType);
43*d90ac83dSHong Zhang 
44*d90ac83dSHong Zhang   PetscFunctionBegin;
45*d90ac83dSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
46*d90ac83dSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftType_C",(void (**)(void))&f);CHKERRQ(ierr);
47*d90ac83dSHong Zhang   if (f) {
48*d90ac83dSHong Zhang     ierr = (*f)(pc,shifttype);CHKERRQ(ierr);
49*d90ac83dSHong Zhang   }
50*d90ac83dSHong Zhang   PetscFunctionReturn(0);
51*d90ac83dSHong Zhang }
52*d90ac83dSHong Zhang 
53*d90ac83dSHong Zhang #undef __FUNCT__
54*d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount"
55*d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
56*d90ac83dSHong Zhang {
57*d90ac83dSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
58*d90ac83dSHong Zhang 
59*d90ac83dSHong Zhang   PetscFunctionBegin;
60*d90ac83dSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
61*d90ac83dSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",(void (**)(void))&f);CHKERRQ(ierr);
62*d90ac83dSHong Zhang   if (f) {
63*d90ac83dSHong Zhang     ierr = (*f)(pc,shiftamount);CHKERRQ(ierr);
64*d90ac83dSHong Zhang   }
65*d90ac83dSHong Zhang   PetscFunctionReturn(0);
66*d90ac83dSHong Zhang }
67*d90ac83dSHong Zhang 
68*d90ac83dSHong Zhang #undef __FUNCT__
695e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero"
705e8efad8SHong Zhang /*@
715e8efad8SHong Zhang    PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during
725e8efad8SHong Zhang      numerical factorization, thus the matrix has nonzero pivots
735e8efad8SHong Zhang 
745e8efad8SHong Zhang    Collective on PC
755e8efad8SHong Zhang 
765e8efad8SHong Zhang    Input Parameters:
77afaefe49SHong Zhang +  pc - the preconditioner context
78afaefe49SHong Zhang -  shift - amount of shift
795e8efad8SHong Zhang 
805e8efad8SHong Zhang    Options Database Key:
819f95998fSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
825e8efad8SHong Zhang 
835e8efad8SHong Zhang    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
845e8efad8SHong Zhang          pivot, then the shift is doubled until this is alleviated.
855e8efad8SHong Zhang 
865e8efad8SHong Zhang    Level: intermediate
875e8efad8SHong Zhang 
885e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill
895e8efad8SHong Zhang 
90d6e5152cSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks()
915e8efad8SHong Zhang @*/
92dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift)
935e8efad8SHong Zhang {
94afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
95afaefe49SHong Zhang 
965e8efad8SHong Zhang   PetscFunctionBegin;
97afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
98afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr);
99afaefe49SHong Zhang   if (f) {
100afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
1015e8efad8SHong Zhang   }
1025e8efad8SHong Zhang   PetscFunctionReturn(0);
1035e8efad8SHong Zhang }
1040a29876aSHong Zhang 
1050a29876aSHong Zhang #undef __FUNCT__
1060a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd"
1070a29876aSHong Zhang /*@
1080a29876aSHong Zhang    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
1090a29876aSHong Zhang    If a matrix factorisation breaks down because of nonpositive pivots,
1100a29876aSHong Zhang    adding sufficient identity to the diagonal will remedy this.
1110a29876aSHong Zhang    Setting this causes a bisection method to find the minimum shift that
1120a29876aSHong Zhang    will lead to a well-defined matrix factor.
1130a29876aSHong Zhang 
1140a29876aSHong Zhang    Collective on PC
1150a29876aSHong Zhang 
1160a29876aSHong Zhang    Input parameters:
117afaefe49SHong Zhang +  pc - the preconditioner context
118afaefe49SHong Zhang -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
1190a29876aSHong Zhang 
1200a29876aSHong Zhang    Options Database Key:
12153ae6812SBarry Smith .  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
12253ae6812SBarry Smith    is optional with false being the default
1230a29876aSHong Zhang 
1240a29876aSHong Zhang    Level: intermediate
1250a29876aSHong Zhang 
1260a29876aSHong Zhang .keywords: PC, indefinite, factorization
1270a29876aSHong Zhang 
128d6e5152cSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftInBlocks()
1290a29876aSHong Zhang @*/
130dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift)
1310a29876aSHong Zhang {
132afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
133afaefe49SHong Zhang 
1340a29876aSHong Zhang   PetscFunctionBegin;
135afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
136afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr);
137afaefe49SHong Zhang   if (f) {
138afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
139afaefe49SHong Zhang   }
1400a29876aSHong Zhang   PetscFunctionReturn(0);
1410a29876aSHong Zhang }
14285317021SBarry Smith 
14385317021SBarry Smith #undef __FUNCT__
144b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance"
14585317021SBarry Smith /*@
146b7c853c4SBarry Smith    PCFactorSetDropTolerance - The preconditioner will use an ILU
14785317021SBarry Smith    based on a drop tolerance.
14885317021SBarry Smith 
14985317021SBarry Smith    Collective on PC
15085317021SBarry Smith 
15185317021SBarry Smith    Input Parameters:
15285317021SBarry Smith +  pc - the preconditioner context
15385317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
15485317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
15585317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
15685317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
15785317021SBarry Smith 
15885317021SBarry Smith    Options Database Key:
159b7c853c4SBarry Smith .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
16085317021SBarry Smith 
16185317021SBarry Smith    Level: intermediate
16285317021SBarry Smith 
16385317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
16485317021SBarry Smith       matrix. We don't know how to compute reasonable values.
16585317021SBarry Smith 
16685317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU
16785317021SBarry Smith @*/
168b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
16985317021SBarry Smith {
17085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
17185317021SBarry Smith 
17285317021SBarry Smith   PetscFunctionBegin;
17385317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
174b7c853c4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
17585317021SBarry Smith   if (f) {
17685317021SBarry Smith     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
17785317021SBarry Smith   }
17885317021SBarry Smith   PetscFunctionReturn(0);
17985317021SBarry Smith }
18085317021SBarry Smith 
18185317021SBarry Smith #undef __FUNCT__
18285317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels"
18385317021SBarry Smith /*@
18485317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
18585317021SBarry Smith 
18685317021SBarry Smith    Collective on PC
18785317021SBarry Smith 
18885317021SBarry Smith    Input Parameters:
18985317021SBarry Smith +  pc - the preconditioner context
19085317021SBarry Smith -  levels - number of levels of fill
19185317021SBarry Smith 
19285317021SBarry Smith    Options Database Key:
19385317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
19485317021SBarry Smith 
19585317021SBarry Smith    Level: intermediate
19685317021SBarry Smith 
19785317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
19885317021SBarry Smith @*/
19985317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
20085317021SBarry Smith {
20185317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
20285317021SBarry Smith 
20385317021SBarry Smith   PetscFunctionBegin;
20485317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
20585317021SBarry Smith   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
20685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
20785317021SBarry Smith   if (f) {
20885317021SBarry Smith     ierr = (*f)(pc,levels);CHKERRQ(ierr);
20985317021SBarry Smith   }
21085317021SBarry Smith   PetscFunctionReturn(0);
21185317021SBarry Smith }
21285317021SBarry Smith 
21385317021SBarry Smith #undef __FUNCT__
21485317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
21585317021SBarry Smith /*@
21685317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
21785317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
21885317021SBarry Smith 
21985317021SBarry Smith    Collective on PC
22085317021SBarry Smith 
22185317021SBarry Smith    Input Parameters:
22285317021SBarry Smith +  pc - the preconditioner context
22385317021SBarry Smith 
22485317021SBarry Smith    Options Database Key:
22585317021SBarry Smith .  -pc_factor_diagonal_fill
22685317021SBarry Smith 
22785317021SBarry Smith    Notes:
22885317021SBarry Smith    Does not apply with 0 fill.
22985317021SBarry Smith 
23085317021SBarry Smith    Level: intermediate
23185317021SBarry Smith 
23285317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
23385317021SBarry Smith @*/
23485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
23585317021SBarry Smith {
23685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
23785317021SBarry Smith 
23885317021SBarry Smith   PetscFunctionBegin;
23985317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
24085317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
24185317021SBarry Smith   if (f) {
24285317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
24385317021SBarry Smith   }
24485317021SBarry Smith   PetscFunctionReturn(0);
24585317021SBarry Smith }
24685317021SBarry Smith 
24785317021SBarry Smith #undef __FUNCT__
24885317021SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks"
24985317021SBarry Smith /*@
250d6e5152cSHong Zhang     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if matrix factor on that block detects a zero pivot.
25185317021SBarry Smith 
25285317021SBarry Smith     Collective on PC
25385317021SBarry Smith 
25485317021SBarry Smith     Input Parameters:
25585317021SBarry Smith +   pc - the preconditioner context
25685317021SBarry Smith -   shift - amount to shift or PETSC_DECIDE
25785317021SBarry Smith 
25885317021SBarry Smith     Options Database Key:
25985317021SBarry Smith .   -pc_factor_shift_in_blocks <shift>
26085317021SBarry Smith 
26185317021SBarry Smith     Level: intermediate
26285317021SBarry Smith 
2638ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot(), PCFactorSetShiftInBlocks()
26485317021SBarry Smith @*/
26585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
26685317021SBarry Smith {
26785317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
26885317021SBarry Smith 
26985317021SBarry Smith   PetscFunctionBegin;
27085317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
27185317021SBarry Smith   if (f) {
27285317021SBarry Smith     ierr = (*f)(pc,shift);CHKERRQ(ierr);
27385317021SBarry Smith   }
27485317021SBarry Smith   PetscFunctionReturn(0);
27585317021SBarry Smith }
27685317021SBarry Smith 
27785317021SBarry Smith #undef __FUNCT__
27885317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
27985317021SBarry Smith /*@
28085317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
28185317021SBarry Smith 
28285317021SBarry Smith    Collective on PC
28385317021SBarry Smith 
28485317021SBarry Smith    Input Parameters:
28585317021SBarry Smith +  pc - the preconditioner context
28685317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
28785317021SBarry Smith 
28885317021SBarry Smith    Options Database Key:
28985317021SBarry Smith .  -pc_factor_nonzeros_along_diagonal
29085317021SBarry Smith 
29185317021SBarry Smith    Level: intermediate
29285317021SBarry Smith 
29385317021SBarry Smith .keywords: PC, set, factorization, direct, fill
29485317021SBarry Smith 
29585317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
29685317021SBarry Smith @*/
29785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
29885317021SBarry Smith {
29985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
30085317021SBarry Smith 
30185317021SBarry Smith   PetscFunctionBegin;
30285317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
30385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
30485317021SBarry Smith   if (f) {
30585317021SBarry Smith     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
30685317021SBarry Smith   }
30785317021SBarry Smith   PetscFunctionReturn(0);
30885317021SBarry Smith }
30985317021SBarry Smith 
31085317021SBarry Smith #undef __FUNCT__
31185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
312bf6011e8SBarry Smith /*@C
31385317021SBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
31485317021SBarry Smith 
31585317021SBarry Smith    Collective on PC
31685317021SBarry Smith 
31785317021SBarry Smith    Input Parameters:
31885317021SBarry Smith +  pc - the preconditioner context
31985317021SBarry Smith -  stype - for example, spooles, superlu, superlu_dist
32085317021SBarry Smith 
32185317021SBarry Smith    Options Database Key:
32285317021SBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
32385317021SBarry Smith 
32485317021SBarry Smith    Level: intermediate
32585317021SBarry Smith 
32685317021SBarry Smith    Note:
32785317021SBarry Smith      By default this will use the PETSc factorization if it exists
32885317021SBarry Smith 
32985317021SBarry Smith 
33085317021SBarry Smith .keywords: PC, set, factorization, direct, fill
33185317021SBarry Smith 
3327112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
33385317021SBarry Smith 
33485317021SBarry Smith @*/
33585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
33685317021SBarry Smith {
33785317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
33885317021SBarry Smith 
33985317021SBarry Smith   PetscFunctionBegin;
34085317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
34185317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
34285317021SBarry Smith   if (f) {
34385317021SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
34485317021SBarry Smith   }
34585317021SBarry Smith   PetscFunctionReturn(0);
34685317021SBarry Smith }
34785317021SBarry Smith 
34885317021SBarry Smith #undef __FUNCT__
3497112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage"
350bf6011e8SBarry Smith /*@C
3517112b564SBarry Smith    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
3527112b564SBarry Smith 
3537112b564SBarry Smith    Collective on PC
3547112b564SBarry Smith 
3557112b564SBarry Smith    Input Parameter:
3567112b564SBarry Smith .  pc - the preconditioner context
3577112b564SBarry Smith 
3587112b564SBarry Smith    Output Parameter:
3597112b564SBarry Smith .   stype - for example, spooles, superlu, superlu_dist
3607112b564SBarry Smith 
3617112b564SBarry Smith    Level: intermediate
3627112b564SBarry Smith 
3637112b564SBarry Smith 
3647112b564SBarry Smith .keywords: PC, set, factorization, direct, fill
3657112b564SBarry Smith 
3667112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
3677112b564SBarry Smith 
3687112b564SBarry Smith @*/
3697112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
3707112b564SBarry Smith {
3717112b564SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
3727112b564SBarry Smith 
3737112b564SBarry Smith   PetscFunctionBegin;
3747112b564SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3757112b564SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
3767112b564SBarry Smith   if (f) {
3777112b564SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
3787112b564SBarry Smith   }
3797112b564SBarry Smith   PetscFunctionReturn(0);
3807112b564SBarry Smith }
3817112b564SBarry Smith 
3827112b564SBarry Smith #undef __FUNCT__
38385317021SBarry Smith #define __FUNCT__ "PCFactorSetFill"
38485317021SBarry Smith /*@
38585317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
38685317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
38785317021SBarry Smith 
38885317021SBarry Smith    Collective on PC
38985317021SBarry Smith 
39085317021SBarry Smith    Input Parameters:
39185317021SBarry Smith +  pc - the preconditioner context
39285317021SBarry Smith -  fill - amount of expected fill
39385317021SBarry Smith 
39485317021SBarry Smith    Options Database Key:
39585317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
39685317021SBarry Smith 
39785317021SBarry Smith    Level: intermediate
39885317021SBarry Smith 
39985317021SBarry Smith    Note:
40085317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
40185317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
40285317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
40385317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
40485317021SBarry Smith 
40585317021SBarry Smith .keywords: PC, set, factorization, direct, fill
40685317021SBarry Smith 
40785317021SBarry Smith @*/
40885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
40985317021SBarry Smith {
41085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
41185317021SBarry Smith 
41285317021SBarry Smith   PetscFunctionBegin;
41385317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
41485317021SBarry Smith   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
41585317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
41685317021SBarry Smith   if (f) {
41785317021SBarry Smith     ierr = (*f)(pc,fill);CHKERRQ(ierr);
41885317021SBarry Smith   }
41985317021SBarry Smith   PetscFunctionReturn(0);
42085317021SBarry Smith }
42185317021SBarry Smith 
42285317021SBarry Smith #undef __FUNCT__
42385317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
42485317021SBarry Smith /*@
42585317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
42685317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
42785317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
42885317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
42985317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
43085317021SBarry Smith    space.
43185317021SBarry Smith 
43285317021SBarry Smith    Collective on PC
43385317021SBarry Smith 
43485317021SBarry Smith    Input Parameters:
43585317021SBarry Smith .  pc - the preconditioner context
43685317021SBarry Smith 
43785317021SBarry Smith    Options Database Key:
43885317021SBarry Smith .  -pc_factor_in_place - Activates in-place factorization
43985317021SBarry Smith 
44085317021SBarry Smith    Notes:
44185317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
44285317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
44385317021SBarry Smith    a call to KSPSetOperators().
44485317021SBarry Smith    This is because the Krylov space methods require an application of the
44585317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
44685317021SBarry Smith    been factored in-place, replacing the original matrix.
44785317021SBarry Smith 
44885317021SBarry Smith    Level: intermediate
44985317021SBarry Smith 
45085317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
45185317021SBarry Smith 
45285317021SBarry Smith .seealso: PCILUSetUseInPlace()
45385317021SBarry Smith @*/
45485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
45585317021SBarry Smith {
45685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
45785317021SBarry Smith 
45885317021SBarry Smith   PetscFunctionBegin;
45985317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
46085317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
46185317021SBarry Smith   if (f) {
46285317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
46385317021SBarry Smith   }
46485317021SBarry Smith   PetscFunctionReturn(0);
46585317021SBarry Smith }
46685317021SBarry Smith 
46785317021SBarry Smith #undef __FUNCT__
46885317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
46985317021SBarry Smith /*@C
47085317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
47185317021SBarry Smith     be used in the LU factorization.
47285317021SBarry Smith 
47385317021SBarry Smith     Collective on PC
47485317021SBarry Smith 
47585317021SBarry Smith     Input Parameters:
47685317021SBarry Smith +   pc - the preconditioner context
47785317021SBarry Smith -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
47885317021SBarry Smith 
47985317021SBarry Smith     Options Database Key:
48085317021SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
48185317021SBarry Smith 
48285317021SBarry Smith     Level: intermediate
48385317021SBarry Smith 
48485317021SBarry Smith     Notes: nested dissection is used by default
48585317021SBarry Smith 
48685317021SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
48785317021SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
48885317021SBarry Smith     SeqAIJ format in this case to get reorderings.
48985317021SBarry Smith 
49085317021SBarry Smith @*/
49185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
49285317021SBarry Smith {
49385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
49485317021SBarry Smith 
49585317021SBarry Smith   PetscFunctionBegin;
49685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
49785317021SBarry Smith   if (f) {
49885317021SBarry Smith     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
49985317021SBarry Smith   }
50085317021SBarry Smith   PetscFunctionReturn(0);
50185317021SBarry Smith }
50285317021SBarry Smith 
50385317021SBarry Smith #undef __FUNCT__
5048ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot"
50585317021SBarry Smith /*@
5068ff23777SHong Zhang     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
50785317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
50885317021SBarry Smith       it is never done. For the Matlab and SuperLU factorization this is used.
50985317021SBarry Smith 
51085317021SBarry Smith     Collective on PC
51185317021SBarry Smith 
51285317021SBarry Smith     Input Parameters:
51385317021SBarry Smith +   pc - the preconditioner context
51485317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
51585317021SBarry Smith 
51685317021SBarry Smith     Options Database Key:
51785317021SBarry Smith .   -pc_factor_pivoting <dtcol>
51885317021SBarry Smith 
51985317021SBarry Smith     Level: intermediate
52085317021SBarry Smith 
52185317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
52285317021SBarry Smith @*/
5238ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
52485317021SBarry Smith {
52585317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
52685317021SBarry Smith 
52785317021SBarry Smith   PetscFunctionBegin;
5288ff23777SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
52985317021SBarry Smith   if (f) {
53085317021SBarry Smith     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
53185317021SBarry Smith   }
53285317021SBarry Smith   PetscFunctionReturn(0);
53385317021SBarry Smith }
53485317021SBarry Smith 
53585317021SBarry Smith #undef __FUNCT__
53685317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
53785317021SBarry Smith /*@
53885317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
53985317021SBarry Smith       with BAIJ or SBAIJ matrices
54085317021SBarry Smith 
54185317021SBarry Smith     Collective on PC
54285317021SBarry Smith 
54385317021SBarry Smith     Input Parameters:
54485317021SBarry Smith +   pc - the preconditioner context
54585317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
54685317021SBarry Smith 
54785317021SBarry Smith     Options Database Key:
54885317021SBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
54985317021SBarry Smith 
55085317021SBarry Smith     Level: intermediate
55185317021SBarry Smith 
5528ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
55385317021SBarry Smith @*/
55485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
55585317021SBarry Smith {
55685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
55785317021SBarry Smith 
55885317021SBarry Smith   PetscFunctionBegin;
55985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
56085317021SBarry Smith   if (f) {
56185317021SBarry Smith     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
56285317021SBarry Smith   }
56385317021SBarry Smith   PetscFunctionReturn(0);
56485317021SBarry Smith }
56585317021SBarry Smith 
56685317021SBarry Smith #undef __FUNCT__
56785317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill"
56885317021SBarry Smith /*@
56985317021SBarry Smith    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
57085317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
57185317021SBarry Smith 
57285317021SBarry Smith    Collective on PC
57385317021SBarry Smith 
57485317021SBarry Smith    Input Parameters:
57585317021SBarry Smith +  pc - the preconditioner context
57685317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
57785317021SBarry Smith 
57885317021SBarry Smith    Options Database Key:
57985317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
58085317021SBarry Smith 
58185317021SBarry Smith    Level: intermediate
58285317021SBarry Smith 
58385317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
58485317021SBarry Smith 
58585317021SBarry Smith .seealso: PCFactorSetReuseOrdering()
58685317021SBarry Smith @*/
58785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
58885317021SBarry Smith {
58985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
59085317021SBarry Smith 
59185317021SBarry Smith   PetscFunctionBegin;
59285317021SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
59385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
59485317021SBarry Smith   if (f) {
59585317021SBarry Smith     ierr = (*f)(pc,flag);CHKERRQ(ierr);
59685317021SBarry Smith   }
59785317021SBarry Smith   PetscFunctionReturn(0);
59885317021SBarry Smith }
599