xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 288e7d534a092f6d31bbec1ea31d944215d65538)
15e8efad8SHong Zhang 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h>  /*I "petscpc.h" I*/
35e8efad8SHong Zhang 
4ee45ca4aSHong Zhang #undef __FUNCT__
53d1c1ea0SBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Factor"
63d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag)
73d1c1ea0SBarry Smith {
83d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
93d1c1ea0SBarry Smith 
103d1c1ea0SBarry Smith   PetscFunctionBegin;
113d1c1ea0SBarry Smith   lu->reuseordering = flag;
123d1c1ea0SBarry Smith   PetscFunctionReturn(0);
133d1c1ea0SBarry Smith }
143d1c1ea0SBarry Smith 
153d1c1ea0SBarry Smith #undef __FUNCT__
163d1c1ea0SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Factor"
173d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc,PetscBool flag)
183d1c1ea0SBarry Smith {
193d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
203d1c1ea0SBarry Smith 
213d1c1ea0SBarry Smith   PetscFunctionBegin;
223d1c1ea0SBarry Smith   lu->reusefill = flag;
233d1c1ea0SBarry Smith   PetscFunctionReturn(0);
243d1c1ea0SBarry Smith }
253d1c1ea0SBarry Smith 
263d1c1ea0SBarry Smith #undef __FUNCT__
273d1c1ea0SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Factor"
283d1c1ea0SBarry Smith static PetscErrorCode  PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg)
293d1c1ea0SBarry Smith {
303d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
313d1c1ea0SBarry Smith 
323d1c1ea0SBarry Smith   PetscFunctionBegin;
333d1c1ea0SBarry Smith   dir->inplace = flg;
343d1c1ea0SBarry Smith   PetscFunctionReturn(0);
353d1c1ea0SBarry Smith }
363d1c1ea0SBarry Smith 
373d1c1ea0SBarry Smith #undef __FUNCT__
383d1c1ea0SBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_Factor"
393d1c1ea0SBarry Smith static PetscErrorCode  PCFactorGetUseInPlace_Factor(PC pc,PetscBool *flg)
403d1c1ea0SBarry Smith {
413d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
423d1c1ea0SBarry Smith 
433d1c1ea0SBarry Smith   PetscFunctionBegin;
443d1c1ea0SBarry Smith   *flg = dir->inplace;
453d1c1ea0SBarry Smith   PetscFunctionReturn(0);
463d1c1ea0SBarry Smith }
473d1c1ea0SBarry Smith 
483d1c1ea0SBarry Smith #undef __FUNCT__
49f8260c8fSBarry Smith #define __FUNCT__ "PCFactorSetUpMatSolverPackage"
50f8260c8fSBarry Smith /*@
51f8260c8fSBarry Smith     PCFactorSetUpMatSolverPackage - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
52f8260c8fSBarry Smith        set the options for that particular factorization object.
53f8260c8fSBarry Smith 
54f8260c8fSBarry Smith   Input Parameter:
55f8260c8fSBarry Smith .  pc  - the preconditioner context
56f8260c8fSBarry Smith 
57f8260c8fSBarry Smith   Notes: After you have called this function (which has to be after the KSPSetOperators() or PCSetOperators()) you can call PCFactorGetMatrix() and then set factor options on that matrix.
58f8260c8fSBarry Smith 
59f8260c8fSBarry Smith .seealso: PCFactorSetMatSolverPackage(), PCFactorGetMatrix()
60f8260c8fSBarry Smith 
612bd2b0e6SSatish Balay   Level: intermediate
622bd2b0e6SSatish Balay 
63f8260c8fSBarry Smith @*/
64f8260c8fSBarry Smith PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
65f8260c8fSBarry Smith {
66f8260c8fSBarry Smith   PetscErrorCode ierr;
67f8260c8fSBarry Smith 
68f8260c8fSBarry Smith   PetscFunctionBegin;
69f8260c8fSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
70f8260c8fSBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetUpMatSolverPackage_C",(PC),(pc));CHKERRQ(ierr);
71b3a44c85SBarry Smith   PetscFunctionReturn(0);
72f8260c8fSBarry Smith }
73f8260c8fSBarry Smith 
74f8260c8fSBarry Smith #undef __FUNCT__
75ee45ca4aSHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot"
76ee45ca4aSHong Zhang /*@
77ee45ca4aSHong Zhang    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
78ee45ca4aSHong Zhang 
79ad4df100SBarry Smith    Logically Collective on PC
80ee45ca4aSHong Zhang 
81ee45ca4aSHong Zhang    Input Parameters:
82afaefe49SHong Zhang +  pc - the preconditioner context
83afaefe49SHong Zhang -  zero - all pivots smaller than this will be considered zero
84ee45ca4aSHong Zhang 
85ee45ca4aSHong Zhang    Options Database Key:
86ee45ca4aSHong Zhang .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
87ee45ca4aSHong Zhang 
88ee45ca4aSHong Zhang    Level: intermediate
89ee45ca4aSHong Zhang 
90ee45ca4aSHong Zhang .keywords: PC, set, factorization, direct, fill
91ee45ca4aSHong Zhang 
92daa17b54SHong Zhang .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
93ee45ca4aSHong Zhang @*/
947087cfbeSBarry Smith PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
95ee45ca4aSHong Zhang {
964ac538c5SBarry Smith   PetscErrorCode ierr;
97afaefe49SHong Zhang 
98ee45ca4aSHong Zhang   PetscFunctionBegin;
990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
100c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,zero,2);
1014ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));CHKERRQ(ierr);
102ee45ca4aSHong Zhang   PetscFunctionReturn(0);
103ee45ca4aSHong Zhang }
104ee45ca4aSHong Zhang 
1055e8efad8SHong Zhang #undef __FUNCT__
106d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType"
107915743fcSHong Zhang /*@
108915743fcSHong Zhang    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
109915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
110915743fcSHong Zhang 
111ad4df100SBarry Smith    Logically Collective on PC
112915743fcSHong Zhang 
113915743fcSHong Zhang    Input Parameters:
114915743fcSHong Zhang +  pc - the preconditioner context
115915743fcSHong Zhang -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
116915743fcSHong Zhang 
117915743fcSHong Zhang    Options Database Key:
118915743fcSHong Zhang .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
119915743fcSHong Zhang 
120915743fcSHong Zhang    Level: intermediate
121915743fcSHong Zhang 
122915743fcSHong Zhang .keywords: PC, set, factorization,
123915743fcSHong Zhang 
124915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
125915743fcSHong Zhang @*/
1267087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
127d90ac83dSHong Zhang {
1284ac538c5SBarry Smith   PetscErrorCode ierr;
129d90ac83dSHong Zhang 
130d90ac83dSHong Zhang   PetscFunctionBegin;
1310700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
132c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,shifttype,2);
1334ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));CHKERRQ(ierr);
134d90ac83dSHong Zhang   PetscFunctionReturn(0);
135d90ac83dSHong Zhang }
136d90ac83dSHong Zhang 
137d90ac83dSHong Zhang #undef __FUNCT__
138d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount"
139915743fcSHong Zhang /*@
140915743fcSHong Zhang    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
141915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
142915743fcSHong Zhang 
143ad4df100SBarry Smith    Logically Collective on PC
144915743fcSHong Zhang 
145915743fcSHong Zhang    Input Parameters:
146915743fcSHong Zhang +  pc - the preconditioner context
147915743fcSHong Zhang -  shiftamount - amount of shift
148915743fcSHong Zhang 
149915743fcSHong Zhang    Options Database Key:
150915743fcSHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
151915743fcSHong Zhang 
152915743fcSHong Zhang    Level: intermediate
153915743fcSHong Zhang 
154915743fcSHong Zhang .keywords: PC, set, factorization,
155915743fcSHong Zhang 
156915743fcSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
157915743fcSHong Zhang @*/
1587087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
159d90ac83dSHong Zhang {
1604ac538c5SBarry Smith   PetscErrorCode ierr;
161d90ac83dSHong Zhang 
162d90ac83dSHong Zhang   PetscFunctionBegin;
1630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
164c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,shiftamount,2);
1654ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));CHKERRQ(ierr);
166d90ac83dSHong Zhang   PetscFunctionReturn(0);
167d90ac83dSHong Zhang }
168d90ac83dSHong Zhang 
169d90ac83dSHong Zhang #undef __FUNCT__
170b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance"
17178fc6b22SHong Zhang /*
172b7c853c4SBarry Smith    PCFactorSetDropTolerance - The preconditioner will use an ILU
17378fc6b22SHong Zhang    based on a drop tolerance. (Under development)
17485317021SBarry Smith 
175ad4df100SBarry Smith    Logically Collective on PC
17685317021SBarry Smith 
17785317021SBarry Smith    Input Parameters:
17885317021SBarry Smith +  pc - the preconditioner context
17985317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
18085317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
18185317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
18285317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
18385317021SBarry Smith 
18485317021SBarry Smith    Options Database Key:
185b7c853c4SBarry Smith .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
18685317021SBarry Smith 
18785317021SBarry Smith    Level: intermediate
18885317021SBarry Smith 
18985317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
19085317021SBarry Smith       matrix. We don't know how to compute reasonable values.
19185317021SBarry Smith 
19285317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, ILU
19378fc6b22SHong Zhang */
1947087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
19585317021SBarry Smith {
1964ac538c5SBarry Smith   PetscErrorCode ierr;
19785317021SBarry Smith 
19885317021SBarry Smith   PetscFunctionBegin;
1990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
200c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,dtcol,2);
201c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,maxrowcount,3);
2024ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));CHKERRQ(ierr);
20385317021SBarry Smith   PetscFunctionReturn(0);
20485317021SBarry Smith }
20585317021SBarry Smith 
20685317021SBarry Smith #undef __FUNCT__
207c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetZeroPivot"
208c7f610a1SBarry Smith /*@
209c7f610a1SBarry Smith    PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
210c7f610a1SBarry Smith 
211c7f610a1SBarry Smith    Not Collective
212c7f610a1SBarry Smith 
213c7f610a1SBarry Smith    Input Parameters:
214c7f610a1SBarry Smith .  pc - the preconditioner context
215c7f610a1SBarry Smith 
216c7f610a1SBarry Smith    Output Parameter:
217c7f610a1SBarry Smith .  pivot - the tolerance
218c7f610a1SBarry Smith 
219c7f610a1SBarry Smith    Level: intermediate
220c7f610a1SBarry Smith 
221c7f610a1SBarry Smith 
222c7f610a1SBarry Smith .seealso: PCFactorSetZeroPivot()
223c7f610a1SBarry Smith @*/
224c7f610a1SBarry Smith PetscErrorCode  PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
225c7f610a1SBarry Smith {
226c7f610a1SBarry Smith   PetscErrorCode ierr;
227c7f610a1SBarry Smith 
228c7f610a1SBarry Smith   PetscFunctionBegin;
229c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
230c7f610a1SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));CHKERRQ(ierr);
231c7f610a1SBarry Smith   PetscFunctionReturn(0);
232c7f610a1SBarry Smith }
233c7f610a1SBarry Smith 
234c7f610a1SBarry Smith #undef __FUNCT__
235c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetShiftAmount"
236c7f610a1SBarry Smith /*@
237c7f610a1SBarry Smith    PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
238c7f610a1SBarry Smith 
239c7f610a1SBarry Smith    Not Collective
240c7f610a1SBarry Smith 
241c7f610a1SBarry Smith    Input Parameters:
242c7f610a1SBarry Smith .  pc - the preconditioner context
243c7f610a1SBarry Smith 
244c7f610a1SBarry Smith    Output Parameter:
245c7f610a1SBarry Smith .  shift - how much to shift the diagonal entry
246c7f610a1SBarry Smith 
247c7f610a1SBarry Smith    Level: intermediate
248c7f610a1SBarry Smith 
249c7f610a1SBarry Smith 
250c7f610a1SBarry Smith .seealso: PCFactorSetShiftAmount(), PCFactorSetShiftType(), PCFactorGetShiftType()
251c7f610a1SBarry Smith @*/
252c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftAmount(PC pc,PetscReal *shift)
253c7f610a1SBarry Smith {
254c7f610a1SBarry Smith   PetscErrorCode ierr;
255c7f610a1SBarry Smith 
256c7f610a1SBarry Smith   PetscFunctionBegin;
257c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
258c7f610a1SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));CHKERRQ(ierr);
259c7f610a1SBarry Smith   PetscFunctionReturn(0);
260c7f610a1SBarry Smith }
261c7f610a1SBarry Smith 
262c7f610a1SBarry Smith #undef __FUNCT__
263c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetShiftType"
264c7f610a1SBarry Smith /*@
265c7f610a1SBarry Smith    PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
266c7f610a1SBarry Smith 
267c7f610a1SBarry Smith    Not Collective
268c7f610a1SBarry Smith 
269c7f610a1SBarry Smith    Input Parameters:
270c7f610a1SBarry Smith .  pc - the preconditioner context
271c7f610a1SBarry Smith 
272c7f610a1SBarry Smith    Output Parameter:
273c7f610a1SBarry Smith .  type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS
274c7f610a1SBarry Smith 
275c7f610a1SBarry Smith    Level: intermediate
276c7f610a1SBarry Smith 
277c7f610a1SBarry Smith 
278c7f610a1SBarry Smith .seealso: PCFactorSetShiftType(), MatFactorShiftType, PCFactorSetShiftAmount(), PCFactorGetShiftAmount()
279c7f610a1SBarry Smith @*/
280c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
281c7f610a1SBarry Smith {
282c7f610a1SBarry Smith   PetscErrorCode ierr;
283c7f610a1SBarry Smith 
284c7f610a1SBarry Smith   PetscFunctionBegin;
285c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
286c7f610a1SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));CHKERRQ(ierr);
287c7f610a1SBarry Smith   PetscFunctionReturn(0);
288c7f610a1SBarry Smith }
289c7f610a1SBarry Smith 
290c7f610a1SBarry Smith #undef __FUNCT__
2912591b318SToby Isaac #define __FUNCT__ "PCFactorGetLevels"
2922591b318SToby Isaac /*@
2932591b318SToby Isaac    PCFactorGetLevels - Gets the number of levels of fill to use.
2942591b318SToby Isaac 
2952591b318SToby Isaac    Logically Collective on PC
2962591b318SToby Isaac 
2972591b318SToby Isaac    Input Parameters:
2982591b318SToby Isaac .  pc - the preconditioner context
2992591b318SToby Isaac 
3002591b318SToby Isaac    Output Parameter:
3012591b318SToby Isaac .  levels - number of levels of fill
3022591b318SToby Isaac 
3032591b318SToby Isaac    Level: intermediate
3042591b318SToby Isaac 
3052591b318SToby Isaac .keywords: PC, levels, fill, factorization, incomplete, ILU
3062591b318SToby Isaac @*/
3072591b318SToby Isaac PetscErrorCode  PCFactorGetLevels(PC pc,PetscInt *levels)
3082591b318SToby Isaac {
3092591b318SToby Isaac   PetscErrorCode ierr;
3102591b318SToby Isaac 
3112591b318SToby Isaac   PetscFunctionBegin;
3122591b318SToby Isaac   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
313c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
3142591b318SToby Isaac   PetscFunctionReturn(0);
3152591b318SToby Isaac }
3162591b318SToby Isaac 
3172591b318SToby Isaac #undef __FUNCT__
31885317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels"
31985317021SBarry Smith /*@
32085317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
32185317021SBarry Smith 
322ad4df100SBarry Smith    Logically Collective on PC
32385317021SBarry Smith 
32485317021SBarry Smith    Input Parameters:
32585317021SBarry Smith +  pc - the preconditioner context
32685317021SBarry Smith -  levels - number of levels of fill
32785317021SBarry Smith 
32885317021SBarry Smith    Options Database Key:
32985317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
33085317021SBarry Smith 
33185317021SBarry Smith    Level: intermediate
33285317021SBarry Smith 
33385317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
33485317021SBarry Smith @*/
3357087cfbeSBarry Smith PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
33685317021SBarry Smith {
3374ac538c5SBarry Smith   PetscErrorCode ierr;
33885317021SBarry Smith 
33985317021SBarry Smith   PetscFunctionBegin;
3400700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
341ce94432eSBarry Smith   if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
342c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,levels,2);
3434ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
34485317021SBarry Smith   PetscFunctionReturn(0);
34585317021SBarry Smith }
34685317021SBarry Smith 
34785317021SBarry Smith #undef __FUNCT__
34885317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
34985317021SBarry Smith /*@
35085317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
35185317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
35285317021SBarry Smith 
353ad4df100SBarry Smith    Logically Collective on PC
35485317021SBarry Smith 
35585317021SBarry Smith    Input Parameters:
35685317021SBarry Smith +  pc - the preconditioner context
35792e9c092SBarry Smith -  flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
35885317021SBarry Smith 
35985317021SBarry Smith    Options Database Key:
36085317021SBarry Smith .  -pc_factor_diagonal_fill
36185317021SBarry Smith 
36285317021SBarry Smith    Notes:
36385317021SBarry Smith    Does not apply with 0 fill.
36485317021SBarry Smith 
36585317021SBarry Smith    Level: intermediate
36685317021SBarry Smith 
36785317021SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
36892e9c092SBarry Smith 
36992e9c092SBarry Smith .seealso: PCFactorGetAllowDiagonalFill()
37092e9c092SBarry Smith 
37185317021SBarry Smith @*/
37292e9c092SBarry Smith PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
37385317021SBarry Smith {
3744ac538c5SBarry Smith   PetscErrorCode ierr;
37585317021SBarry Smith 
37685317021SBarry Smith   PetscFunctionBegin;
3770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
37892e9c092SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
37992e9c092SBarry Smith   PetscFunctionReturn(0);
38092e9c092SBarry Smith }
38192e9c092SBarry Smith 
38292e9c092SBarry Smith #undef __FUNCT__
38392e9c092SBarry Smith #define __FUNCT__ "PCFactorGetAllowDiagonalFill"
38492e9c092SBarry Smith /*@
38592e9c092SBarry Smith    PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
38692e9c092SBarry Smith        treated as level 0 fill even if there is no non-zero location.
38792e9c092SBarry Smith 
38892e9c092SBarry Smith    Logically Collective on PC
38992e9c092SBarry Smith 
39092e9c092SBarry Smith    Input Parameter:
39192e9c092SBarry Smith .  pc - the preconditioner context
39292e9c092SBarry Smith 
39392e9c092SBarry Smith    Output Parameter:
39492e9c092SBarry Smith .   flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
39592e9c092SBarry Smith 
39692e9c092SBarry Smith    Options Database Key:
39792e9c092SBarry Smith .  -pc_factor_diagonal_fill
39892e9c092SBarry Smith 
39992e9c092SBarry Smith    Notes:
40092e9c092SBarry Smith    Does not apply with 0 fill.
40192e9c092SBarry Smith 
40292e9c092SBarry Smith    Level: intermediate
40392e9c092SBarry Smith 
40492e9c092SBarry Smith .keywords: PC, levels, fill, factorization, incomplete, ILU
40592e9c092SBarry Smith 
40692e9c092SBarry Smith .seealso: PCFactorSetAllowDiagonalFill()
40792e9c092SBarry Smith 
40892e9c092SBarry Smith @*/
40992e9c092SBarry Smith PetscErrorCode  PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
41092e9c092SBarry Smith {
41192e9c092SBarry Smith   PetscErrorCode ierr;
41292e9c092SBarry Smith 
41392e9c092SBarry Smith   PetscFunctionBegin;
41492e9c092SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
415c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
41685317021SBarry Smith   PetscFunctionReturn(0);
41785317021SBarry Smith }
41885317021SBarry Smith 
41985317021SBarry Smith #undef __FUNCT__
42085317021SBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
42185317021SBarry Smith /*@
42285317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
42385317021SBarry Smith 
424ad4df100SBarry Smith    Logically Collective on PC
42585317021SBarry Smith 
42685317021SBarry Smith    Input Parameters:
42785317021SBarry Smith +  pc - the preconditioner context
42885317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
42985317021SBarry Smith 
43085317021SBarry Smith    Options Database Key:
43192e9c092SBarry Smith .  -pc_factor_nonzeros_along_diagonal <tol>
43285317021SBarry Smith 
43385317021SBarry Smith    Level: intermediate
43485317021SBarry Smith 
43585317021SBarry Smith .keywords: PC, set, factorization, direct, fill
43685317021SBarry Smith 
43785317021SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
43885317021SBarry Smith @*/
4397087cfbeSBarry Smith PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
44085317021SBarry Smith {
4414ac538c5SBarry Smith   PetscErrorCode ierr;
44285317021SBarry Smith 
44385317021SBarry Smith   PetscFunctionBegin;
4440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
445c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,rtol,2);
4464ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));CHKERRQ(ierr);
44785317021SBarry Smith   PetscFunctionReturn(0);
44885317021SBarry Smith }
44985317021SBarry Smith 
45085317021SBarry Smith #undef __FUNCT__
45185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage"
452bf6011e8SBarry Smith /*@C
45385317021SBarry Smith    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
45485317021SBarry Smith 
455ad4df100SBarry Smith    Logically Collective on PC
45685317021SBarry Smith 
45785317021SBarry Smith    Input Parameters:
45885317021SBarry Smith +  pc - the preconditioner context
459f60c3dc2SHong Zhang -  stype - for example, superlu, superlu_dist
46085317021SBarry Smith 
46185317021SBarry Smith    Options Database Key:
462f60c3dc2SHong Zhang .  -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse
46385317021SBarry Smith 
46485317021SBarry Smith    Level: intermediate
46585317021SBarry Smith 
46685317021SBarry Smith    Note:
46785317021SBarry Smith      By default this will use the PETSc factorization if it exists
46885317021SBarry Smith 
46985317021SBarry Smith 
47085317021SBarry Smith .keywords: PC, set, factorization, direct, fill
47185317021SBarry Smith 
4727112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
47385317021SBarry Smith 
47485317021SBarry Smith @*/
4757087cfbeSBarry Smith PetscErrorCode  PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
47685317021SBarry Smith {
4774ac538c5SBarry Smith   PetscErrorCode ierr;
47885317021SBarry Smith 
47985317021SBarry Smith   PetscFunctionBegin;
4800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4814ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));CHKERRQ(ierr);
48285317021SBarry Smith   PetscFunctionReturn(0);
48385317021SBarry Smith }
48485317021SBarry Smith 
48585317021SBarry Smith #undef __FUNCT__
4867112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage"
487bf6011e8SBarry Smith /*@C
4887112b564SBarry Smith    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
4897112b564SBarry Smith 
490c5eb9154SBarry Smith    Not Collective
4917112b564SBarry Smith 
4927112b564SBarry Smith    Input Parameter:
4937112b564SBarry Smith .  pc - the preconditioner context
4947112b564SBarry Smith 
4957112b564SBarry Smith    Output Parameter:
4960298fd71SBarry Smith .   stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
4977112b564SBarry Smith 
4987112b564SBarry Smith    Level: intermediate
4997112b564SBarry Smith 
5007112b564SBarry Smith 
5017112b564SBarry Smith .keywords: PC, set, factorization, direct, fill
5027112b564SBarry Smith 
5037112b564SBarry Smith .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
5047112b564SBarry Smith 
5057112b564SBarry Smith @*/
5067087cfbeSBarry Smith PetscErrorCode  PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
5077112b564SBarry Smith {
5088b5c83b4SJed Brown   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
5097112b564SBarry Smith 
5107112b564SBarry Smith   PetscFunctionBegin;
5110700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5120005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);CHKERRQ(ierr);
5138b5c83b4SJed Brown   if (f) {
5148b5c83b4SJed Brown     ierr = (*f)(pc,stype);CHKERRQ(ierr);
5158b5c83b4SJed Brown   } else {
5160298fd71SBarry Smith     *stype = NULL;
5178b5c83b4SJed Brown   }
5187112b564SBarry Smith   PetscFunctionReturn(0);
5197112b564SBarry Smith }
5207112b564SBarry Smith 
5217112b564SBarry Smith #undef __FUNCT__
52285317021SBarry Smith #define __FUNCT__ "PCFactorSetFill"
52385317021SBarry Smith /*@
52485317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
52585317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
52685317021SBarry Smith 
527c5eb9154SBarry Smith    Not Collective, each process can expect a different amount of fill
52885317021SBarry Smith 
52985317021SBarry Smith    Input Parameters:
53085317021SBarry Smith +  pc - the preconditioner context
53185317021SBarry Smith -  fill - amount of expected fill
53285317021SBarry Smith 
53385317021SBarry Smith    Options Database Key:
53485317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
53585317021SBarry Smith 
53685317021SBarry Smith    Level: intermediate
53785317021SBarry Smith 
53885317021SBarry Smith    Note:
53985317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
54085317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
54185317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
54285317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
54385317021SBarry Smith 
54401a79839SBarry Smith    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
54501a79839SBarry Smith 
54601a79839SBarry Smith 
54785317021SBarry Smith .keywords: PC, set, factorization, direct, fill
54885317021SBarry Smith 
54985317021SBarry Smith @*/
5507087cfbeSBarry Smith PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
55185317021SBarry Smith {
5524ac538c5SBarry Smith   PetscErrorCode ierr;
55385317021SBarry Smith 
55485317021SBarry Smith   PetscFunctionBegin;
5550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
556ce94432eSBarry Smith   if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
5574ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));CHKERRQ(ierr);
55885317021SBarry Smith   PetscFunctionReturn(0);
55985317021SBarry Smith }
56085317021SBarry Smith 
56185317021SBarry Smith #undef __FUNCT__
56285317021SBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace"
56385317021SBarry Smith /*@
56485317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
56585317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
56685317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
56785317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
56885317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
56985317021SBarry Smith    space.
57085317021SBarry Smith 
571ad4df100SBarry Smith    Logically Collective on PC
57285317021SBarry Smith 
57385317021SBarry Smith    Input Parameters:
5748e37d05fSBarry Smith +  pc - the preconditioner context
5758e37d05fSBarry Smith -  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
57685317021SBarry Smith 
57785317021SBarry Smith    Options Database Key:
5788e37d05fSBarry Smith .  -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
57985317021SBarry Smith 
58085317021SBarry Smith    Notes:
58185317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
58285317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
58385317021SBarry Smith    a call to KSPSetOperators().
58485317021SBarry Smith    This is because the Krylov space methods require an application of the
58585317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
58685317021SBarry Smith    been factored in-place, replacing the original matrix.
58785317021SBarry Smith 
58885317021SBarry Smith    Level: intermediate
58985317021SBarry Smith 
59085317021SBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
59185317021SBarry Smith 
5928e37d05fSBarry Smith .seealso: PCFactorGetUseInPlace()
59385317021SBarry Smith @*/
5948e37d05fSBarry Smith PetscErrorCode  PCFactorSetUseInPlace(PC pc,PetscBool flg)
59585317021SBarry Smith {
5964ac538c5SBarry Smith   PetscErrorCode ierr;
59785317021SBarry Smith 
59885317021SBarry Smith   PetscFunctionBegin;
5990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
6008e37d05fSBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
6018e37d05fSBarry Smith   PetscFunctionReturn(0);
6028e37d05fSBarry Smith }
6038e37d05fSBarry Smith 
6048e37d05fSBarry Smith #undef __FUNCT__
6058e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace"
6068e37d05fSBarry Smith /*@
6078e37d05fSBarry Smith    PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
6088e37d05fSBarry Smith 
6098e37d05fSBarry Smith    Logically Collective on PC
6108e37d05fSBarry Smith 
6118e37d05fSBarry Smith    Input Parameter:
6128e37d05fSBarry Smith .  pc - the preconditioner context
6138e37d05fSBarry Smith 
6148e37d05fSBarry Smith    Output Parameter:
6158e37d05fSBarry Smith .  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
6168e37d05fSBarry Smith 
6178e37d05fSBarry Smith    Level: intermediate
6188e37d05fSBarry Smith 
6198e37d05fSBarry Smith .keywords: PC, set, factorization, direct, inplace, in-place, LU
6208e37d05fSBarry Smith 
6218e37d05fSBarry Smith .seealso: PCFactorSetUseInPlace()
6228e37d05fSBarry Smith @*/
6238e37d05fSBarry Smith PetscErrorCode  PCFactorGetUseInPlace(PC pc,PetscBool *flg)
6248e37d05fSBarry Smith {
6258e37d05fSBarry Smith   PetscErrorCode ierr;
6268e37d05fSBarry Smith 
6278e37d05fSBarry Smith   PetscFunctionBegin;
6288e37d05fSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
629c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
63085317021SBarry Smith   PetscFunctionReturn(0);
63185317021SBarry Smith }
63285317021SBarry Smith 
63385317021SBarry Smith #undef __FUNCT__
63485317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType"
63585317021SBarry Smith /*@C
63685317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
63785317021SBarry Smith     be used in the LU factorization.
63885317021SBarry Smith 
639ad4df100SBarry Smith     Logically Collective on PC
64085317021SBarry Smith 
64185317021SBarry Smith     Input Parameters:
64285317021SBarry Smith +   pc - the preconditioner context
6432692d6eeSBarry Smith -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
64485317021SBarry Smith 
64585317021SBarry Smith     Options Database Key:
64685317021SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
64785317021SBarry Smith 
64885317021SBarry Smith     Level: intermediate
64985317021SBarry Smith 
65085317021SBarry Smith     Notes: nested dissection is used by default
65185317021SBarry Smith 
65285317021SBarry Smith     For Cholesky and ICC and the SBAIJ format reorderings are not available,
65385317021SBarry Smith     since only the upper triangular part of the matrix is stored. You can use the
65485317021SBarry Smith     SeqAIJ format in this case to get reorderings.
65585317021SBarry Smith 
65685317021SBarry Smith @*/
65719fd82e9SBarry Smith PetscErrorCode  PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
65885317021SBarry Smith {
6594ac538c5SBarry Smith   PetscErrorCode ierr;
66085317021SBarry Smith 
66185317021SBarry Smith   PetscFunctionBegin;
662c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
66319fd82e9SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));CHKERRQ(ierr);
66485317021SBarry Smith   PetscFunctionReturn(0);
66585317021SBarry Smith }
66685317021SBarry Smith 
66785317021SBarry Smith #undef __FUNCT__
6688ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot"
66985317021SBarry Smith /*@
6708ff23777SHong Zhang     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
67185317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
672e3c5b3baSBarry Smith       it is never done. For the MATLAB and SuperLU factorization this is used.
67385317021SBarry Smith 
674ad4df100SBarry Smith     Logically Collective on PC
67585317021SBarry Smith 
67685317021SBarry Smith     Input Parameters:
67785317021SBarry Smith +   pc - the preconditioner context
67885317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
67985317021SBarry Smith 
68085317021SBarry Smith     Options Database Key:
68185317021SBarry Smith .   -pc_factor_pivoting <dtcol>
68285317021SBarry Smith 
68385317021SBarry Smith     Level: intermediate
68485317021SBarry Smith 
68585317021SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
68685317021SBarry Smith @*/
6877087cfbeSBarry Smith PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
68885317021SBarry Smith {
6894ac538c5SBarry Smith   PetscErrorCode ierr;
69085317021SBarry Smith 
69185317021SBarry Smith   PetscFunctionBegin;
692c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
693c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,dtcol,2);
6944ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));CHKERRQ(ierr);
69585317021SBarry Smith   PetscFunctionReturn(0);
69685317021SBarry Smith }
69785317021SBarry Smith 
69885317021SBarry Smith #undef __FUNCT__
69985317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks"
70085317021SBarry Smith /*@
70185317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
70285317021SBarry Smith       with BAIJ or SBAIJ matrices
70385317021SBarry Smith 
704ad4df100SBarry Smith     Logically Collective on PC
70585317021SBarry Smith 
70685317021SBarry Smith     Input Parameters:
70785317021SBarry Smith +   pc - the preconditioner context
70885317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
70985317021SBarry Smith 
71085317021SBarry Smith     Options Database Key:
71185317021SBarry Smith .   -pc_factor_pivot_in_blocks <true,false>
71285317021SBarry Smith 
71385317021SBarry Smith     Level: intermediate
71485317021SBarry Smith 
7158ff23777SHong Zhang .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
71685317021SBarry Smith @*/
7177087cfbeSBarry Smith PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
71885317021SBarry Smith {
7194ac538c5SBarry Smith   PetscErrorCode ierr;
72085317021SBarry Smith 
72185317021SBarry Smith   PetscFunctionBegin;
722c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
723acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,pivot,2);
7244ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));CHKERRQ(ierr);
72585317021SBarry Smith   PetscFunctionReturn(0);
72685317021SBarry Smith }
72785317021SBarry Smith 
72885317021SBarry Smith #undef __FUNCT__
72985317021SBarry Smith #define __FUNCT__ "PCFactorSetReuseFill"
73085317021SBarry Smith /*@
731*288e7d53SBarry Smith    PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
73285317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
73385317021SBarry Smith 
734ad4df100SBarry Smith    Logically Collective on PC
73585317021SBarry Smith 
73685317021SBarry Smith    Input Parameters:
73785317021SBarry Smith +  pc - the preconditioner context
73885317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
73985317021SBarry Smith 
74085317021SBarry Smith    Options Database Key:
74185317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
74285317021SBarry Smith 
74385317021SBarry Smith    Level: intermediate
74485317021SBarry Smith 
74585317021SBarry Smith .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
74685317021SBarry Smith 
74785317021SBarry Smith .seealso: PCFactorSetReuseOrdering()
74885317021SBarry Smith @*/
7497087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool flag)
75085317021SBarry Smith {
7514ac538c5SBarry Smith   PetscErrorCode ierr;
75285317021SBarry Smith 
75385317021SBarry Smith   PetscFunctionBegin;
7540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,2);
755acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
7564ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr);
75785317021SBarry Smith   PetscFunctionReturn(0);
75885317021SBarry Smith }
7593d1c1ea0SBarry Smith 
7603d1c1ea0SBarry Smith #undef __FUNCT__
7613d1c1ea0SBarry Smith #define __FUNCT__ "PCFactorInitialize"
7623d1c1ea0SBarry Smith PetscErrorCode PCFactorInitialize(PC pc)
7633d1c1ea0SBarry Smith {
7643d1c1ea0SBarry Smith   PetscErrorCode ierr;
7653d1c1ea0SBarry Smith   PC_Factor       *fact = (PC_Factor*)pc->data;
7663d1c1ea0SBarry Smith 
7673d1c1ea0SBarry Smith   PetscFunctionBegin;
7683d1c1ea0SBarry Smith   ierr                       = MatFactorInfoInitialize(&fact->info);CHKERRQ(ierr);
7693d1c1ea0SBarry Smith   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
7703d1c1ea0SBarry Smith   fact->info.shiftamount     = 100.0*PETSC_MACHINE_EPSILON;
7713d1c1ea0SBarry Smith   fact->info.zeropivot       = 100.0*PETSC_MACHINE_EPSILON;
7723d1c1ea0SBarry Smith   fact->info.pivotinblocks   = 1.0;
7733d1c1ea0SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
7743d1c1ea0SBarry Smith 
7753d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
7763d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);CHKERRQ(ierr);
7773d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
7783d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);CHKERRQ(ierr);
7793d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
7803d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);CHKERRQ(ierr);
7813d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
7823d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
7833d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
7843d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
7853d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
7863d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
7873d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
7883d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
7893d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);CHKERRQ(ierr);
7903d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
7913d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor);CHKERRQ(ierr);
7923d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor);CHKERRQ(ierr);
7933d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor);CHKERRQ(ierr);
7943d1c1ea0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor);CHKERRQ(ierr);
7953d1c1ea0SBarry Smith   PetscFunctionReturn(0);
7963d1c1ea0SBarry Smith }
797