xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision ad4df1005a67cf9575a5c4294f11fd96229b2bb4)
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 
10*ad4df100SBarry 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);
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"
40915743fcSHong Zhang /*@
41915743fcSHong Zhang    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
42915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
43915743fcSHong Zhang 
44*ad4df100SBarry Smith    Logically Collective on PC
45915743fcSHong Zhang 
46915743fcSHong Zhang    Input Parameters:
47915743fcSHong Zhang +  pc - the preconditioner context
48915743fcSHong Zhang -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
49915743fcSHong Zhang 
50915743fcSHong Zhang    Options Database Key:
51915743fcSHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
52915743fcSHong Zhang 
53915743fcSHong Zhang    Level: intermediate
54915743fcSHong Zhang 
55915743fcSHong Zhang .keywords: PC, set, factorization,
56915743fcSHong Zhang 
57915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
58915743fcSHong 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;
640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,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"
74915743fcSHong Zhang /*@
75915743fcSHong Zhang    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
76915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
77915743fcSHong Zhang 
78*ad4df100SBarry Smith    Logically Collective on PC
79915743fcSHong Zhang 
80915743fcSHong Zhang    Input Parameters:
81915743fcSHong Zhang +  pc - the preconditioner context
82915743fcSHong Zhang -  shiftamount - amount of shift
83915743fcSHong Zhang 
84915743fcSHong Zhang    Options Database Key:
85915743fcSHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
86915743fcSHong Zhang 
87915743fcSHong Zhang    Level: intermediate
88915743fcSHong Zhang 
89915743fcSHong Zhang .keywords: PC, set, factorization,
90915743fcSHong Zhang 
91915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
92915743fcSHong 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;
980700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,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__
107b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance"
10878fc6b22SHong Zhang /*
109b7c853c4SBarry Smith    PCFactorSetDropTolerance - The preconditioner will use an ILU
11078fc6b22SHong Zhang    based on a drop tolerance. (Under development)
11185317021SBarry Smith 
112*ad4df100SBarry Smith    Logically Collective on PC
11385317021SBarry Smith 
11485317021SBarry Smith    Input Parameters:
11585317021SBarry Smith +  pc - the preconditioner context
11685317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
11785317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
11885317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
11985317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
12085317021SBarry Smith 
12185317021SBarry Smith    Options Database Key:
122b7c853c4SBarry Smith .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
12385317021SBarry Smith 
12485317021SBarry Smith    Level: intermediate
12585317021SBarry Smith 
12685317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
12785317021SBarry Smith       matrix. We don't know how to compute reasonable values.
12885317021SBarry Smith 
12985317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU
13078fc6b22SHong Zhang */
131b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
13285317021SBarry Smith {
13385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
13485317021SBarry Smith 
13585317021SBarry Smith   PetscFunctionBegin;
1360700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
137b7c853c4SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
13885317021SBarry Smith   if (f) {
13985317021SBarry Smith     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
14085317021SBarry Smith   }
14185317021SBarry Smith   PetscFunctionReturn(0);
14285317021SBarry Smith }
14385317021SBarry Smith 
14485317021SBarry Smith #undef __FUNCT__
14585317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels"
14685317021SBarry Smith /*@
14785317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
14885317021SBarry Smith 
149*ad4df100SBarry Smith    Logically Collective on PC
15085317021SBarry Smith 
15185317021SBarry Smith    Input Parameters:
15285317021SBarry Smith +  pc - the preconditioner context
15385317021SBarry Smith -  levels - number of levels of fill
15485317021SBarry Smith 
15585317021SBarry Smith    Options Database Key:
15685317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
15785317021SBarry Smith 
15885317021SBarry Smith    Level: intermediate
15985317021SBarry Smith 
16085317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
16185317021SBarry Smith @*/
16285317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
16385317021SBarry Smith {
16485317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
16585317021SBarry Smith 
16685317021SBarry Smith   PetscFunctionBegin;
1670700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
168e7e72b3dSBarry Smith   if (levels < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
16985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
17085317021SBarry Smith   if (f) {
17185317021SBarry Smith     ierr = (*f)(pc,levels);CHKERRQ(ierr);
17285317021SBarry Smith   }
17385317021SBarry Smith   PetscFunctionReturn(0);
17485317021SBarry Smith }
17585317021SBarry Smith 
17685317021SBarry Smith #undef __FUNCT__
17785317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
17885317021SBarry Smith /*@
17985317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
18085317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
18185317021SBarry Smith 
182*ad4df100SBarry Smith    Logically Collective on PC
18385317021SBarry Smith 
18485317021SBarry Smith    Input Parameters:
18585317021SBarry Smith +  pc - the preconditioner context
18685317021SBarry Smith 
18785317021SBarry Smith    Options Database Key:
18885317021SBarry Smith .  -pc_factor_diagonal_fill
18985317021SBarry Smith 
19085317021SBarry Smith    Notes:
19185317021SBarry Smith    Does not apply with 0 fill.
19285317021SBarry Smith 
19385317021SBarry Smith    Level: intermediate
19485317021SBarry Smith 
19585317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
19685317021SBarry Smith @*/
19785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
19885317021SBarry Smith {
19985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
20085317021SBarry Smith 
20185317021SBarry Smith   PetscFunctionBegin;
2020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
20385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
20485317021SBarry Smith   if (f) {
20585317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
20685317021SBarry Smith   }
20785317021SBarry Smith   PetscFunctionReturn(0);
20885317021SBarry Smith }
20985317021SBarry Smith 
21085317021SBarry Smith #undef __FUNCT__
21185317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
21285317021SBarry Smith /*@
21385317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
21485317021SBarry Smith 
215*ad4df100SBarry Smith    Logically Collective on PC
21685317021SBarry Smith 
21785317021SBarry Smith    Input Parameters:
21885317021SBarry Smith +  pc - the preconditioner context
21985317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
22085317021SBarry Smith 
22185317021SBarry Smith    Options Database Key:
22285317021SBarry Smith .  -pc_factor_nonzeros_along_diagonal
22385317021SBarry Smith 
22485317021SBarry Smith    Level: intermediate
22585317021SBarry Smith 
22685317021SBarry Smith .keywords: PC, set, factorization, direct, fill
22785317021SBarry Smith 
22885317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
22985317021SBarry Smith @*/
23085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
23185317021SBarry Smith {
23285317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
23385317021SBarry Smith 
23485317021SBarry Smith   PetscFunctionBegin;
2350700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
23785317021SBarry Smith   if (f) {
23885317021SBarry Smith     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
23985317021SBarry Smith   }
24085317021SBarry Smith   PetscFunctionReturn(0);
24185317021SBarry Smith }
24285317021SBarry Smith 
24385317021SBarry Smith #undef __FUNCT__
24485317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
245bf6011e8SBarry Smith /*@C
24685317021SBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
24785317021SBarry Smith 
248*ad4df100SBarry Smith    Logically Collective on PC
24985317021SBarry Smith 
25085317021SBarry Smith    Input Parameters:
25185317021SBarry Smith +  pc - the preconditioner context
25285317021SBarry Smith -  stype - for example, spooles, superlu, superlu_dist
25385317021SBarry Smith 
25485317021SBarry Smith    Options Database Key:
25585317021SBarry Smith .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
25685317021SBarry Smith 
25785317021SBarry Smith    Level: intermediate
25885317021SBarry Smith 
25985317021SBarry Smith    Note:
26085317021SBarry Smith      By default this will use the PETSc factorization if it exists
26185317021SBarry Smith 
26285317021SBarry Smith 
26385317021SBarry Smith .keywords: PC, set, factorization, direct, fill
26485317021SBarry Smith 
2657112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
26685317021SBarry Smith 
26785317021SBarry Smith @*/
26885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
26985317021SBarry Smith {
27085317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
27185317021SBarry Smith 
27285317021SBarry Smith   PetscFunctionBegin;
2730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
27485317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
27585317021SBarry Smith   if (f) {
27685317021SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
27785317021SBarry Smith   }
27885317021SBarry Smith   PetscFunctionReturn(0);
27985317021SBarry Smith }
28085317021SBarry Smith 
28185317021SBarry Smith #undef __FUNCT__
2827112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage"
283bf6011e8SBarry Smith /*@C
2847112b564SBarry Smith    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
2857112b564SBarry Smith 
286*ad4df100SBarry Smith    Logically Collective on PC
2877112b564SBarry Smith 
2887112b564SBarry Smith    Input Parameter:
2897112b564SBarry Smith .  pc - the preconditioner context
2907112b564SBarry Smith 
2917112b564SBarry Smith    Output Parameter:
2927112b564SBarry Smith .   stype - for example, spooles, superlu, superlu_dist
2937112b564SBarry Smith 
2947112b564SBarry Smith    Level: intermediate
2957112b564SBarry Smith 
2967112b564SBarry Smith 
2977112b564SBarry Smith .keywords: PC, set, factorization, direct, fill
2987112b564SBarry Smith 
2997112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
3007112b564SBarry Smith 
3017112b564SBarry Smith @*/
3027112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
3037112b564SBarry Smith {
3047112b564SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
3057112b564SBarry Smith 
3067112b564SBarry Smith   PetscFunctionBegin;
3070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3087112b564SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
3097112b564SBarry Smith   if (f) {
3107112b564SBarry Smith     ierr = (*f)(pc,stype);CHKERRQ(ierr);
3117112b564SBarry Smith   }
3127112b564SBarry Smith   PetscFunctionReturn(0);
3137112b564SBarry Smith }
3147112b564SBarry Smith 
3157112b564SBarry Smith #undef __FUNCT__
31685317021SBarry Smith #define __FUNCT__ "PCFactorSetFill"
31785317021SBarry Smith /*@
31885317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
31985317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
32085317021SBarry Smith 
321*ad4df100SBarry Smith    Logically Collective on PC
32285317021SBarry Smith 
32385317021SBarry Smith    Input Parameters:
32485317021SBarry Smith +  pc - the preconditioner context
32585317021SBarry Smith -  fill - amount of expected fill
32685317021SBarry Smith 
32785317021SBarry Smith    Options Database Key:
32885317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
32985317021SBarry Smith 
33085317021SBarry Smith    Level: intermediate
33185317021SBarry Smith 
33285317021SBarry Smith    Note:
33385317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
33485317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
33585317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
33685317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
33785317021SBarry Smith 
33885317021SBarry Smith .keywords: PC, set, factorization, direct, fill
33985317021SBarry Smith 
34085317021SBarry Smith @*/
34185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
34285317021SBarry Smith {
34385317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
34485317021SBarry Smith 
34585317021SBarry Smith   PetscFunctionBegin;
3460700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
347e7e72b3dSBarry Smith   if (fill < 1.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
34885317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
34985317021SBarry Smith   if (f) {
35085317021SBarry Smith     ierr = (*f)(pc,fill);CHKERRQ(ierr);
35185317021SBarry Smith   }
35285317021SBarry Smith   PetscFunctionReturn(0);
35385317021SBarry Smith }
35485317021SBarry Smith 
35585317021SBarry Smith #undef __FUNCT__
35685317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
35785317021SBarry Smith /*@
35885317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
35985317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
36085317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
36185317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
36285317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
36385317021SBarry Smith    space.
36485317021SBarry Smith 
365*ad4df100SBarry Smith    Logically Collective on PC
36685317021SBarry Smith 
36785317021SBarry Smith    Input Parameters:
36885317021SBarry Smith .  pc - the preconditioner context
36985317021SBarry Smith 
37085317021SBarry Smith    Options Database Key:
37185317021SBarry Smith .  -pc_factor_in_place - Activates in-place factorization
37285317021SBarry Smith 
37385317021SBarry Smith    Notes:
37485317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
37585317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
37685317021SBarry Smith    a call to KSPSetOperators().
37785317021SBarry Smith    This is because the Krylov space methods require an application of the
37885317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
37985317021SBarry Smith    been factored in-place, replacing the original matrix.
38085317021SBarry Smith 
38185317021SBarry Smith    Level: intermediate
38285317021SBarry Smith 
38385317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
38485317021SBarry Smith 
38585317021SBarry Smith .seealso: PCILUSetUseInPlace()
38685317021SBarry Smith @*/
38785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
38885317021SBarry Smith {
38985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC);
39085317021SBarry Smith 
39185317021SBarry Smith   PetscFunctionBegin;
3920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
39385317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
39485317021SBarry Smith   if (f) {
39585317021SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
39685317021SBarry Smith   }
39785317021SBarry Smith   PetscFunctionReturn(0);
39885317021SBarry Smith }
39985317021SBarry Smith 
40085317021SBarry Smith #undef __FUNCT__
40185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
40285317021SBarry Smith /*@C
40385317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
40485317021SBarry Smith     be used in the LU factorization.
40585317021SBarry Smith 
406*ad4df100SBarry Smith     Logically Collective on PC
40785317021SBarry Smith 
40885317021SBarry Smith     Input Parameters:
40985317021SBarry Smith +   pc - the preconditioner context
4102692d6eeSBarry Smith -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
41185317021SBarry Smith 
41285317021SBarry Smith     Options Database Key:
41385317021SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
41485317021SBarry Smith 
41585317021SBarry Smith     Level: intermediate
41685317021SBarry Smith 
41785317021SBarry Smith     Notes: nested dissection is used by default
41885317021SBarry Smith 
41985317021SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
42085317021SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
42185317021SBarry Smith     SeqAIJ format in this case to get reorderings.
42285317021SBarry Smith 
42385317021SBarry Smith @*/
42485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
42585317021SBarry Smith {
42685317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
42785317021SBarry Smith 
42885317021SBarry Smith   PetscFunctionBegin;
42985317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
43085317021SBarry Smith   if (f) {
43185317021SBarry Smith     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
43285317021SBarry Smith   }
43385317021SBarry Smith   PetscFunctionReturn(0);
43485317021SBarry Smith }
43585317021SBarry Smith 
43685317021SBarry Smith #undef __FUNCT__
4378ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot"
43885317021SBarry Smith /*@
4398ff23777SHong Zhang     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
44085317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
44185317021SBarry Smith       it is never done. For the Matlab and SuperLU factorization this is used.
44285317021SBarry Smith 
443*ad4df100SBarry Smith     Logically Collective on PC
44485317021SBarry Smith 
44585317021SBarry Smith     Input Parameters:
44685317021SBarry Smith +   pc - the preconditioner context
44785317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
44885317021SBarry Smith 
44985317021SBarry Smith     Options Database Key:
45085317021SBarry Smith .   -pc_factor_pivoting <dtcol>
45185317021SBarry Smith 
45285317021SBarry Smith     Level: intermediate
45385317021SBarry Smith 
45485317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
45585317021SBarry Smith @*/
4568ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
45785317021SBarry Smith {
45885317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
45985317021SBarry Smith 
46085317021SBarry Smith   PetscFunctionBegin;
4618ff23777SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
46285317021SBarry Smith   if (f) {
46385317021SBarry Smith     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
46485317021SBarry Smith   }
46585317021SBarry Smith   PetscFunctionReturn(0);
46685317021SBarry Smith }
46785317021SBarry Smith 
46885317021SBarry Smith #undef __FUNCT__
46985317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
47085317021SBarry Smith /*@
47185317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
47285317021SBarry Smith       with BAIJ or SBAIJ matrices
47385317021SBarry Smith 
474*ad4df100SBarry Smith     Logically Collective on PC
47585317021SBarry Smith 
47685317021SBarry Smith     Input Parameters:
47785317021SBarry Smith +   pc - the preconditioner context
47885317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
47985317021SBarry Smith 
48085317021SBarry Smith     Options Database Key:
48185317021SBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
48285317021SBarry Smith 
48385317021SBarry Smith     Level: intermediate
48485317021SBarry Smith 
4858ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
48685317021SBarry Smith @*/
48785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
48885317021SBarry Smith {
48985317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
49085317021SBarry Smith 
49185317021SBarry Smith   PetscFunctionBegin;
49285317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
49385317021SBarry Smith   if (f) {
49485317021SBarry Smith     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
49585317021SBarry Smith   }
49685317021SBarry Smith   PetscFunctionReturn(0);
49785317021SBarry Smith }
49885317021SBarry Smith 
49985317021SBarry Smith #undef __FUNCT__
50085317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill"
50185317021SBarry Smith /*@
50285317021SBarry Smith    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
50385317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
50485317021SBarry Smith 
505*ad4df100SBarry Smith    Logically Collective on PC
50685317021SBarry Smith 
50785317021SBarry Smith    Input Parameters:
50885317021SBarry Smith +  pc - the preconditioner context
50985317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
51085317021SBarry Smith 
51185317021SBarry Smith    Options Database Key:
51285317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
51385317021SBarry Smith 
51485317021SBarry Smith    Level: intermediate
51585317021SBarry Smith 
51685317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
51785317021SBarry Smith 
51885317021SBarry Smith .seealso: PCFactorSetReuseOrdering()
51985317021SBarry Smith @*/
52085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
52185317021SBarry Smith {
52285317021SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscTruth);
52385317021SBarry Smith 
52485317021SBarry Smith   PetscFunctionBegin;
5250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,2);
52685317021SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
52785317021SBarry Smith   if (f) {
52885317021SBarry Smith     ierr = (*f)(pc,flag);CHKERRQ(ierr);
52985317021SBarry Smith   }
53085317021SBarry Smith   PetscFunctionReturn(0);
53185317021SBarry Smith }
532