xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 04c3f3b8f563eff258c1de90d4e02d41e6027585)
15e8efad8SHong Zhang 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
3978beb7fSPierre Jolivet #include <petsc/private/matimpl.h>
45e8efad8SHong Zhang 
54ac6704cSBarry Smith /*
64ac6704cSBarry Smith     If an ordering is not yet set and the matrix is available determine a default ordering
74ac6704cSBarry Smith */
8d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
9d71ae5a4SJacob Faibussowitsch {
107def7855SStefano Zampini   PetscBool   foundmtype, flg, destroy = PETSC_FALSE;
1126cc229bSBarry Smith   const char *prefix;
124ac6704cSBarry Smith 
134ac6704cSBarry Smith   PetscFunctionBegin;
144ac6704cSBarry Smith   if (pc->pmat) {
1526cc229bSBarry Smith     PetscCall(PCGetOptionsPrefix(pc, &prefix));
1626cc229bSBarry Smith     PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
174ac6704cSBarry Smith     PC_Factor *fact = (PC_Factor *)pc->data;
189566063dSJacob Faibussowitsch     PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL));
19978beb7fSPierre Jolivet     if (foundmtype) {
207def7855SStefano Zampini       if (!fact->fact) {
217def7855SStefano Zampini         /* factored matrix is not present at this point, we want to create it during PCSetUp.
227def7855SStefano Zampini            Since this can be called from setfromoptions, we destroy it when we are done with it */
237def7855SStefano Zampini         PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact));
247def7855SStefano Zampini         destroy = PETSC_TRUE;
257def7855SStefano Zampini       }
2603e5aca4SStefano Zampini       if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS);
2703e5aca4SStefano Zampini       if (!fact->fact->assembled) {
289566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg));
29978beb7fSPierre Jolivet         if (!flg) {
3026cc229bSBarry Smith           Mat B;
319566063dSJacob Faibussowitsch           PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B));
329566063dSJacob Faibussowitsch           PetscCall(MatHeaderReplace(fact->fact, &B));
33978beb7fSPierre Jolivet         }
344ac6704cSBarry Smith       }
354ac6704cSBarry Smith       if (!fact->ordering) {
364ac6704cSBarry Smith         PetscBool       canuseordering;
374ac6704cSBarry Smith         MatOrderingType otype;
384ac6704cSBarry Smith 
399566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering));
404ac6704cSBarry Smith         if (canuseordering) {
419566063dSJacob Faibussowitsch           PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype));
424ac6704cSBarry Smith         } else otype = MATORDERINGEXTERNAL;
439566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering));
444ac6704cSBarry Smith       }
457def7855SStefano Zampini       if (destroy) PetscCall(MatDestroy(&fact->fact));
464ac6704cSBarry Smith     }
47978beb7fSPierre Jolivet   }
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
494ac6704cSBarry Smith }
504ac6704cSBarry Smith 
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag)
52d71ae5a4SJacob Faibussowitsch {
533d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
543d1c1ea0SBarry Smith 
553d1c1ea0SBarry Smith   PetscFunctionBegin;
563d1c1ea0SBarry Smith   lu->reuseordering = flag;
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
583d1c1ea0SBarry Smith }
593d1c1ea0SBarry Smith 
60d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag)
61d71ae5a4SJacob Faibussowitsch {
623d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
633d1c1ea0SBarry Smith 
643d1c1ea0SBarry Smith   PetscFunctionBegin;
653d1c1ea0SBarry Smith   lu->reusefill = flag;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
673d1c1ea0SBarry Smith }
683d1c1ea0SBarry Smith 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg)
70d71ae5a4SJacob Faibussowitsch {
713d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
723d1c1ea0SBarry Smith 
733d1c1ea0SBarry Smith   PetscFunctionBegin;
743d1c1ea0SBarry Smith   dir->inplace = flg;
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
763d1c1ea0SBarry Smith }
773d1c1ea0SBarry Smith 
78d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg)
79d71ae5a4SJacob Faibussowitsch {
803d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
813d1c1ea0SBarry Smith 
823d1c1ea0SBarry Smith   PetscFunctionBegin;
833d1c1ea0SBarry Smith   *flg = dir->inplace;
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
853d1c1ea0SBarry Smith }
863d1c1ea0SBarry Smith 
87f8260c8fSBarry Smith /*@
88f1580f4eSBarry Smith   PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may
89f8260c8fSBarry Smith   set the options for that particular factorization object.
90f8260c8fSBarry Smith 
91f8260c8fSBarry Smith   Input Parameter:
92f8260c8fSBarry Smith . pc - the preconditioner context
93f8260c8fSBarry Smith 
94f1580f4eSBarry Smith   Note:
95f1580f4eSBarry Smith   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.
9603e5aca4SStefano Zampini   This function raises an error if the requested combination of solver package and matrix type is not supported.
97f8260c8fSBarry Smith 
982bd2b0e6SSatish Balay   Level: intermediate
992bd2b0e6SSatish Balay 
100f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
101f8260c8fSBarry Smith @*/
102d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
103d71ae5a4SJacob Faibussowitsch {
104f8260c8fSBarry Smith   PetscFunctionBegin;
105f8260c8fSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
106cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108f8260c8fSBarry Smith }
109f8260c8fSBarry Smith 
110ee45ca4aSHong Zhang /*@
111ee45ca4aSHong Zhang   PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
112ee45ca4aSHong Zhang 
113c3339decSBarry Smith   Logically Collective
114ee45ca4aSHong Zhang 
115ee45ca4aSHong Zhang   Input Parameters:
116afaefe49SHong Zhang + pc   - the preconditioner context
117afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero
118ee45ca4aSHong Zhang 
119ee45ca4aSHong Zhang   Options Database Key:
120ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
121ee45ca4aSHong Zhang 
122ee45ca4aSHong Zhang   Level: intermediate
123ee45ca4aSHong Zhang 
124f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
125ee45ca4aSHong Zhang @*/
126d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero)
127d71ae5a4SJacob Faibussowitsch {
128ee45ca4aSHong Zhang   PetscFunctionBegin;
1290700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
130c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, zero, 2);
131cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero));
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133ee45ca4aSHong Zhang }
134ee45ca4aSHong Zhang 
135915743fcSHong Zhang /*@
136915743fcSHong Zhang   PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
137915743fcSHong Zhang   numerical factorization, thus the matrix has nonzero pivots
138915743fcSHong Zhang 
139c3339decSBarry Smith   Logically Collective
140915743fcSHong Zhang 
141915743fcSHong Zhang   Input Parameters:
142915743fcSHong Zhang + pc        - the preconditioner context
143f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS`
144915743fcSHong Zhang 
145915743fcSHong Zhang   Options Database Key:
14628d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
147915743fcSHong Zhang 
148915743fcSHong Zhang   Level: intermediate
149915743fcSHong Zhang 
150f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
151915743fcSHong Zhang @*/
152d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype)
153d71ae5a4SJacob Faibussowitsch {
154d90ac83dSHong Zhang   PetscFunctionBegin;
1550700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
156c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, shifttype, 2);
157cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype));
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159d90ac83dSHong Zhang }
160d90ac83dSHong Zhang 
161915743fcSHong Zhang /*@
162915743fcSHong Zhang   PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
163915743fcSHong Zhang   numerical factorization, thus the matrix has nonzero pivots
164915743fcSHong Zhang 
165c3339decSBarry Smith   Logically Collective
166915743fcSHong Zhang 
167915743fcSHong Zhang   Input Parameters:
168915743fcSHong Zhang + pc          - the preconditioner context
169f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default
170915743fcSHong Zhang 
171915743fcSHong Zhang   Options Database Key:
172f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
173915743fcSHong Zhang 
174915743fcSHong Zhang   Level: intermediate
175915743fcSHong Zhang 
176f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, ``PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
177915743fcSHong Zhang @*/
178d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount)
179d71ae5a4SJacob Faibussowitsch {
180d90ac83dSHong Zhang   PetscFunctionBegin;
1810700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
182c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, shiftamount, 2);
183cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount));
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
185d90ac83dSHong Zhang }
186d90ac83dSHong Zhang 
1876d33885cSprj- /*@
188f1580f4eSBarry Smith   PCFactorSetDropTolerance - The preconditioner will use an `PCILU`
189f1580f4eSBarry Smith   based on a drop tolerance.
19085317021SBarry Smith 
191c3339decSBarry Smith   Logically Collective
19285317021SBarry Smith 
19385317021SBarry Smith   Input Parameters:
19485317021SBarry Smith + pc          - the preconditioner context
19585317021SBarry Smith . dt          - the drop tolerance, try from 1.e-10 to .1
19685317021SBarry Smith . dtcol       - tolerance for column pivot, good values [0.1 to 0.01]
19785317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value
19885317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
19985317021SBarry Smith 
20085317021SBarry Smith   Options Database Key:
201b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
20285317021SBarry Smith 
20385317021SBarry Smith   Level: intermediate
20485317021SBarry Smith 
205f1580f4eSBarry Smith   Note:
20685317021SBarry Smith   There are NO default values for the 3 parameters, you must set them with reasonable values for your
20785317021SBarry Smith   matrix. We don't know how to compute reasonable values.
20885317021SBarry Smith 
209f1580f4eSBarry Smith .seealso: `PCILU`
2106d33885cSprj- @*/
211d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount)
212d71ae5a4SJacob Faibussowitsch {
21385317021SBarry Smith   PetscFunctionBegin;
2140700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
215064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveReal(pc, dtcol, 3);
216064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(pc, maxrowcount, 4);
217cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount));
2183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21985317021SBarry Smith }
22085317021SBarry Smith 
221c7f610a1SBarry Smith /*@
222c7f610a1SBarry Smith   PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
223c7f610a1SBarry Smith 
224c7f610a1SBarry Smith   Not Collective
225c7f610a1SBarry Smith 
2262fe279fdSBarry Smith   Input Parameter:
227c7f610a1SBarry Smith . pc - the preconditioner context
228c7f610a1SBarry Smith 
229c7f610a1SBarry Smith   Output Parameter:
230c7f610a1SBarry Smith . pivot - the tolerance
231c7f610a1SBarry Smith 
232c7f610a1SBarry Smith   Level: intermediate
233c7f610a1SBarry Smith 
234f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()`
235c7f610a1SBarry Smith @*/
236d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot)
237d71ae5a4SJacob Faibussowitsch {
238c7f610a1SBarry Smith   PetscFunctionBegin;
239c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
240cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot));
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
242c7f610a1SBarry Smith }
243c7f610a1SBarry Smith 
244c7f610a1SBarry Smith /*@
245c7f610a1SBarry Smith   PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
246c7f610a1SBarry Smith 
247c7f610a1SBarry Smith   Not Collective
248c7f610a1SBarry Smith 
2492fe279fdSBarry Smith   Input Parameter:
250c7f610a1SBarry Smith . pc - the preconditioner context
251c7f610a1SBarry Smith 
252c7f610a1SBarry Smith   Output Parameter:
253c7f610a1SBarry Smith . shift - how much to shift the diagonal entry
254c7f610a1SBarry Smith 
255c7f610a1SBarry Smith   Level: intermediate
256c7f610a1SBarry Smith 
257f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
258c7f610a1SBarry Smith @*/
259d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift)
260d71ae5a4SJacob Faibussowitsch {
261c7f610a1SBarry Smith   PetscFunctionBegin;
262c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
263cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift));
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
265c7f610a1SBarry Smith }
266c7f610a1SBarry Smith 
267c7f610a1SBarry Smith /*@
268c7f610a1SBarry Smith   PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
269c7f610a1SBarry Smith 
270c7f610a1SBarry Smith   Not Collective
271c7f610a1SBarry Smith 
272f1580f4eSBarry Smith   Input Parameter:
273c7f610a1SBarry Smith . pc - the preconditioner context
274c7f610a1SBarry Smith 
275c7f610a1SBarry Smith   Output Parameter:
276f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS`
277c7f610a1SBarry Smith 
278c7f610a1SBarry Smith   Level: intermediate
279c7f610a1SBarry Smith 
280f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
281c7f610a1SBarry Smith @*/
282d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type)
283d71ae5a4SJacob Faibussowitsch {
284c7f610a1SBarry Smith   PetscFunctionBegin;
285c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
286cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type));
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288c7f610a1SBarry Smith }
289c7f610a1SBarry Smith 
2902591b318SToby Isaac /*@
2912591b318SToby Isaac   PCFactorGetLevels - Gets the number of levels of fill to use.
2922591b318SToby Isaac 
293c3339decSBarry Smith   Logically Collective
2942591b318SToby Isaac 
2952fe279fdSBarry Smith   Input Parameter:
2962591b318SToby Isaac . pc - the preconditioner context
2972591b318SToby Isaac 
2982591b318SToby Isaac   Output Parameter:
2992591b318SToby Isaac . levels - number of levels of fill
3002591b318SToby Isaac 
3012591b318SToby Isaac   Level: intermediate
3022591b318SToby Isaac 
303f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetLevels()`
3042591b318SToby Isaac @*/
305d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels)
306d71ae5a4SJacob Faibussowitsch {
3072591b318SToby Isaac   PetscFunctionBegin;
3082591b318SToby Isaac   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
309cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3112591b318SToby Isaac }
3122591b318SToby Isaac 
31385317021SBarry Smith /*@
31485317021SBarry Smith   PCFactorSetLevels - Sets the number of levels of fill to use.
31585317021SBarry Smith 
316c3339decSBarry Smith   Logically Collective
31785317021SBarry Smith 
31885317021SBarry Smith   Input Parameters:
31985317021SBarry Smith + pc     - the preconditioner context
32085317021SBarry Smith - levels - number of levels of fill
32185317021SBarry Smith 
32285317021SBarry Smith   Options Database Key:
32385317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level
32485317021SBarry Smith 
32585317021SBarry Smith   Level: intermediate
32685317021SBarry Smith 
327f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetLevels()`
32885317021SBarry Smith @*/
329d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels)
330d71ae5a4SJacob Faibussowitsch {
33185317021SBarry Smith   PetscFunctionBegin;
3320700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3335f80ce2aSJacob Faibussowitsch   PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels");
334c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, levels, 2);
335cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels));
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33785317021SBarry Smith }
33885317021SBarry Smith 
33985317021SBarry Smith /*@
34085317021SBarry Smith   PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
34185317021SBarry Smith   treated as level 0 fill even if there is no non-zero location.
34285317021SBarry Smith 
343c3339decSBarry Smith   Logically Collective
34485317021SBarry Smith 
34585317021SBarry Smith   Input Parameters:
34685317021SBarry Smith + pc  - the preconditioner context
347f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
34885317021SBarry Smith 
34985317021SBarry Smith   Options Database Key:
35067b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill
35185317021SBarry Smith 
352f1580f4eSBarry Smith   Note:
35385317021SBarry Smith   Does not apply with 0 fill.
35485317021SBarry Smith 
35585317021SBarry Smith   Level: intermediate
35685317021SBarry Smith 
357f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()`
35885317021SBarry Smith @*/
359d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg)
360d71ae5a4SJacob Faibussowitsch {
36185317021SBarry Smith   PetscFunctionBegin;
3620700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
363cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg));
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36592e9c092SBarry Smith }
36692e9c092SBarry Smith 
36792e9c092SBarry Smith /*@
36892e9c092SBarry Smith   PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
36992e9c092SBarry Smith   treated as level 0 fill even if there is no non-zero location.
37092e9c092SBarry Smith 
371c3339decSBarry Smith   Logically Collective
37292e9c092SBarry Smith 
37392e9c092SBarry Smith   Input Parameter:
37492e9c092SBarry Smith . pc - the preconditioner context
37592e9c092SBarry Smith 
37692e9c092SBarry Smith   Output Parameter:
377f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
37892e9c092SBarry Smith 
379f1580f4eSBarry Smith   Note:
38092e9c092SBarry Smith   Does not apply with 0 fill.
38192e9c092SBarry Smith 
38292e9c092SBarry Smith   Level: intermediate
38392e9c092SBarry Smith 
384f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()`
38592e9c092SBarry Smith @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg)
387d71ae5a4SJacob Faibussowitsch {
38892e9c092SBarry Smith   PetscFunctionBegin;
38992e9c092SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
390cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg));
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39285317021SBarry Smith }
39385317021SBarry Smith 
39485317021SBarry Smith /*@
39585317021SBarry Smith   PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
39685317021SBarry Smith 
397c3339decSBarry Smith   Logically Collective
39885317021SBarry Smith 
39985317021SBarry Smith   Input Parameters:
40085317021SBarry Smith + pc   - the preconditioner context
401feefa0e1SJacob Faibussowitsch - rtol - diagonal entries smaller than this in absolute value are considered zero
40285317021SBarry Smith 
40385317021SBarry Smith   Options Database Key:
404147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance
40585317021SBarry Smith 
40685317021SBarry Smith   Level: intermediate
40785317021SBarry Smith 
408f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
40985317021SBarry Smith @*/
410d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol)
411d71ae5a4SJacob Faibussowitsch {
41285317021SBarry Smith   PetscFunctionBegin;
4130700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
414c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, rtol, 2);
415cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol));
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41785317021SBarry Smith }
41885317021SBarry Smith 
419bf6011e8SBarry Smith /*@C
420f1580f4eSBarry Smith   PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization
42185317021SBarry Smith 
422c3339decSBarry Smith   Logically Collective
42385317021SBarry Smith 
42485317021SBarry Smith   Input Parameters:
42585317021SBarry Smith + pc    - the preconditioner context
426f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
42785317021SBarry Smith 
42885317021SBarry Smith   Options Database Key:
4293ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
43085317021SBarry Smith 
43185317021SBarry Smith   Level: intermediate
43285317021SBarry Smith 
43385317021SBarry Smith   Note:
43485317021SBarry Smith   By default this will use the PETSc factorization if it exists
43585317021SBarry Smith 
436feefa0e1SJacob Faibussowitsch .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`,
437f1580f4eSBarry Smith           `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
43885317021SBarry Smith @*/
439d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype)
440d71ae5a4SJacob Faibussowitsch {
44185317021SBarry Smith   PetscFunctionBegin;
4420700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
443cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype));
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44585317021SBarry Smith }
44685317021SBarry Smith 
447bf6011e8SBarry Smith /*@C
448f1580f4eSBarry Smith   PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization
4497112b564SBarry Smith 
450c5eb9154SBarry Smith   Not Collective
4517112b564SBarry Smith 
4527112b564SBarry Smith   Input Parameter:
4537112b564SBarry Smith . pc - the preconditioner context
4547112b564SBarry Smith 
4557112b564SBarry Smith   Output Parameter:
456f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
4577112b564SBarry Smith 
4587112b564SBarry Smith   Level: intermediate
4597112b564SBarry Smith 
46042747ad1SJacob Faibussowitsch .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `MATSOLVERSUPERLU`,
46142747ad1SJacob Faibussowitsch `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
4627112b564SBarry Smith @*/
463d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype)
464d71ae5a4SJacob Faibussowitsch {
4655f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(PC, MatSolverType *);
4667112b564SBarry Smith 
4677112b564SBarry Smith   PetscFunctionBegin;
4680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4694f572ea9SToby Isaac   PetscAssertPointer(stype, 2);
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f));
4719566063dSJacob Faibussowitsch   if (f) PetscCall((*f)(pc, stype));
4725f80ce2aSJacob Faibussowitsch   else *stype = NULL;
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4747112b564SBarry Smith }
4757112b564SBarry Smith 
47685317021SBarry Smith /*@
47785317021SBarry Smith   PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
47885317021SBarry Smith   fill = number nonzeros in factor/number nonzeros in original matrix.
47985317021SBarry Smith 
480c5eb9154SBarry Smith   Not Collective, each process can expect a different amount of fill
48185317021SBarry Smith 
48285317021SBarry Smith   Input Parameters:
48385317021SBarry Smith + pc   - the preconditioner context
48485317021SBarry Smith - fill - amount of expected fill
48585317021SBarry Smith 
48685317021SBarry Smith   Options Database Key:
48785317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount
48885317021SBarry Smith 
48985317021SBarry Smith   Level: intermediate
49085317021SBarry Smith 
491f1580f4eSBarry Smith   Notes:
49285317021SBarry Smith   For sparse matrix factorizations it is difficult to predict how much
49385317021SBarry Smith   fill to expect. By running with the option -info PETSc will print the
49485317021SBarry Smith   actual amount of fill used; allowing you to set the value accurately for
49585317021SBarry Smith   future runs. Default PETSc uses a value of 5.0
49685317021SBarry Smith 
497f1580f4eSBarry Smith   This is ignored for most solver packages
49801a79839SBarry Smith 
499f1580f4eSBarry Smith   This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels.
500f1580f4eSBarry Smith 
501f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()`
50285317021SBarry Smith @*/
503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill)
504d71ae5a4SJacob Faibussowitsch {
50585317021SBarry Smith   PetscFunctionBegin;
5060700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5075f80ce2aSJacob Faibussowitsch   PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0");
508cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51085317021SBarry Smith }
51185317021SBarry Smith 
51285317021SBarry Smith /*@
513*04c3f3b8SBarry Smith   PCFactorSetUseInPlace - Tells the preconditioner to do an in-place factorization.
51485317021SBarry Smith 
515c3339decSBarry Smith   Logically Collective
51685317021SBarry Smith 
51785317021SBarry Smith   Input Parameters:
5188e37d05fSBarry Smith + pc  - the preconditioner context
519f1580f4eSBarry Smith - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
52085317021SBarry Smith 
52185317021SBarry Smith   Options Database Key:
5228e37d05fSBarry Smith . -pc_factor_in_place <true,false> - Activate/deactivate in-place factorization
52385317021SBarry Smith 
524f1580f4eSBarry Smith   Note:
525*04c3f3b8SBarry Smith   For dense matrices, this enables the solution of much larger problems.
526*04c3f3b8SBarry Smith   For sparse matrices the factorization cannot be done truly in-place
527*04c3f3b8SBarry Smith   so this does not save memory during the factorization, but after the matrix
528*04c3f3b8SBarry Smith   is factored, the original unfactored matrix is freed, thus recovering that
529*04c3f3b8SBarry Smith   space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
530*04c3f3b8SBarry Smith 
531f1580f4eSBarry Smith   `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when
53285317021SBarry Smith   a different matrix is provided for the multiply and the preconditioner in
533f1580f4eSBarry Smith   a call to `KSPSetOperators()`.
53485317021SBarry Smith   This is because the Krylov space methods require an application of the
53585317021SBarry Smith   matrix multiplication, which is not possible here because the matrix has
53685317021SBarry Smith   been factored in-place, replacing the original matrix.
53785317021SBarry Smith 
53885317021SBarry Smith   Level: intermediate
53985317021SBarry Smith 
540*04c3f3b8SBarry Smith .seealso: `PC`, `Mat`, `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()`
54185317021SBarry Smith @*/
542d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg)
543d71ae5a4SJacob Faibussowitsch {
54485317021SBarry Smith   PetscFunctionBegin;
5450700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
546cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg));
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5488e37d05fSBarry Smith }
5498e37d05fSBarry Smith 
5508e37d05fSBarry Smith /*@
5518e37d05fSBarry Smith   PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
5528e37d05fSBarry Smith 
553c3339decSBarry Smith   Logically Collective
5548e37d05fSBarry Smith 
5558e37d05fSBarry Smith   Input Parameter:
5568e37d05fSBarry Smith . pc - the preconditioner context
5578e37d05fSBarry Smith 
5588e37d05fSBarry Smith   Output Parameter:
559f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
5608e37d05fSBarry Smith 
5618e37d05fSBarry Smith   Level: intermediate
5628e37d05fSBarry Smith 
563f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()`
5648e37d05fSBarry Smith @*/
565d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg)
566d71ae5a4SJacob Faibussowitsch {
5678e37d05fSBarry Smith   PetscFunctionBegin;
5688e37d05fSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
569cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg));
5703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57185317021SBarry Smith }
57285317021SBarry Smith 
57385317021SBarry Smith /*@C
57485317021SBarry Smith   PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
575f1580f4eSBarry Smith   be used in the `PCLU`, `PCCHOLESKY`, `PCILU`,  or `PCICC` preconditioners
57685317021SBarry Smith 
577c3339decSBarry Smith   Logically Collective
57885317021SBarry Smith 
57985317021SBarry Smith   Input Parameters:
58085317021SBarry Smith + pc       - the preconditioner context
581f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM`
58285317021SBarry Smith 
58385317021SBarry Smith   Options Database Key:
5842c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine
58585317021SBarry Smith 
58685317021SBarry Smith   Level: intermediate
58785317021SBarry Smith 
58895452b02SPatrick Sanan   Notes:
5894ac6704cSBarry Smith   Nested dissection is used by default for some of PETSc's sparse matrix formats
59085317021SBarry Smith 
591f1580f4eSBarry Smith   For `PCCHOLESKY` and `PCICC` and the `MATSBAIJ` format the only reordering available is natural since only the upper half of the matrix is stored
5929bd791bbSBarry Smith   and reordering this matrix is very expensive.
5939bd791bbSBarry Smith 
594f1580f4eSBarry Smith   You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering.
5959bd791bbSBarry Smith 
596f1580f4eSBarry Smith   `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others.
5972c7c0729SBarry Smith 
598f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM`
59985317021SBarry Smith @*/
600d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering)
601d71ae5a4SJacob Faibussowitsch {
60285317021SBarry Smith   PetscFunctionBegin;
603c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
604cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering));
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60685317021SBarry Smith }
60785317021SBarry Smith 
60885317021SBarry Smith /*@
6098ff23777SHong Zhang   PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
61085317021SBarry Smith   For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
611f1580f4eSBarry Smith   it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used.
61285317021SBarry Smith 
613c3339decSBarry Smith   Logically Collective
61485317021SBarry Smith 
61585317021SBarry Smith   Input Parameters:
61685317021SBarry Smith + pc    - the preconditioner context
61785317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
61885317021SBarry Smith 
61985317021SBarry Smith   Options Database Key:
620147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance
62185317021SBarry Smith 
62285317021SBarry Smith   Level: intermediate
62385317021SBarry Smith 
624f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
62585317021SBarry Smith @*/
626d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol)
627d71ae5a4SJacob Faibussowitsch {
62885317021SBarry Smith   PetscFunctionBegin;
629c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
630c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, dtcol, 2);
631cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol));
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63385317021SBarry Smith }
63485317021SBarry Smith 
63585317021SBarry Smith /*@
63685317021SBarry Smith   PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
637f1580f4eSBarry Smith   with `MATBAIJ` or `MATSBAIJ` matrices
63885317021SBarry Smith 
639c3339decSBarry Smith   Logically Collective
64085317021SBarry Smith 
64185317021SBarry Smith   Input Parameters:
64285317021SBarry Smith + pc    - the preconditioner context
643f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE`
64485317021SBarry Smith 
64585317021SBarry Smith   Options Database Key:
646f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ`
64785317021SBarry Smith 
64885317021SBarry Smith   Level: intermediate
64985317021SBarry Smith 
650f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
65185317021SBarry Smith @*/
652d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot)
653d71ae5a4SJacob Faibussowitsch {
65485317021SBarry Smith   PetscFunctionBegin;
655c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
656acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, pivot, 2);
657cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot));
6583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65985317021SBarry Smith }
66085317021SBarry Smith 
66185317021SBarry Smith /*@
662288e7d53SBarry Smith   PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
66385317021SBarry Smith   this causes later ones to use the fill ratio computed in the initial factorization.
66485317021SBarry Smith 
665c3339decSBarry Smith   Logically Collective
66685317021SBarry Smith 
66785317021SBarry Smith   Input Parameters:
66885317021SBarry Smith + pc   - the preconditioner context
669f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
67085317021SBarry Smith 
67185317021SBarry Smith   Options Database Key:
672f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
67385317021SBarry Smith 
67485317021SBarry Smith   Level: intermediate
67585317021SBarry Smith 
676f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()`
67785317021SBarry Smith @*/
678d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag)
679d71ae5a4SJacob Faibussowitsch {
68085317021SBarry Smith   PetscFunctionBegin;
681064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
682acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, flag, 2);
683cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag));
6843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68585317021SBarry Smith }
6863d1c1ea0SBarry Smith 
687d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype)
688d71ae5a4SJacob Faibussowitsch {
6893d1c1ea0SBarry Smith   PC_Factor *fact = (PC_Factor *)pc->data;
6903d1c1ea0SBarry Smith 
6913d1c1ea0SBarry Smith   PetscFunctionBegin;
6929566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&fact->info));
6934ac6704cSBarry Smith   fact->factortype           = ftype;
6943d1c1ea0SBarry Smith   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
6953d1c1ea0SBarry Smith   fact->info.shiftamount     = 100.0 * PETSC_MACHINE_EPSILON;
6963d1c1ea0SBarry Smith   fact->info.zeropivot       = 100.0 * PETSC_MACHINE_EPSILON;
6973d1c1ea0SBarry Smith   fact->info.pivotinblocks   = 1.0;
6983d1c1ea0SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
6993d1c1ea0SBarry Smith 
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor));
7019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor));
7029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor));
7039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor));
7049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor));
7059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor));
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor));
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor));
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor));
7099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor));
7109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor));
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor));
7129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor));
7139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor));
7149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor));
7159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor));
7169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor));
7179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor));
7189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor));
7199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor));
7203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7213d1c1ea0SBarry Smith }
7222e956fe4SStefano Zampini 
723d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc)
724d71ae5a4SJacob Faibussowitsch {
7252e956fe4SStefano Zampini   PetscFunctionBegin;
7262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL));
7272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL));
7282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL));
7292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL));
7302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL));
7312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL));
7322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL));
7332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL));
7342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL));
7352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL));
7362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL));
7372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL));
7382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL));
7392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL));
7402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL));
7412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL));
7422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL));
7432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL));
7442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL));
7452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL));
7462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL));
7472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL));
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7492e956fe4SStefano Zampini }
750