xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 901f93825bbaa60582af604c6700caf57884a2e1)
1c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
2978beb7fSPierre Jolivet #include <petsc/private/matimpl.h>
35e8efad8SHong Zhang 
44ac6704cSBarry Smith /*
54ac6704cSBarry Smith     If an ordering is not yet set and the matrix is available determine a default ordering
64ac6704cSBarry Smith */
PCFactorSetDefaultOrdering_Factor(PC pc)7d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
8d71ae5a4SJacob Faibussowitsch {
97def7855SStefano Zampini   PetscBool   foundmtype, flg, destroy = PETSC_FALSE;
1026cc229bSBarry Smith   const char *prefix;
114ac6704cSBarry Smith 
124ac6704cSBarry Smith   PetscFunctionBegin;
134ac6704cSBarry Smith   if (pc->pmat) {
1426cc229bSBarry Smith     PetscCall(PCGetOptionsPrefix(pc, &prefix));
1526cc229bSBarry Smith     PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
164ac6704cSBarry Smith     PC_Factor *fact = (PC_Factor *)pc->data;
179566063dSJacob Faibussowitsch     PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL));
18978beb7fSPierre Jolivet     if (foundmtype) {
197def7855SStefano Zampini       if (!fact->fact) {
207def7855SStefano Zampini         /* factored matrix is not present at this point, we want to create it during PCSetUp.
217def7855SStefano Zampini            Since this can be called from setfromoptions, we destroy it when we are done with it */
227def7855SStefano Zampini         PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact));
237def7855SStefano Zampini         destroy = PETSC_TRUE;
247def7855SStefano Zampini       }
2503e5aca4SStefano Zampini       if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS);
2603e5aca4SStefano Zampini       if (!fact->fact->assembled) {
279566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg));
28978beb7fSPierre Jolivet         if (!flg) {
2926cc229bSBarry Smith           Mat B;
309566063dSJacob Faibussowitsch           PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B));
319566063dSJacob Faibussowitsch           PetscCall(MatHeaderReplace(fact->fact, &B));
32978beb7fSPierre Jolivet         }
334ac6704cSBarry Smith       }
344ac6704cSBarry Smith       if (!fact->ordering) {
354ac6704cSBarry Smith         PetscBool       canuseordering;
364ac6704cSBarry Smith         MatOrderingType otype;
374ac6704cSBarry Smith 
389566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering));
39ac530a7eSPierre Jolivet         if (canuseordering) PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype));
40ac530a7eSPierre Jolivet         else otype = MATORDERINGEXTERNAL;
419566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering));
424ac6704cSBarry Smith       }
437def7855SStefano Zampini       if (destroy) PetscCall(MatDestroy(&fact->fact));
444ac6704cSBarry Smith     }
45978beb7fSPierre Jolivet   }
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
474ac6704cSBarry Smith }
484ac6704cSBarry Smith 
PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag)49d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag)
50d71ae5a4SJacob Faibussowitsch {
513d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
523d1c1ea0SBarry Smith 
533d1c1ea0SBarry Smith   PetscFunctionBegin;
543d1c1ea0SBarry Smith   lu->reuseordering = flag;
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563d1c1ea0SBarry Smith }
573d1c1ea0SBarry Smith 
PCFactorSetReuseFill_Factor(PC pc,PetscBool flag)58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag)
59d71ae5a4SJacob Faibussowitsch {
603d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
613d1c1ea0SBarry Smith 
623d1c1ea0SBarry Smith   PetscFunctionBegin;
633d1c1ea0SBarry Smith   lu->reusefill = flag;
643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653d1c1ea0SBarry Smith }
663d1c1ea0SBarry Smith 
PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg)67d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg)
68d71ae5a4SJacob Faibussowitsch {
693d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
703d1c1ea0SBarry Smith 
713d1c1ea0SBarry Smith   PetscFunctionBegin;
723d1c1ea0SBarry Smith   dir->inplace = flg;
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
743d1c1ea0SBarry Smith }
753d1c1ea0SBarry Smith 
PCFactorGetUseInPlace_Factor(PC pc,PetscBool * flg)76d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg)
77d71ae5a4SJacob Faibussowitsch {
783d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
793d1c1ea0SBarry Smith 
803d1c1ea0SBarry Smith   PetscFunctionBegin;
813d1c1ea0SBarry Smith   *flg = dir->inplace;
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
833d1c1ea0SBarry Smith }
843d1c1ea0SBarry Smith 
85f8260c8fSBarry Smith /*@
86f1580f4eSBarry Smith   PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may
87f8260c8fSBarry Smith   set the options for that particular factorization object.
88f8260c8fSBarry Smith 
89f8260c8fSBarry Smith   Input Parameter:
90f8260c8fSBarry Smith . pc - the preconditioner context
91f8260c8fSBarry Smith 
92*da34b7cdSBarry Smith   Level: intermediate
93*da34b7cdSBarry 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 
98*da34b7cdSBarry Smith   Developer Note:
99*da34b7cdSBarry Smith   This function should have a name that clearly indicates that this calls `MatGetFactor()` and thus populates the `->factor` field. The `MatSetSolverType` portion of the name
100*da34b7cdSBarry Smith   may not add value to the clarity of the purpose of the function and could be removed.
1012bd2b0e6SSatish Balay 
102562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
103f8260c8fSBarry Smith @*/
PCFactorSetUpMatSolverType(PC pc)104d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
105d71ae5a4SJacob Faibussowitsch {
106f8260c8fSBarry Smith   PetscFunctionBegin;
107f8260c8fSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
108cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc));
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110f8260c8fSBarry Smith }
111f8260c8fSBarry Smith 
112ee45ca4aSHong Zhang /*@
113ee45ca4aSHong Zhang   PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
114ee45ca4aSHong Zhang 
115c3339decSBarry Smith   Logically Collective
116ee45ca4aSHong Zhang 
117ee45ca4aSHong Zhang   Input Parameters:
118afaefe49SHong Zhang + pc   - the preconditioner context
119afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero
120ee45ca4aSHong Zhang 
121ee45ca4aSHong Zhang   Options Database Key:
122ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
123ee45ca4aSHong Zhang 
124ee45ca4aSHong Zhang   Level: intermediate
125ee45ca4aSHong Zhang 
126562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
127ee45ca4aSHong Zhang @*/
PCFactorSetZeroPivot(PC pc,PetscReal zero)128d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero)
129d71ae5a4SJacob Faibussowitsch {
130ee45ca4aSHong Zhang   PetscFunctionBegin;
1310700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
132c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, zero, 2);
133cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135ee45ca4aSHong Zhang }
136ee45ca4aSHong Zhang 
137915743fcSHong Zhang /*@
138915743fcSHong Zhang   PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
139915743fcSHong Zhang   numerical factorization, thus the matrix has nonzero pivots
140915743fcSHong Zhang 
141c3339decSBarry Smith   Logically Collective
142915743fcSHong Zhang 
143915743fcSHong Zhang   Input Parameters:
144915743fcSHong Zhang + pc        - the preconditioner context
145f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS`
146915743fcSHong Zhang 
147915743fcSHong Zhang   Options Database Key:
14828d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
149915743fcSHong Zhang 
150915743fcSHong Zhang   Level: intermediate
151915743fcSHong Zhang 
152562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
153915743fcSHong Zhang @*/
PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)154d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype)
155d71ae5a4SJacob Faibussowitsch {
156d90ac83dSHong Zhang   PetscFunctionBegin;
1570700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
158c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, shifttype, 2);
159cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype));
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161d90ac83dSHong Zhang }
162d90ac83dSHong Zhang 
163915743fcSHong Zhang /*@
164915743fcSHong Zhang   PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
165915743fcSHong Zhang   numerical factorization, thus the matrix has nonzero pivots
166915743fcSHong Zhang 
167c3339decSBarry Smith   Logically Collective
168915743fcSHong Zhang 
169915743fcSHong Zhang   Input Parameters:
170915743fcSHong Zhang + pc          - the preconditioner context
171f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default
172915743fcSHong Zhang 
173915743fcSHong Zhang   Options Database Key:
174f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
175915743fcSHong Zhang 
176915743fcSHong Zhang   Level: intermediate
177915743fcSHong Zhang 
178a94f484eSPierre Jolivet .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
179915743fcSHong Zhang @*/
PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)180d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount)
181d71ae5a4SJacob Faibussowitsch {
182d90ac83dSHong Zhang   PetscFunctionBegin;
1830700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
184c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, shiftamount, 2);
185cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount));
1863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187d90ac83dSHong Zhang }
188d90ac83dSHong Zhang 
1896d33885cSprj- /*@
190f1580f4eSBarry Smith   PCFactorSetDropTolerance - The preconditioner will use an `PCILU`
191f1580f4eSBarry Smith   based on a drop tolerance.
19285317021SBarry Smith 
193c3339decSBarry Smith   Logically Collective
19485317021SBarry Smith 
19585317021SBarry Smith   Input Parameters:
19685317021SBarry Smith + pc          - the preconditioner context
19785317021SBarry Smith . dt          - the drop tolerance, try from 1.e-10 to .1
19885317021SBarry Smith . dtcol       - tolerance for column pivot, good values [0.1 to 0.01]
19985317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value
20085317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
20185317021SBarry Smith 
20285317021SBarry Smith   Options Database Key:
203b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
20485317021SBarry Smith 
20585317021SBarry Smith   Level: intermediate
20685317021SBarry Smith 
207f1580f4eSBarry Smith   Note:
20885317021SBarry Smith   There are NO default values for the 3 parameters, you must set them with reasonable values for your
20985317021SBarry Smith   matrix. We don't know how to compute reasonable values.
21085317021SBarry Smith 
211562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`
2126d33885cSprj- @*/
PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)213d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount)
214d71ae5a4SJacob Faibussowitsch {
21585317021SBarry Smith   PetscFunctionBegin;
2160700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
217064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveReal(pc, dtcol, 3);
218064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(pc, maxrowcount, 4);
219cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount));
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22185317021SBarry Smith }
22285317021SBarry Smith 
223c7f610a1SBarry Smith /*@
224c7f610a1SBarry Smith   PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
225c7f610a1SBarry Smith 
226c7f610a1SBarry Smith   Not Collective
227c7f610a1SBarry Smith 
2282fe279fdSBarry Smith   Input Parameter:
229c7f610a1SBarry Smith . pc - the preconditioner context
230c7f610a1SBarry Smith 
231c7f610a1SBarry Smith   Output Parameter:
232c7f610a1SBarry Smith . pivot - the tolerance
233c7f610a1SBarry Smith 
234c7f610a1SBarry Smith   Level: intermediate
235c7f610a1SBarry Smith 
236562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()`
237c7f610a1SBarry Smith @*/
PCFactorGetZeroPivot(PC pc,PetscReal * pivot)238d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot)
239d71ae5a4SJacob Faibussowitsch {
240c7f610a1SBarry Smith   PetscFunctionBegin;
241c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
242cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot));
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244c7f610a1SBarry Smith }
245c7f610a1SBarry Smith 
246c7f610a1SBarry Smith /*@
247c7f610a1SBarry Smith   PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
248c7f610a1SBarry Smith 
249c7f610a1SBarry Smith   Not Collective
250c7f610a1SBarry Smith 
2512fe279fdSBarry Smith   Input Parameter:
252c7f610a1SBarry Smith . pc - the preconditioner context
253c7f610a1SBarry Smith 
254c7f610a1SBarry Smith   Output Parameter:
255c7f610a1SBarry Smith . shift - how much to shift the diagonal entry
256c7f610a1SBarry Smith 
257c7f610a1SBarry Smith   Level: intermediate
258c7f610a1SBarry Smith 
259562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
260c7f610a1SBarry Smith @*/
PCFactorGetShiftAmount(PC pc,PetscReal * shift)261d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift)
262d71ae5a4SJacob Faibussowitsch {
263c7f610a1SBarry Smith   PetscFunctionBegin;
264c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
265cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267c7f610a1SBarry Smith }
268c7f610a1SBarry Smith 
269c7f610a1SBarry Smith /*@
270c7f610a1SBarry Smith   PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
271c7f610a1SBarry Smith 
272c7f610a1SBarry Smith   Not Collective
273c7f610a1SBarry Smith 
274f1580f4eSBarry Smith   Input Parameter:
275c7f610a1SBarry Smith . pc - the preconditioner context
276c7f610a1SBarry Smith 
277c7f610a1SBarry Smith   Output Parameter:
278f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS`
279c7f610a1SBarry Smith 
280c7f610a1SBarry Smith   Level: intermediate
281c7f610a1SBarry Smith 
282562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
283c7f610a1SBarry Smith @*/
PCFactorGetShiftType(PC pc,MatFactorShiftType * type)284d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type)
285d71ae5a4SJacob Faibussowitsch {
286c7f610a1SBarry Smith   PetscFunctionBegin;
287c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
288cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type));
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
290c7f610a1SBarry Smith }
291c7f610a1SBarry Smith 
2922591b318SToby Isaac /*@
2932591b318SToby Isaac   PCFactorGetLevels - Gets the number of levels of fill to use.
2942591b318SToby Isaac 
295c3339decSBarry Smith   Logically Collective
2962591b318SToby Isaac 
2972fe279fdSBarry Smith   Input Parameter:
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 
305562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetLevels()`
3062591b318SToby Isaac @*/
PCFactorGetLevels(PC pc,PetscInt * levels)307d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels)
308d71ae5a4SJacob Faibussowitsch {
3092591b318SToby Isaac   PetscFunctionBegin;
3102591b318SToby Isaac   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
311cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3132591b318SToby Isaac }
3142591b318SToby Isaac 
31585317021SBarry Smith /*@
31685317021SBarry Smith   PCFactorSetLevels - Sets the number of levels of fill to use.
31785317021SBarry Smith 
318c3339decSBarry Smith   Logically Collective
31985317021SBarry Smith 
32085317021SBarry Smith   Input Parameters:
32185317021SBarry Smith + pc     - the preconditioner context
32285317021SBarry Smith - levels - number of levels of fill
32385317021SBarry Smith 
32485317021SBarry Smith   Options Database Key:
32585317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level
32685317021SBarry Smith 
32785317021SBarry Smith   Level: intermediate
32885317021SBarry Smith 
329562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetLevels()`
33085317021SBarry Smith @*/
PCFactorSetLevels(PC pc,PetscInt levels)331d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels)
332d71ae5a4SJacob Faibussowitsch {
33385317021SBarry Smith   PetscFunctionBegin;
3340700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3355f80ce2aSJacob Faibussowitsch   PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels");
336c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, levels, 2);
337cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33985317021SBarry Smith }
34085317021SBarry Smith 
34185317021SBarry Smith /*@
34285317021SBarry Smith   PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
34385317021SBarry Smith   treated as level 0 fill even if there is no non-zero location.
34485317021SBarry Smith 
345c3339decSBarry Smith   Logically Collective
34685317021SBarry Smith 
34785317021SBarry Smith   Input Parameters:
34885317021SBarry Smith + pc  - the preconditioner context
349f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
35085317021SBarry Smith 
35185317021SBarry Smith   Options Database Key:
35267b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill
35385317021SBarry Smith 
354f1580f4eSBarry Smith   Note:
35585317021SBarry Smith   Does not apply with 0 fill.
35685317021SBarry Smith 
35785317021SBarry Smith   Level: intermediate
35885317021SBarry Smith 
359562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()`
36085317021SBarry Smith @*/
PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)361d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg)
362d71ae5a4SJacob Faibussowitsch {
36385317021SBarry Smith   PetscFunctionBegin;
3640700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
365cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg));
3663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36792e9c092SBarry Smith }
36892e9c092SBarry Smith 
36992e9c092SBarry Smith /*@
37092e9c092SBarry Smith   PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
37192e9c092SBarry Smith   treated as level 0 fill even if there is no non-zero location.
37292e9c092SBarry Smith 
373c3339decSBarry Smith   Logically Collective
37492e9c092SBarry Smith 
37592e9c092SBarry Smith   Input Parameter:
37692e9c092SBarry Smith . pc - the preconditioner context
37792e9c092SBarry Smith 
37892e9c092SBarry Smith   Output Parameter:
379f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off
38092e9c092SBarry Smith 
381f1580f4eSBarry Smith   Note:
38292e9c092SBarry Smith   Does not apply with 0 fill.
38392e9c092SBarry Smith 
38492e9c092SBarry Smith   Level: intermediate
38592e9c092SBarry Smith 
386562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()`
38792e9c092SBarry Smith @*/
PCFactorGetAllowDiagonalFill(PC pc,PetscBool * flg)388d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg)
389d71ae5a4SJacob Faibussowitsch {
39092e9c092SBarry Smith   PetscFunctionBegin;
39192e9c092SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
392cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg));
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39485317021SBarry Smith }
39585317021SBarry Smith 
39685317021SBarry Smith /*@
39785317021SBarry Smith   PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
39885317021SBarry Smith 
399c3339decSBarry Smith   Logically Collective
40085317021SBarry Smith 
40185317021SBarry Smith   Input Parameters:
40285317021SBarry Smith + pc   - the preconditioner context
403feefa0e1SJacob Faibussowitsch - rtol - diagonal entries smaller than this in absolute value are considered zero
40485317021SBarry Smith 
40585317021SBarry Smith   Options Database Key:
406147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance
40785317021SBarry Smith 
40885317021SBarry Smith   Level: intermediate
40985317021SBarry Smith 
410ef959800SPierre Jolivet .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftAmount()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
41185317021SBarry Smith @*/
PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)412d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol)
413d71ae5a4SJacob Faibussowitsch {
41485317021SBarry Smith   PetscFunctionBegin;
4150700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
416c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, rtol, 2);
417cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol));
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41985317021SBarry Smith }
42085317021SBarry Smith 
421cc4c1da9SBarry Smith /*@
422f1580f4eSBarry Smith   PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization
42385317021SBarry Smith 
424c3339decSBarry Smith   Logically Collective
42585317021SBarry Smith 
42685317021SBarry Smith   Input Parameters:
42785317021SBarry Smith + pc    - the preconditioner context
428f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
42985317021SBarry Smith 
43085317021SBarry Smith   Options Database Key:
4313ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
43285317021SBarry Smith 
43385317021SBarry Smith   Level: intermediate
43485317021SBarry Smith 
43585317021SBarry Smith   Note:
436998e4596SBarry Smith   The default type is set by searching for available types based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`.
437998e4596SBarry Smith   Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver.
438998e4596SBarry Smith   For example if one configuration had --download-mumps while a different one had --download-superlu_dist.
43985317021SBarry Smith 
440998e4596SBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverTypeRegister()`,
441998e4596SBarry Smith           `MatInitializePackage()`, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `MatSolverTypeGet()`
44285317021SBarry Smith @*/
PCFactorSetMatSolverType(PC pc,MatSolverType stype)443d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype)
444d71ae5a4SJacob Faibussowitsch {
44585317021SBarry Smith   PetscFunctionBegin;
4460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
447cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype));
4483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44985317021SBarry Smith }
45085317021SBarry Smith 
451cc4c1da9SBarry Smith /*@
452f1580f4eSBarry Smith   PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization
4537112b564SBarry Smith 
454c5eb9154SBarry Smith   Not Collective
4557112b564SBarry Smith 
4567112b564SBarry Smith   Input Parameter:
4577112b564SBarry Smith . pc - the preconditioner context
4587112b564SBarry Smith 
4597112b564SBarry Smith   Output Parameter:
460f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
4617112b564SBarry Smith 
4627112b564SBarry Smith   Level: intermediate
4637112b564SBarry Smith 
464562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `MATSOLVERSUPERLU`,
46542747ad1SJacob Faibussowitsch `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`
4667112b564SBarry Smith @*/
PCFactorGetMatSolverType(PC pc,MatSolverType * stype)467d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype)
468d71ae5a4SJacob Faibussowitsch {
4695f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(PC, MatSolverType *);
4707112b564SBarry Smith 
4717112b564SBarry Smith   PetscFunctionBegin;
4720700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4734f572ea9SToby Isaac   PetscAssertPointer(stype, 2);
4749566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f));
4759566063dSJacob Faibussowitsch   if (f) PetscCall((*f)(pc, stype));
4765f80ce2aSJacob Faibussowitsch   else *stype = NULL;
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4787112b564SBarry Smith }
4797112b564SBarry Smith 
48085317021SBarry Smith /*@
48185317021SBarry Smith   PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
48285317021SBarry Smith   fill = number nonzeros in factor/number nonzeros in original matrix.
48385317021SBarry Smith 
484c5eb9154SBarry Smith   Not Collective, each process can expect a different amount of fill
48585317021SBarry Smith 
48685317021SBarry Smith   Input Parameters:
48785317021SBarry Smith + pc   - the preconditioner context
48885317021SBarry Smith - fill - amount of expected fill
48985317021SBarry Smith 
49085317021SBarry Smith   Options Database Key:
49185317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount
49285317021SBarry Smith 
49385317021SBarry Smith   Level: intermediate
49485317021SBarry Smith 
495f1580f4eSBarry Smith   Notes:
49685317021SBarry Smith   For sparse matrix factorizations it is difficult to predict how much
49785317021SBarry Smith   fill to expect. By running with the option -info PETSc will print the
49885317021SBarry Smith   actual amount of fill used; allowing you to set the value accurately for
49985317021SBarry Smith   future runs. Default PETSc uses a value of 5.0
50085317021SBarry Smith 
501f1580f4eSBarry Smith   This is ignored for most solver packages
50201a79839SBarry Smith 
503f1580f4eSBarry Smith   This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels.
504f1580f4eSBarry Smith 
505562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()`
50685317021SBarry Smith @*/
PCFactorSetFill(PC pc,PetscReal fill)507d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill)
508d71ae5a4SJacob Faibussowitsch {
50985317021SBarry Smith   PetscFunctionBegin;
5100700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
511467446fbSPierre Jolivet   PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less than 1.0");
512cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill));
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51485317021SBarry Smith }
51585317021SBarry Smith 
51685317021SBarry Smith /*@
51704c3f3b8SBarry Smith   PCFactorSetUseInPlace - Tells the preconditioner to do an in-place factorization.
51885317021SBarry Smith 
519c3339decSBarry Smith   Logically Collective
52085317021SBarry Smith 
52185317021SBarry Smith   Input Parameters:
5228e37d05fSBarry Smith + pc  - the preconditioner context
523f1580f4eSBarry Smith - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
52485317021SBarry Smith 
52585317021SBarry Smith   Options Database Key:
5268e37d05fSBarry Smith . -pc_factor_in_place <true,false> - Activate/deactivate in-place factorization
52785317021SBarry Smith 
528f1580f4eSBarry Smith   Note:
52904c3f3b8SBarry Smith   For dense matrices, this enables the solution of much larger problems.
53004c3f3b8SBarry Smith   For sparse matrices the factorization cannot be done truly in-place
53104c3f3b8SBarry Smith   so this does not save memory during the factorization, but after the matrix
53204c3f3b8SBarry Smith   is factored, the original unfactored matrix is freed, thus recovering that
53304c3f3b8SBarry Smith   space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
53404c3f3b8SBarry Smith 
535f1580f4eSBarry Smith   `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when
53685317021SBarry Smith   a different matrix is provided for the multiply and the preconditioner in
537f1580f4eSBarry Smith   a call to `KSPSetOperators()`.
53885317021SBarry Smith   This is because the Krylov space methods require an application of the
53985317021SBarry Smith   matrix multiplication, which is not possible here because the matrix has
54085317021SBarry Smith   been factored in-place, replacing the original matrix.
54185317021SBarry Smith 
54285317021SBarry Smith   Level: intermediate
54385317021SBarry Smith 
544562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `Mat`, `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()`
54585317021SBarry Smith @*/
PCFactorSetUseInPlace(PC pc,PetscBool flg)546d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg)
547d71ae5a4SJacob Faibussowitsch {
54885317021SBarry Smith   PetscFunctionBegin;
5490700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
550cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg));
5513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5528e37d05fSBarry Smith }
5538e37d05fSBarry Smith 
5548e37d05fSBarry Smith /*@
5558e37d05fSBarry Smith   PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
5568e37d05fSBarry Smith 
557c3339decSBarry Smith   Logically Collective
5588e37d05fSBarry Smith 
5598e37d05fSBarry Smith   Input Parameter:
5608e37d05fSBarry Smith . pc - the preconditioner context
5618e37d05fSBarry Smith 
5628e37d05fSBarry Smith   Output Parameter:
563f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable
5648e37d05fSBarry Smith 
5658e37d05fSBarry Smith   Level: intermediate
5668e37d05fSBarry Smith 
567562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()`
5688e37d05fSBarry Smith @*/
PCFactorGetUseInPlace(PC pc,PetscBool * flg)569d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg)
570d71ae5a4SJacob Faibussowitsch {
5718e37d05fSBarry Smith   PetscFunctionBegin;
5728e37d05fSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
573cac4c232SBarry Smith   PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg));
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57585317021SBarry Smith }
57685317021SBarry Smith 
5775d83a8b1SBarry Smith /*@
57885317021SBarry Smith   PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
579f1580f4eSBarry Smith   be used in the `PCLU`, `PCCHOLESKY`, `PCILU`,  or `PCICC` preconditioners
58085317021SBarry Smith 
581c3339decSBarry Smith   Logically Collective
58285317021SBarry Smith 
58385317021SBarry Smith   Input Parameters:
58485317021SBarry Smith + pc       - the preconditioner context
585f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM`
58685317021SBarry Smith 
58785317021SBarry Smith   Options Database Key:
5882c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine
58985317021SBarry Smith 
59085317021SBarry Smith   Level: intermediate
59185317021SBarry Smith 
59295452b02SPatrick Sanan   Notes:
5934ac6704cSBarry Smith   Nested dissection is used by default for some of PETSc's sparse matrix formats
59485317021SBarry Smith 
595f1580f4eSBarry 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
5969bd791bbSBarry Smith   and reordering this matrix is very expensive.
5979bd791bbSBarry Smith 
598f1580f4eSBarry Smith   You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering.
5999bd791bbSBarry Smith 
600f1580f4eSBarry Smith   `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others.
6012c7c0729SBarry Smith 
602562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM`
60385317021SBarry Smith @*/
PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)604d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering)
605d71ae5a4SJacob Faibussowitsch {
60685317021SBarry Smith   PetscFunctionBegin;
607c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
608cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering));
6093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61085317021SBarry Smith }
61185317021SBarry Smith 
61285317021SBarry Smith /*@
6138ff23777SHong Zhang   PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
61485317021SBarry Smith   For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
615f1580f4eSBarry Smith   it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used.
61685317021SBarry Smith 
617c3339decSBarry Smith   Logically Collective
61885317021SBarry Smith 
61985317021SBarry Smith   Input Parameters:
62085317021SBarry Smith + pc    - the preconditioner context
62185317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
62285317021SBarry Smith 
62385317021SBarry Smith   Options Database Key:
624147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance
62585317021SBarry Smith 
62685317021SBarry Smith   Level: intermediate
62785317021SBarry Smith 
628562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
62985317021SBarry Smith @*/
PCFactorSetColumnPivot(PC pc,PetscReal dtcol)630d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol)
631d71ae5a4SJacob Faibussowitsch {
63285317021SBarry Smith   PetscFunctionBegin;
633c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
634c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, dtcol, 2);
635cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol));
6363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63785317021SBarry Smith }
63885317021SBarry Smith 
63985317021SBarry Smith /*@
64085317021SBarry Smith   PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
641f1580f4eSBarry Smith   with `MATBAIJ` or `MATSBAIJ` matrices
64285317021SBarry Smith 
643c3339decSBarry Smith   Logically Collective
64485317021SBarry Smith 
64585317021SBarry Smith   Input Parameters:
64685317021SBarry Smith + pc    - the preconditioner context
647f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE`
64885317021SBarry Smith 
64985317021SBarry Smith   Options Database Key:
650f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ`
65185317021SBarry Smith 
65285317021SBarry Smith   Level: intermediate
65385317021SBarry Smith 
654562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
65585317021SBarry Smith @*/
PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)656d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot)
657d71ae5a4SJacob Faibussowitsch {
65885317021SBarry Smith   PetscFunctionBegin;
659c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
660acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, pivot, 2);
661cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot));
6623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66385317021SBarry Smith }
66485317021SBarry Smith 
66585317021SBarry Smith /*@
666288e7d53SBarry Smith   PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
66785317021SBarry Smith   this causes later ones to use the fill ratio computed in the initial factorization.
66885317021SBarry Smith 
669c3339decSBarry Smith   Logically Collective
67085317021SBarry Smith 
67185317021SBarry Smith   Input Parameters:
67285317021SBarry Smith + pc   - the preconditioner context
673f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
67485317021SBarry Smith 
67585317021SBarry Smith   Options Database Key:
676f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
67785317021SBarry Smith 
67885317021SBarry Smith   Level: intermediate
67985317021SBarry Smith 
680562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()`
68185317021SBarry Smith @*/
PCFactorSetReuseFill(PC pc,PetscBool flag)682d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag)
683d71ae5a4SJacob Faibussowitsch {
68485317021SBarry Smith   PetscFunctionBegin;
685064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
686acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, flag, 2);
687cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68985317021SBarry Smith }
6903d1c1ea0SBarry Smith 
PCFactorInitialize(PC pc,MatFactorType ftype)691d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype)
692d71ae5a4SJacob Faibussowitsch {
6933d1c1ea0SBarry Smith   PC_Factor *fact = (PC_Factor *)pc->data;
6943d1c1ea0SBarry Smith 
6953d1c1ea0SBarry Smith   PetscFunctionBegin;
6969566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&fact->info));
6974ac6704cSBarry Smith   fact->factortype           = ftype;
6983d1c1ea0SBarry Smith   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
6993d1c1ea0SBarry Smith   fact->info.shiftamount     = 100.0 * PETSC_MACHINE_EPSILON;
7003d1c1ea0SBarry Smith   fact->info.zeropivot       = 100.0 * PETSC_MACHINE_EPSILON;
7013d1c1ea0SBarry Smith   fact->info.pivotinblocks   = 1.0;
7023d1c1ea0SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
7033d1c1ea0SBarry Smith 
7049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor));
7059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor));
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor));
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor));
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor));
7099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor));
7109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor));
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor));
7129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor));
7139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor));
7149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor));
7159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor));
7169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor));
7179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor));
7189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor));
7199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor));
7209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor));
7219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor));
7229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor));
7239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor));
7243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7253d1c1ea0SBarry Smith }
7262e956fe4SStefano Zampini 
PCFactorClearComposedFunctions(PC pc)727d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc)
728d71ae5a4SJacob Faibussowitsch {
7292e956fe4SStefano Zampini   PetscFunctionBegin;
7302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL));
7312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL));
7322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL));
7332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL));
7342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL));
7352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL));
7362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL));
7372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL));
7382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL));
7392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL));
7402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL));
7412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL));
7422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL));
7432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL));
7442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL));
7452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL));
7462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL));
7472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL));
7482e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL));
7492e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL));
7502e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL));
7512e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL));
7523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7532e956fe4SStefano Zampini }
754