xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision c5eb91543d2ee8daf880d93389b892228ddada03)
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 
10ad4df100SBarry Smith    Logically 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 
23daa17b54SHong Zhang .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
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;
300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
31*c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,zero,2);
32afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
33afaefe49SHong Zhang   if (f) {
34afaefe49SHong Zhang     ierr = (*f)(pc,zero);CHKERRQ(ierr);
35afaefe49SHong Zhang   }
36ee45ca4aSHong Zhang   PetscFunctionReturn(0);
37ee45ca4aSHong Zhang }
38ee45ca4aSHong Zhang 
395e8efad8SHong Zhang #undef __FUNCT__
40d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType"
41915743fcSHong Zhang /*@
42915743fcSHong Zhang    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
43915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
44915743fcSHong Zhang 
45ad4df100SBarry Smith    Logically Collective on PC
46915743fcSHong Zhang 
47915743fcSHong Zhang    Input Parameters:
48915743fcSHong Zhang +  pc - the preconditioner context
49915743fcSHong Zhang -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
50915743fcSHong Zhang 
51915743fcSHong Zhang    Options Database Key:
52915743fcSHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
53915743fcSHong Zhang 
54915743fcSHong Zhang    Level: intermediate
55915743fcSHong Zhang 
56915743fcSHong Zhang .keywords: PC, set, factorization,
57915743fcSHong Zhang 
58915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
59915743fcSHong Zhang @*/
60d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
61d90ac83dSHong Zhang {
62d90ac83dSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatFactorShiftType);
63d90ac83dSHong Zhang 
64d90ac83dSHong Zhang   PetscFunctionBegin;
650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
66*c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,shifttype,2);
67d90ac83dSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftType_C",(void (**)(void))&f);CHKERRQ(ierr);
68d90ac83dSHong Zhang   if (f) {
69d90ac83dSHong Zhang     ierr = (*f)(pc,shifttype);CHKERRQ(ierr);
70d90ac83dSHong Zhang   }
71d90ac83dSHong Zhang   PetscFunctionReturn(0);
72d90ac83dSHong Zhang }
73d90ac83dSHong Zhang 
74d90ac83dSHong Zhang #undef __FUNCT__
75d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount"
76915743fcSHong Zhang /*@
77915743fcSHong Zhang    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
78915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
79915743fcSHong Zhang 
80ad4df100SBarry Smith    Logically Collective on PC
81915743fcSHong Zhang 
82915743fcSHong Zhang    Input Parameters:
83915743fcSHong Zhang +  pc - the preconditioner context
84915743fcSHong Zhang -  shiftamount - amount of shift
85915743fcSHong Zhang 
86915743fcSHong Zhang    Options Database Key:
87915743fcSHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
88915743fcSHong Zhang 
89915743fcSHong Zhang    Level: intermediate
90915743fcSHong Zhang 
91915743fcSHong Zhang .keywords: PC, set, factorization,
92915743fcSHong Zhang 
93915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
94915743fcSHong Zhang @*/
95d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
96d90ac83dSHong Zhang {
97d90ac83dSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
98d90ac83dSHong Zhang 
99d90ac83dSHong Zhang   PetscFunctionBegin;
1000700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
101*c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,shiftamount,2);
102d90ac83dSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",(void (**)(void))&f);CHKERRQ(ierr);
103d90ac83dSHong Zhang   if (f) {
104d90ac83dSHong Zhang     ierr = (*f)(pc,shiftamount);CHKERRQ(ierr);
105d90ac83dSHong Zhang   }
106d90ac83dSHong Zhang   PetscFunctionReturn(0);
107d90ac83dSHong Zhang }
108d90ac83dSHong Zhang 
109d90ac83dSHong Zhang #undef __FUNCT__
110b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance"
11178fc6b22SHong Zhang /*
112b7c853c4SBarry Smith    PCFactorSetDropTolerance - The preconditioner will use an ILU
11378fc6b22SHong Zhang    based on a drop tolerance. (Under development)
11485317021SBarry Smith 
115ad4df100SBarry Smith    Logically Collective on PC
11685317021SBarry Smith 
11785317021SBarry Smith    Input Parameters:
11885317021SBarry Smith +  pc - the preconditioner context
11985317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
12085317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
12185317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
12285317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
12385317021SBarry Smith 
12485317021SBarry Smith    Options Database Key:
125b7c853c4SBarry Smith .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
12685317021SBarry Smith 
12785317021SBarry Smith    Level: intermediate
12885317021SBarry Smith 
12985317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
13085317021SBarry Smith       matrix. We don't know how to compute reasonable values.
13185317021SBarry Smith 
13285317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU
13378fc6b22SHong Zhang */
134b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
13585317021SBarry Smith {
13685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
13785317021SBarry Smith 
13885317021SBarry Smith   PetscFunctionBegin;
1390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
140*c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,dtcol,2);
141*c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,maxrowcount,3);
142b7c853c4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
14385317021SBarry Smith   if (f) {
14485317021SBarry Smith     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
14585317021SBarry Smith   }
14685317021SBarry Smith   PetscFunctionReturn(0);
14785317021SBarry Smith }
14885317021SBarry Smith 
14985317021SBarry Smith #undef __FUNCT__
15085317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels"
15185317021SBarry Smith /*@
15285317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
15385317021SBarry Smith 
154ad4df100SBarry Smith    Logically Collective on PC
15585317021SBarry Smith 
15685317021SBarry Smith    Input Parameters:
15785317021SBarry Smith +  pc - the preconditioner context
15885317021SBarry Smith -  levels - number of levels of fill
15985317021SBarry Smith 
16085317021SBarry Smith    Options Database Key:
16185317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
16285317021SBarry Smith 
16385317021SBarry Smith    Level: intermediate
16485317021SBarry Smith 
16585317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
16685317021SBarry Smith @*/
16785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
16885317021SBarry Smith {
16985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
17085317021SBarry Smith 
17185317021SBarry Smith   PetscFunctionBegin;
1720700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
173e7e72b3dSBarry Smith   if (levels < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
174*c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,levels,2);
17585317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
17685317021SBarry Smith   if (f) {
17785317021SBarry Smith     ierr = (*f)(pc,levels);CHKERRQ(ierr);
17885317021SBarry Smith   }
17985317021SBarry Smith   PetscFunctionReturn(0);
18085317021SBarry Smith }
18185317021SBarry Smith 
18285317021SBarry Smith #undef __FUNCT__
18385317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
18485317021SBarry Smith /*@
18585317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
18685317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
18785317021SBarry Smith 
188ad4df100SBarry Smith    Logically Collective on PC
18985317021SBarry Smith 
19085317021SBarry Smith    Input Parameters:
19185317021SBarry Smith +  pc - the preconditioner context
19285317021SBarry Smith 
19385317021SBarry Smith    Options Database Key:
19485317021SBarry Smith .  -pc_factor_diagonal_fill
19585317021SBarry Smith 
19685317021SBarry Smith    Notes:
19785317021SBarry Smith    Does not apply with 0 fill.
19885317021SBarry Smith 
19985317021SBarry Smith    Level: intermediate
20085317021SBarry Smith 
20185317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
20285317021SBarry Smith @*/
20385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
20485317021SBarry Smith {
20585317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
20685317021SBarry Smith 
20785317021SBarry Smith   PetscFunctionBegin;
2080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
20985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
21085317021SBarry Smith   if (f) {
21185317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
21285317021SBarry Smith   }
21385317021SBarry Smith   PetscFunctionReturn(0);
21485317021SBarry Smith }
21585317021SBarry Smith 
21685317021SBarry Smith #undef __FUNCT__
21785317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
21885317021SBarry Smith /*@
21985317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
22085317021SBarry Smith 
221ad4df100SBarry Smith    Logically Collective on PC
22285317021SBarry Smith 
22385317021SBarry Smith    Input Parameters:
22485317021SBarry Smith +  pc - the preconditioner context
22585317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
22685317021SBarry Smith 
22785317021SBarry Smith    Options Database Key:
22885317021SBarry Smith .  -pc_factor_nonzeros_along_diagonal
22985317021SBarry Smith 
23085317021SBarry Smith    Level: intermediate
23185317021SBarry Smith 
23285317021SBarry Smith .keywords: PC, set, factorization, direct, fill
23385317021SBarry Smith 
23485317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
23585317021SBarry Smith @*/
23685317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
23785317021SBarry Smith {
23885317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
23985317021SBarry Smith 
24085317021SBarry Smith   PetscFunctionBegin;
2410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
242*c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,rtol,2);
24385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
24485317021SBarry Smith   if (f) {
24585317021SBarry Smith     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
24685317021SBarry Smith   }
24785317021SBarry Smith   PetscFunctionReturn(0);
24885317021SBarry Smith }
24985317021SBarry Smith 
25085317021SBarry Smith #undef __FUNCT__
25185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
252bf6011e8SBarry Smith /*@C
25385317021SBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
25485317021SBarry Smith 
255ad4df100SBarry Smith    Logically Collective on PC
25685317021SBarry Smith 
25785317021SBarry Smith    Input Parameters:
25885317021SBarry Smith +  pc - the preconditioner context
25985317021SBarry Smith -  stype - for example, spooles, superlu, superlu_dist
26085317021SBarry Smith 
26185317021SBarry Smith    Options Database Key:
26285317021SBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
26385317021SBarry Smith 
26485317021SBarry Smith    Level: intermediate
26585317021SBarry Smith 
26685317021SBarry Smith    Note:
26785317021SBarry Smith      By default this will use the PETSc factorization if it exists
26885317021SBarry Smith 
26985317021SBarry Smith 
27085317021SBarry Smith .keywords: PC, set, factorization, direct, fill
27185317021SBarry Smith 
2727112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
27385317021SBarry Smith 
27485317021SBarry Smith @*/
27585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
27685317021SBarry Smith {
27785317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
27885317021SBarry Smith 
27985317021SBarry Smith   PetscFunctionBegin;
2800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
28185317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
28285317021SBarry Smith   if (f) {
28385317021SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
28485317021SBarry Smith   }
28585317021SBarry Smith   PetscFunctionReturn(0);
28685317021SBarry Smith }
28785317021SBarry Smith 
28885317021SBarry Smith #undef __FUNCT__
2897112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage"
290bf6011e8SBarry Smith /*@C
2917112b564SBarry Smith    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
2927112b564SBarry Smith 
293*c5eb9154SBarry Smith    Not Collective
2947112b564SBarry Smith 
2957112b564SBarry Smith    Input Parameter:
2967112b564SBarry Smith .  pc - the preconditioner context
2977112b564SBarry Smith 
2987112b564SBarry Smith    Output Parameter:
2997112b564SBarry Smith .   stype - for example, spooles, superlu, superlu_dist
3007112b564SBarry Smith 
3017112b564SBarry Smith    Level: intermediate
3027112b564SBarry Smith 
3037112b564SBarry Smith 
3047112b564SBarry Smith .keywords: PC, set, factorization, direct, fill
3057112b564SBarry Smith 
3067112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
3077112b564SBarry Smith 
3087112b564SBarry Smith @*/
3097112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
3107112b564SBarry Smith {
3117112b564SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
3127112b564SBarry Smith 
3137112b564SBarry Smith   PetscFunctionBegin;
3140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3157112b564SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
3167112b564SBarry Smith   if (f) {
3177112b564SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
3187112b564SBarry Smith   }
3197112b564SBarry Smith   PetscFunctionReturn(0);
3207112b564SBarry Smith }
3217112b564SBarry Smith 
3227112b564SBarry Smith #undef __FUNCT__
32385317021SBarry Smith #define __FUNCT__ "PCFactorSetFill"
32485317021SBarry Smith /*@
32585317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
32685317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
32785317021SBarry Smith 
328*c5eb9154SBarry Smith    Not Collective, each process can expect a different amount of fill
32985317021SBarry Smith 
33085317021SBarry Smith    Input Parameters:
33185317021SBarry Smith +  pc - the preconditioner context
33285317021SBarry Smith -  fill - amount of expected fill
33385317021SBarry Smith 
33485317021SBarry Smith    Options Database Key:
33585317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
33685317021SBarry Smith 
33785317021SBarry Smith    Level: intermediate
33885317021SBarry Smith 
33985317021SBarry Smith    Note:
34085317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
34185317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
34285317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
34385317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
34485317021SBarry Smith 
34585317021SBarry Smith .keywords: PC, set, factorization, direct, fill
34685317021SBarry Smith 
34785317021SBarry Smith @*/
34885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
34985317021SBarry Smith {
35085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
35185317021SBarry Smith 
35285317021SBarry Smith   PetscFunctionBegin;
3530700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
354e7e72b3dSBarry Smith   if (fill < 1.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
35585317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
35685317021SBarry Smith   if (f) {
35785317021SBarry Smith     ierr = (*f)(pc,fill);CHKERRQ(ierr);
35885317021SBarry Smith   }
35985317021SBarry Smith   PetscFunctionReturn(0);
36085317021SBarry Smith }
36185317021SBarry Smith 
36285317021SBarry Smith #undef __FUNCT__
36385317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
36485317021SBarry Smith /*@
36585317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
36685317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
36785317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
36885317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
36985317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
37085317021SBarry Smith    space.
37185317021SBarry Smith 
372ad4df100SBarry Smith    Logically Collective on PC
37385317021SBarry Smith 
37485317021SBarry Smith    Input Parameters:
37585317021SBarry Smith .  pc - the preconditioner context
37685317021SBarry Smith 
37785317021SBarry Smith    Options Database Key:
37885317021SBarry Smith .  -pc_factor_in_place - Activates in-place factorization
37985317021SBarry Smith 
38085317021SBarry Smith    Notes:
38185317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
38285317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
38385317021SBarry Smith    a call to KSPSetOperators().
38485317021SBarry Smith    This is because the Krylov space methods require an application of the
38585317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
38685317021SBarry Smith    been factored in-place, replacing the original matrix.
38785317021SBarry Smith 
38885317021SBarry Smith    Level: intermediate
38985317021SBarry Smith 
39085317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
39185317021SBarry Smith 
39285317021SBarry Smith .seealso: PCILUSetUseInPlace()
39385317021SBarry Smith @*/
39485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
39585317021SBarry Smith {
39685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
39785317021SBarry Smith 
39885317021SBarry Smith   PetscFunctionBegin;
3990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
40085317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
40185317021SBarry Smith   if (f) {
40285317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
40385317021SBarry Smith   }
40485317021SBarry Smith   PetscFunctionReturn(0);
40585317021SBarry Smith }
40685317021SBarry Smith 
40785317021SBarry Smith #undef __FUNCT__
40885317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
40985317021SBarry Smith /*@C
41085317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
41185317021SBarry Smith     be used in the LU factorization.
41285317021SBarry Smith 
413ad4df100SBarry Smith     Logically Collective on PC
41485317021SBarry Smith 
41585317021SBarry Smith     Input Parameters:
41685317021SBarry Smith +   pc - the preconditioner context
4172692d6eeSBarry Smith -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
41885317021SBarry Smith 
41985317021SBarry Smith     Options Database Key:
42085317021SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
42185317021SBarry Smith 
42285317021SBarry Smith     Level: intermediate
42385317021SBarry Smith 
42485317021SBarry Smith     Notes: nested dissection is used by default
42585317021SBarry Smith 
42685317021SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
42785317021SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
42885317021SBarry Smith     SeqAIJ format in this case to get reorderings.
42985317021SBarry Smith 
43085317021SBarry Smith @*/
43185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
43285317021SBarry Smith {
43385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
43485317021SBarry Smith 
43585317021SBarry Smith   PetscFunctionBegin;
436*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
437*c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,ordering,2);
43885317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
43985317021SBarry Smith   if (f) {
44085317021SBarry Smith     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
44185317021SBarry Smith   }
44285317021SBarry Smith   PetscFunctionReturn(0);
44385317021SBarry Smith }
44485317021SBarry Smith 
44585317021SBarry Smith #undef __FUNCT__
4468ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot"
44785317021SBarry Smith /*@
4488ff23777SHong Zhang     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
44985317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
45085317021SBarry Smith       it is never done. For the Matlab and SuperLU factorization this is used.
45185317021SBarry Smith 
452ad4df100SBarry Smith     Logically Collective on PC
45385317021SBarry Smith 
45485317021SBarry Smith     Input Parameters:
45585317021SBarry Smith +   pc - the preconditioner context
45685317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
45785317021SBarry Smith 
45885317021SBarry Smith     Options Database Key:
45985317021SBarry Smith .   -pc_factor_pivoting <dtcol>
46085317021SBarry Smith 
46185317021SBarry Smith     Level: intermediate
46285317021SBarry Smith 
46385317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
46485317021SBarry Smith @*/
4658ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
46685317021SBarry Smith {
46785317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
46885317021SBarry Smith 
46985317021SBarry Smith   PetscFunctionBegin;
470*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
471*c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,dtcol,2);
4728ff23777SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
47385317021SBarry Smith   if (f) {
47485317021SBarry Smith     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
47585317021SBarry Smith   }
47685317021SBarry Smith   PetscFunctionReturn(0);
47785317021SBarry Smith }
47885317021SBarry Smith 
47985317021SBarry Smith #undef __FUNCT__
48085317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
48185317021SBarry Smith /*@
48285317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
48385317021SBarry Smith       with BAIJ or SBAIJ matrices
48485317021SBarry Smith 
485ad4df100SBarry Smith     Logically Collective on PC
48685317021SBarry Smith 
48785317021SBarry Smith     Input Parameters:
48885317021SBarry Smith +   pc - the preconditioner context
48985317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
49085317021SBarry Smith 
49185317021SBarry Smith     Options Database Key:
49285317021SBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
49385317021SBarry Smith 
49485317021SBarry Smith     Level: intermediate
49585317021SBarry Smith 
4968ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
49785317021SBarry Smith @*/
49885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
49985317021SBarry Smith {
50085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
50185317021SBarry Smith 
50285317021SBarry Smith   PetscFunctionBegin;
503*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
504*c5eb9154SBarry Smith   PetscValidLogicalCollectiveTruth(pc,pivot,2);
50585317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
50685317021SBarry Smith   if (f) {
50785317021SBarry Smith     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
50885317021SBarry Smith   }
50985317021SBarry Smith   PetscFunctionReturn(0);
51085317021SBarry Smith }
51185317021SBarry Smith 
51285317021SBarry Smith #undef __FUNCT__
51385317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill"
51485317021SBarry Smith /*@
51585317021SBarry Smith    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
51685317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
51785317021SBarry Smith 
518ad4df100SBarry Smith    Logically Collective on PC
51985317021SBarry Smith 
52085317021SBarry Smith    Input Parameters:
52185317021SBarry Smith +  pc - the preconditioner context
52285317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
52385317021SBarry Smith 
52485317021SBarry Smith    Options Database Key:
52585317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
52685317021SBarry Smith 
52785317021SBarry Smith    Level: intermediate
52885317021SBarry Smith 
52985317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
53085317021SBarry Smith 
53185317021SBarry Smith .seealso: PCFactorSetReuseOrdering()
53285317021SBarry Smith @*/
53385317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
53485317021SBarry Smith {
53585317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
53685317021SBarry Smith 
53785317021SBarry Smith   PetscFunctionBegin;
5380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,2);
539*c5eb9154SBarry Smith   PetscValidLogicalCollectiveTruth(pc,flag,2);
54085317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
54185317021SBarry Smith   if (f) {
54285317021SBarry Smith     ierr = (*f)(pc,flag);CHKERRQ(ierr);
54385317021SBarry Smith   }
54485317021SBarry Smith   PetscFunctionReturn(0);
54585317021SBarry Smith }
546