xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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 */
84ac6704cSBarry Smith PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
94ac6704cSBarry Smith {
10978beb7fSPierre Jolivet   Mat            B;
11978beb7fSPierre Jolivet   PetscBool      foundmtype,flg;
124ac6704cSBarry Smith 
134ac6704cSBarry Smith   PetscFunctionBegin;
144ac6704cSBarry Smith   if (pc->pmat) {
154ac6704cSBarry Smith     PC_Factor *fact = (PC_Factor*)pc->data;
169566063dSJacob Faibussowitsch     PetscCall(MatSolverTypeGet(fact->solvertype,((PetscObject)pc->pmat)->type_name,fact->factortype,NULL,&foundmtype,NULL));
17978beb7fSPierre Jolivet     if (foundmtype) {
184ac6704cSBarry Smith       if (!fact->fact) {
199566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat,fact->solvertype,fact->factortype,&fact->fact));
20978beb7fSPierre Jolivet       } else if (!fact->fact->assembled) {
219566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(fact->solvertype,fact->fact->solvertype,&flg));
22978beb7fSPierre Jolivet         if (!flg) {
239566063dSJacob Faibussowitsch           PetscCall(MatGetFactor(pc->pmat,fact->solvertype,fact->factortype,&B));
249566063dSJacob Faibussowitsch           PetscCall(MatHeaderReplace(fact->fact,&B));
25978beb7fSPierre Jolivet         }
264ac6704cSBarry Smith       }
274ac6704cSBarry Smith       if (!fact->ordering) {
284ac6704cSBarry Smith         PetscBool       canuseordering;
294ac6704cSBarry Smith         MatOrderingType otype;
304ac6704cSBarry Smith 
319566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(fact->fact,&canuseordering));
324ac6704cSBarry Smith         if (canuseordering) {
339566063dSJacob Faibussowitsch           PetscCall(MatFactorGetPreferredOrdering(fact->fact,fact->factortype,&otype));
344ac6704cSBarry Smith         } else otype = MATORDERINGEXTERNAL;
359566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(otype,(char **)&fact->ordering));
364ac6704cSBarry Smith       }
374ac6704cSBarry Smith     }
38978beb7fSPierre Jolivet   }
394ac6704cSBarry Smith   PetscFunctionReturn(0);
404ac6704cSBarry Smith }
414ac6704cSBarry Smith 
423d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag)
433d1c1ea0SBarry Smith {
443d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
453d1c1ea0SBarry Smith 
463d1c1ea0SBarry Smith   PetscFunctionBegin;
473d1c1ea0SBarry Smith   lu->reuseordering = flag;
483d1c1ea0SBarry Smith   PetscFunctionReturn(0);
493d1c1ea0SBarry Smith }
503d1c1ea0SBarry Smith 
513d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc,PetscBool flag)
523d1c1ea0SBarry Smith {
533d1c1ea0SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
543d1c1ea0SBarry Smith 
553d1c1ea0SBarry Smith   PetscFunctionBegin;
563d1c1ea0SBarry Smith   lu->reusefill = flag;
573d1c1ea0SBarry Smith   PetscFunctionReturn(0);
583d1c1ea0SBarry Smith }
593d1c1ea0SBarry Smith 
603d1c1ea0SBarry Smith static PetscErrorCode  PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg)
613d1c1ea0SBarry Smith {
623d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
633d1c1ea0SBarry Smith 
643d1c1ea0SBarry Smith   PetscFunctionBegin;
653d1c1ea0SBarry Smith   dir->inplace = flg;
663d1c1ea0SBarry Smith   PetscFunctionReturn(0);
673d1c1ea0SBarry Smith }
683d1c1ea0SBarry Smith 
693d1c1ea0SBarry Smith static PetscErrorCode  PCFactorGetUseInPlace_Factor(PC pc,PetscBool *flg)
703d1c1ea0SBarry Smith {
713d1c1ea0SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
723d1c1ea0SBarry Smith 
733d1c1ea0SBarry Smith   PetscFunctionBegin;
743d1c1ea0SBarry Smith   *flg = dir->inplace;
753d1c1ea0SBarry Smith   PetscFunctionReturn(0);
763d1c1ea0SBarry Smith }
773d1c1ea0SBarry Smith 
78f8260c8fSBarry Smith /*@
793ca39a21SBarry Smith     PCFactorSetUpMatSolverType - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
80f8260c8fSBarry Smith        set the options for that particular factorization object.
81f8260c8fSBarry Smith 
82f8260c8fSBarry Smith   Input Parameter:
83f8260c8fSBarry Smith .  pc  - the preconditioner context
84f8260c8fSBarry Smith 
8595452b02SPatrick Sanan   Notes:
8695452b02SPatrick Sanan     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.
87f8260c8fSBarry Smith 
882bd2b0e6SSatish Balay   Level: intermediate
892bd2b0e6SSatish Balay 
90db781477SPatrick Sanan .seealso: `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
91f8260c8fSBarry Smith @*/
923ca39a21SBarry Smith PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
93f8260c8fSBarry Smith {
94f8260c8fSBarry Smith   PetscFunctionBegin;
95f8260c8fSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
96cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetUpMatSolverType_C",(PC),(pc));
97b3a44c85SBarry Smith   PetscFunctionReturn(0);
98f8260c8fSBarry Smith }
99f8260c8fSBarry Smith 
100ee45ca4aSHong Zhang /*@
101ee45ca4aSHong Zhang    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
102ee45ca4aSHong Zhang 
103ad4df100SBarry Smith    Logically Collective on PC
104ee45ca4aSHong Zhang 
105ee45ca4aSHong Zhang    Input Parameters:
106afaefe49SHong Zhang +  pc - the preconditioner context
107afaefe49SHong Zhang -  zero - all pivots smaller than this will be considered zero
108ee45ca4aSHong Zhang 
109ee45ca4aSHong Zhang    Options Database Key:
110ee45ca4aSHong Zhang .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
111ee45ca4aSHong Zhang 
112ee45ca4aSHong Zhang    Level: intermediate
113ee45ca4aSHong Zhang 
114db781477SPatrick Sanan .seealso: `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
115ee45ca4aSHong Zhang @*/
1167087cfbeSBarry Smith PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
117ee45ca4aSHong Zhang {
118ee45ca4aSHong Zhang   PetscFunctionBegin;
1190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,zero,2);
121cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));
122ee45ca4aSHong Zhang   PetscFunctionReturn(0);
123ee45ca4aSHong Zhang }
124ee45ca4aSHong Zhang 
125915743fcSHong Zhang /*@
126915743fcSHong Zhang    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
127915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
128915743fcSHong Zhang 
129ad4df100SBarry Smith    Logically Collective on PC
130915743fcSHong Zhang 
131915743fcSHong Zhang    Input Parameters:
132915743fcSHong Zhang +  pc - the preconditioner context
133915743fcSHong Zhang -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
134915743fcSHong Zhang 
135915743fcSHong Zhang    Options Database Key:
13628d58a37SPierre Jolivet .  -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
137915743fcSHong Zhang 
138915743fcSHong Zhang    Level: intermediate
139915743fcSHong Zhang 
140db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
141915743fcSHong Zhang @*/
1427087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
143d90ac83dSHong Zhang {
144d90ac83dSHong Zhang   PetscFunctionBegin;
1450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
146c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,shifttype,2);
147cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));
148d90ac83dSHong Zhang   PetscFunctionReturn(0);
149d90ac83dSHong Zhang }
150d90ac83dSHong Zhang 
151915743fcSHong Zhang /*@
152915743fcSHong Zhang    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
153915743fcSHong Zhang      numerical factorization, thus the matrix has nonzero pivots
154915743fcSHong Zhang 
155ad4df100SBarry Smith    Logically Collective on PC
156915743fcSHong Zhang 
157915743fcSHong Zhang    Input Parameters:
158915743fcSHong Zhang +  pc - the preconditioner context
159915743fcSHong Zhang -  shiftamount - amount of shift
160915743fcSHong Zhang 
161915743fcSHong Zhang    Options Database Key:
162915743fcSHong Zhang .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
163915743fcSHong Zhang 
164915743fcSHong Zhang    Level: intermediate
165915743fcSHong Zhang 
166db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
167915743fcSHong Zhang @*/
1687087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
169d90ac83dSHong Zhang {
170d90ac83dSHong Zhang   PetscFunctionBegin;
1710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
172c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,shiftamount,2);
173cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));
174d90ac83dSHong Zhang   PetscFunctionReturn(0);
175d90ac83dSHong Zhang }
176d90ac83dSHong Zhang 
1776d33885cSprj- /*@
178b7c853c4SBarry Smith    PCFactorSetDropTolerance - The preconditioner will use an ILU
17978fc6b22SHong Zhang    based on a drop tolerance. (Under development)
18085317021SBarry Smith 
181ad4df100SBarry Smith    Logically Collective on PC
18285317021SBarry Smith 
18385317021SBarry Smith    Input Parameters:
18485317021SBarry Smith +  pc - the preconditioner context
18585317021SBarry Smith .  dt - the drop tolerance, try from 1.e-10 to .1
18685317021SBarry Smith .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
18785317021SBarry Smith -  maxrowcount - the max number of nonzeros allowed in a row, best value
18885317021SBarry Smith                  depends on the number of nonzeros in row of original matrix
18985317021SBarry Smith 
19085317021SBarry Smith    Options Database Key:
191b7c853c4SBarry Smith .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
19285317021SBarry Smith 
19385317021SBarry Smith    Level: intermediate
19485317021SBarry Smith 
19585317021SBarry Smith       There are NO default values for the 3 parameters, you must set them with reasonable values for your
19685317021SBarry Smith       matrix. We don't know how to compute reasonable values.
19785317021SBarry Smith 
1986d33885cSprj- @*/
1997087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
20085317021SBarry Smith {
20185317021SBarry Smith   PetscFunctionBegin;
2020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
203064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveReal(pc,dtcol,3);
204064a246eSJacob Faibussowitsch   PetscValidLogicalCollectiveInt(pc,maxrowcount,4);
205cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));
20685317021SBarry Smith   PetscFunctionReturn(0);
20785317021SBarry Smith }
20885317021SBarry Smith 
209c7f610a1SBarry Smith /*@
210c7f610a1SBarry Smith    PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
211c7f610a1SBarry Smith 
212c7f610a1SBarry Smith    Not Collective
213c7f610a1SBarry Smith 
214c7f610a1SBarry Smith    Input Parameters:
215c7f610a1SBarry Smith .  pc - the preconditioner context
216c7f610a1SBarry Smith 
217c7f610a1SBarry Smith    Output Parameter:
218c7f610a1SBarry Smith .  pivot - the tolerance
219c7f610a1SBarry Smith 
220c7f610a1SBarry Smith    Level: intermediate
221c7f610a1SBarry Smith 
222db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()`
223c7f610a1SBarry Smith @*/
224c7f610a1SBarry Smith PetscErrorCode  PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
225c7f610a1SBarry Smith {
226c7f610a1SBarry Smith   PetscFunctionBegin;
227c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
228cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));
229c7f610a1SBarry Smith   PetscFunctionReturn(0);
230c7f610a1SBarry Smith }
231c7f610a1SBarry Smith 
232c7f610a1SBarry Smith /*@
233c7f610a1SBarry Smith    PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
234c7f610a1SBarry Smith 
235c7f610a1SBarry Smith    Not Collective
236c7f610a1SBarry Smith 
237c7f610a1SBarry Smith    Input Parameters:
238c7f610a1SBarry Smith .  pc - the preconditioner context
239c7f610a1SBarry Smith 
240c7f610a1SBarry Smith    Output Parameter:
241c7f610a1SBarry Smith .  shift - how much to shift the diagonal entry
242c7f610a1SBarry Smith 
243c7f610a1SBarry Smith    Level: intermediate
244c7f610a1SBarry Smith 
245db781477SPatrick Sanan .seealso: `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
246c7f610a1SBarry Smith @*/
247c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftAmount(PC pc,PetscReal *shift)
248c7f610a1SBarry Smith {
249c7f610a1SBarry Smith   PetscFunctionBegin;
250c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
251cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));
252c7f610a1SBarry Smith   PetscFunctionReturn(0);
253c7f610a1SBarry Smith }
254c7f610a1SBarry Smith 
255c7f610a1SBarry Smith /*@
256c7f610a1SBarry Smith    PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
257c7f610a1SBarry Smith 
258c7f610a1SBarry Smith    Not Collective
259c7f610a1SBarry Smith 
260c7f610a1SBarry Smith    Input Parameters:
261c7f610a1SBarry Smith .  pc - the preconditioner context
262c7f610a1SBarry Smith 
263c7f610a1SBarry Smith    Output Parameter:
264c7f610a1SBarry Smith .  type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS
265c7f610a1SBarry Smith 
266c7f610a1SBarry Smith    Level: intermediate
267c7f610a1SBarry Smith 
268db781477SPatrick Sanan .seealso: `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
269c7f610a1SBarry Smith @*/
270c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
271c7f610a1SBarry Smith {
272c7f610a1SBarry Smith   PetscFunctionBegin;
273c7f610a1SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
274cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));
275c7f610a1SBarry Smith   PetscFunctionReturn(0);
276c7f610a1SBarry Smith }
277c7f610a1SBarry Smith 
2782591b318SToby Isaac /*@
2792591b318SToby Isaac    PCFactorGetLevels - Gets the number of levels of fill to use.
2802591b318SToby Isaac 
2812591b318SToby Isaac    Logically Collective on PC
2822591b318SToby Isaac 
2832591b318SToby Isaac    Input Parameters:
2842591b318SToby Isaac .  pc - the preconditioner context
2852591b318SToby Isaac 
2862591b318SToby Isaac    Output Parameter:
2872591b318SToby Isaac .  levels - number of levels of fill
2882591b318SToby Isaac 
2892591b318SToby Isaac    Level: intermediate
2902591b318SToby Isaac 
2912591b318SToby Isaac @*/
2922591b318SToby Isaac PetscErrorCode  PCFactorGetLevels(PC pc,PetscInt *levels)
2932591b318SToby Isaac {
2942591b318SToby Isaac   PetscFunctionBegin;
2952591b318SToby Isaac   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
296cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
2972591b318SToby Isaac   PetscFunctionReturn(0);
2982591b318SToby Isaac }
2992591b318SToby Isaac 
30085317021SBarry Smith /*@
30185317021SBarry Smith    PCFactorSetLevels - Sets the number of levels of fill to use.
30285317021SBarry Smith 
303ad4df100SBarry Smith    Logically Collective on PC
30485317021SBarry Smith 
30585317021SBarry Smith    Input Parameters:
30685317021SBarry Smith +  pc - the preconditioner context
30785317021SBarry Smith -  levels - number of levels of fill
30885317021SBarry Smith 
30985317021SBarry Smith    Options Database Key:
31085317021SBarry Smith .  -pc_factor_levels <levels> - Sets fill level
31185317021SBarry Smith 
31285317021SBarry Smith    Level: intermediate
31385317021SBarry Smith 
31485317021SBarry Smith @*/
3157087cfbeSBarry Smith PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
31685317021SBarry Smith {
31785317021SBarry Smith   PetscFunctionBegin;
3180700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3195f80ce2aSJacob Faibussowitsch   PetscCheck(levels >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
320c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,levels,2);
321cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
32285317021SBarry Smith   PetscFunctionReturn(0);
32385317021SBarry Smith }
32485317021SBarry Smith 
32585317021SBarry Smith /*@
32685317021SBarry Smith    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
32785317021SBarry Smith    treated as level 0 fill even if there is no non-zero location.
32885317021SBarry Smith 
329ad4df100SBarry Smith    Logically Collective on PC
33085317021SBarry Smith 
33185317021SBarry Smith    Input Parameters:
33285317021SBarry Smith +  pc - the preconditioner context
33392e9c092SBarry Smith -  flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
33485317021SBarry Smith 
33585317021SBarry Smith    Options Database Key:
33667b8a455SSatish Balay .  -pc_factor_diagonal_fill <bool> - allow the diagonal fill
33785317021SBarry Smith 
33885317021SBarry Smith    Notes:
33985317021SBarry Smith    Does not apply with 0 fill.
34085317021SBarry Smith 
34185317021SBarry Smith    Level: intermediate
34285317021SBarry Smith 
343db781477SPatrick Sanan .seealso: `PCFactorGetAllowDiagonalFill()`
34485317021SBarry Smith @*/
34592e9c092SBarry Smith PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
34685317021SBarry Smith {
34785317021SBarry Smith   PetscFunctionBegin;
3480700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
349cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
35092e9c092SBarry Smith   PetscFunctionReturn(0);
35192e9c092SBarry Smith }
35292e9c092SBarry Smith 
35392e9c092SBarry Smith /*@
35492e9c092SBarry Smith    PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
35592e9c092SBarry Smith        treated as level 0 fill even if there is no non-zero location.
35692e9c092SBarry Smith 
35792e9c092SBarry Smith    Logically Collective on PC
35892e9c092SBarry Smith 
35992e9c092SBarry Smith    Input Parameter:
36092e9c092SBarry Smith .  pc - the preconditioner context
36192e9c092SBarry Smith 
36292e9c092SBarry Smith    Output Parameter:
36392e9c092SBarry Smith .   flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
36492e9c092SBarry Smith 
36592e9c092SBarry Smith    Notes:
36692e9c092SBarry Smith    Does not apply with 0 fill.
36792e9c092SBarry Smith 
36892e9c092SBarry Smith    Level: intermediate
36992e9c092SBarry Smith 
370db781477SPatrick Sanan .seealso: `PCFactorSetAllowDiagonalFill()`
37192e9c092SBarry Smith @*/
37292e9c092SBarry Smith PetscErrorCode  PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
37392e9c092SBarry Smith {
37492e9c092SBarry Smith   PetscFunctionBegin;
37592e9c092SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
376cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
37785317021SBarry Smith   PetscFunctionReturn(0);
37885317021SBarry Smith }
37985317021SBarry Smith 
38085317021SBarry Smith /*@
38185317021SBarry Smith    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
38285317021SBarry Smith 
383ad4df100SBarry Smith    Logically Collective on PC
38485317021SBarry Smith 
38585317021SBarry Smith    Input Parameters:
38685317021SBarry Smith +  pc - the preconditioner context
38785317021SBarry Smith -  tol - diagonal entries smaller than this in absolute value are considered zero
38885317021SBarry Smith 
38985317021SBarry Smith    Options Database Key:
390147403d9SBarry Smith .  -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance
39185317021SBarry Smith 
39285317021SBarry Smith    Level: intermediate
39385317021SBarry Smith 
394db781477SPatrick Sanan .seealso: `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
39585317021SBarry Smith @*/
3967087cfbeSBarry Smith PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
39785317021SBarry Smith {
39885317021SBarry Smith   PetscFunctionBegin;
3990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
400c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,rtol,2);
401cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
40285317021SBarry Smith   PetscFunctionReturn(0);
40385317021SBarry Smith }
40485317021SBarry Smith 
405bf6011e8SBarry Smith /*@C
4063ca39a21SBarry Smith    PCFactorSetMatSolverType - sets the software that is used to perform the factorization
40785317021SBarry Smith 
408ad4df100SBarry Smith    Logically Collective on PC
40985317021SBarry Smith 
41085317021SBarry Smith    Input Parameters:
41185317021SBarry Smith +  pc - the preconditioner context
412f60c3dc2SHong Zhang -  stype - for example, superlu, superlu_dist
41385317021SBarry Smith 
41485317021SBarry Smith    Options Database Key:
4153ca39a21SBarry Smith .  -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
41685317021SBarry Smith 
41785317021SBarry Smith    Level: intermediate
41885317021SBarry Smith 
41985317021SBarry Smith    Note:
42085317021SBarry Smith      By default this will use the PETSc factorization if it exists
42185317021SBarry Smith 
422db781477SPatrick Sanan .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`
42385317021SBarry Smith @*/
424ea799195SBarry Smith PetscErrorCode  PCFactorSetMatSolverType(PC pc,MatSolverType stype)
42585317021SBarry Smith {
42685317021SBarry Smith   PetscFunctionBegin;
4270700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
428cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetMatSolverType_C",(PC,MatSolverType),(pc,stype));
42985317021SBarry Smith   PetscFunctionReturn(0);
43085317021SBarry Smith }
43185317021SBarry Smith 
432bf6011e8SBarry Smith /*@C
4333ca39a21SBarry Smith    PCFactorGetMatSolverType - gets the software that is used to perform the factorization
4347112b564SBarry Smith 
435c5eb9154SBarry Smith    Not Collective
4367112b564SBarry Smith 
4377112b564SBarry Smith    Input Parameter:
4387112b564SBarry Smith .  pc - the preconditioner context
4397112b564SBarry Smith 
4407112b564SBarry Smith    Output Parameter:
4410298fd71SBarry Smith .   stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
4427112b564SBarry Smith 
4437112b564SBarry Smith    Level: intermediate
4447112b564SBarry Smith 
445db781477SPatrick Sanan .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`
4467112b564SBarry Smith @*/
447ea799195SBarry Smith PetscErrorCode  PCFactorGetMatSolverType(PC pc,MatSolverType *stype)
4487112b564SBarry Smith {
4495f80ce2aSJacob Faibussowitsch   PetscErrorCode (*f)(PC,MatSolverType*);
4507112b564SBarry Smith 
4517112b564SBarry Smith   PetscFunctionBegin;
4520700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4535f80ce2aSJacob Faibussowitsch   PetscValidPointer(stype,2);
4549566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",&f));
4559566063dSJacob Faibussowitsch   if (f) PetscCall((*f)(pc,stype));
4565f80ce2aSJacob Faibussowitsch   else *stype = NULL;
4577112b564SBarry Smith   PetscFunctionReturn(0);
4587112b564SBarry Smith }
4597112b564SBarry Smith 
46085317021SBarry Smith /*@
46185317021SBarry Smith    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
46285317021SBarry Smith    fill = number nonzeros in factor/number nonzeros in original matrix.
46385317021SBarry Smith 
464c5eb9154SBarry Smith    Not Collective, each process can expect a different amount of fill
46585317021SBarry Smith 
46685317021SBarry Smith    Input Parameters:
46785317021SBarry Smith +  pc - the preconditioner context
46885317021SBarry Smith -  fill - amount of expected fill
46985317021SBarry Smith 
47085317021SBarry Smith    Options Database Key:
47185317021SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
47285317021SBarry Smith 
47385317021SBarry Smith    Level: intermediate
47485317021SBarry Smith 
47585317021SBarry Smith    Note:
47685317021SBarry Smith    For sparse matrix factorizations it is difficult to predict how much
47785317021SBarry Smith    fill to expect. By running with the option -info PETSc will print the
47885317021SBarry Smith    actual amount of fill used; allowing you to set the value accurately for
47985317021SBarry Smith    future runs. Default PETSc uses a value of 5.0
48085317021SBarry Smith 
48101a79839SBarry Smith    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
48201a79839SBarry Smith 
48385317021SBarry Smith @*/
4847087cfbeSBarry Smith PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
48585317021SBarry Smith {
48685317021SBarry Smith   PetscFunctionBegin;
4870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4885f80ce2aSJacob Faibussowitsch   PetscCheck(fill >= 1.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
489cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
49085317021SBarry Smith   PetscFunctionReturn(0);
49185317021SBarry Smith }
49285317021SBarry Smith 
49385317021SBarry Smith /*@
49485317021SBarry Smith    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
49585317021SBarry Smith    For dense matrices, this enables the solution of much larger problems.
49685317021SBarry Smith    For sparse matrices the factorization cannot be done truly in-place
49785317021SBarry Smith    so this does not save memory during the factorization, but after the matrix
49885317021SBarry Smith    is factored, the original unfactored matrix is freed, thus recovering that
499ec5066bdSBarry Smith    space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
50085317021SBarry Smith 
501ad4df100SBarry Smith    Logically Collective on PC
50285317021SBarry Smith 
50385317021SBarry Smith    Input Parameters:
5048e37d05fSBarry Smith +  pc - the preconditioner context
5058e37d05fSBarry Smith -  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
50685317021SBarry Smith 
50785317021SBarry Smith    Options Database Key:
5088e37d05fSBarry Smith .  -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
50985317021SBarry Smith 
51085317021SBarry Smith    Notes:
51185317021SBarry Smith    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
51285317021SBarry Smith    a different matrix is provided for the multiply and the preconditioner in
51385317021SBarry Smith    a call to KSPSetOperators().
51485317021SBarry Smith    This is because the Krylov space methods require an application of the
51585317021SBarry Smith    matrix multiplication, which is not possible here because the matrix has
51685317021SBarry Smith    been factored in-place, replacing the original matrix.
51785317021SBarry Smith 
51885317021SBarry Smith    Level: intermediate
51985317021SBarry Smith 
520db781477SPatrick Sanan .seealso: `PCFactorGetUseInPlace()`
52185317021SBarry Smith @*/
5228e37d05fSBarry Smith PetscErrorCode  PCFactorSetUseInPlace(PC pc,PetscBool flg)
52385317021SBarry Smith {
52485317021SBarry Smith   PetscFunctionBegin;
5250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
526cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
5278e37d05fSBarry Smith   PetscFunctionReturn(0);
5288e37d05fSBarry Smith }
5298e37d05fSBarry Smith 
5308e37d05fSBarry Smith /*@
5318e37d05fSBarry Smith    PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
5328e37d05fSBarry Smith 
5338e37d05fSBarry Smith    Logically Collective on PC
5348e37d05fSBarry Smith 
5358e37d05fSBarry Smith    Input Parameter:
5368e37d05fSBarry Smith .  pc - the preconditioner context
5378e37d05fSBarry Smith 
5388e37d05fSBarry Smith    Output Parameter:
5398e37d05fSBarry Smith .  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
5408e37d05fSBarry Smith 
5418e37d05fSBarry Smith    Level: intermediate
5428e37d05fSBarry Smith 
543db781477SPatrick Sanan .seealso: `PCFactorSetUseInPlace()`
5448e37d05fSBarry Smith @*/
5458e37d05fSBarry Smith PetscErrorCode  PCFactorGetUseInPlace(PC pc,PetscBool *flg)
5468e37d05fSBarry Smith {
5478e37d05fSBarry Smith   PetscFunctionBegin;
5488e37d05fSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
549cac4c232SBarry Smith   PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
55085317021SBarry Smith   PetscFunctionReturn(0);
55185317021SBarry Smith }
55285317021SBarry Smith 
55385317021SBarry Smith /*@C
55485317021SBarry Smith     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
5552c7c0729SBarry Smith     be used in the LU, ILU, Cholesky, and ICC factorizations.
55685317021SBarry Smith 
557ad4df100SBarry Smith     Logically Collective on PC
55885317021SBarry Smith 
55985317021SBarry Smith     Input Parameters:
56085317021SBarry Smith +   pc - the preconditioner context
5612692d6eeSBarry Smith -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
56285317021SBarry Smith 
56385317021SBarry Smith     Options Database Key:
5642c7c0729SBarry Smith .   -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine
56585317021SBarry Smith 
56685317021SBarry Smith     Level: intermediate
56785317021SBarry Smith 
56895452b02SPatrick Sanan     Notes:
5694ac6704cSBarry Smith       Nested dissection is used by default for some of PETSc's sparse matrix formats
57085317021SBarry Smith 
5719bd791bbSBarry Smith      For Cholesky and ICC and the SBAIJ format the only reordering available is natural since only the upper half of the matrix is stored
5729bd791bbSBarry Smith      and reordering this matrix is very expensive.
5739bd791bbSBarry Smith 
5744ac6704cSBarry Smith       You can use a SeqAIJ matrix with Cholesky and ICC and use any ordering.
5759bd791bbSBarry Smith 
5764ac6704cSBarry Smith       MATORDERINGEXTERNAL means PETSc will not compute an ordering and the package will use its own ordering, usable with MATSOLVERCHOLMOD, MATSOLVERUMFPACK, and others.
5772c7c0729SBarry Smith 
578db781477SPatrick Sanan .seealso: `MatOrderingType`
57985317021SBarry Smith 
58085317021SBarry Smith @*/
58119fd82e9SBarry Smith PetscErrorCode  PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
58285317021SBarry Smith {
58385317021SBarry Smith   PetscFunctionBegin;
584c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
585cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));
58685317021SBarry Smith   PetscFunctionReturn(0);
58785317021SBarry Smith }
58885317021SBarry Smith 
58985317021SBarry Smith /*@
5908ff23777SHong Zhang     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
59185317021SBarry Smith       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
592e3c5b3baSBarry Smith       it is never done. For the MATLAB and SuperLU factorization this is used.
59385317021SBarry Smith 
594ad4df100SBarry Smith     Logically Collective on PC
59585317021SBarry Smith 
59685317021SBarry Smith     Input Parameters:
59785317021SBarry Smith +   pc - the preconditioner context
59885317021SBarry Smith -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
59985317021SBarry Smith 
60085317021SBarry Smith     Options Database Key:
601147403d9SBarry Smith .   -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance
60285317021SBarry Smith 
60385317021SBarry Smith     Level: intermediate
60485317021SBarry Smith 
605db781477SPatrick Sanan .seealso: `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
60685317021SBarry Smith @*/
6077087cfbeSBarry Smith PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
60885317021SBarry Smith {
60985317021SBarry Smith   PetscFunctionBegin;
610c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
611c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,dtcol,2);
612cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
61385317021SBarry Smith   PetscFunctionReturn(0);
61485317021SBarry Smith }
61585317021SBarry Smith 
61685317021SBarry Smith /*@
61785317021SBarry Smith     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
61885317021SBarry Smith       with BAIJ or SBAIJ matrices
61985317021SBarry Smith 
620ad4df100SBarry Smith     Logically Collective on PC
62185317021SBarry Smith 
62285317021SBarry Smith     Input Parameters:
62385317021SBarry Smith +   pc - the preconditioner context
62485317021SBarry Smith -   pivot - PETSC_TRUE or PETSC_FALSE
62585317021SBarry Smith 
62685317021SBarry Smith     Options Database Key:
62767b8a455SSatish Balay .   -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for BAIJ and SBAIJ
62885317021SBarry Smith 
62985317021SBarry Smith     Level: intermediate
63085317021SBarry Smith 
631db781477SPatrick Sanan .seealso: `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
63285317021SBarry Smith @*/
6337087cfbeSBarry Smith PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
63485317021SBarry Smith {
63585317021SBarry Smith   PetscFunctionBegin;
636c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
637acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,pivot,2);
638cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
63985317021SBarry Smith   PetscFunctionReturn(0);
64085317021SBarry Smith }
64185317021SBarry Smith 
64285317021SBarry Smith /*@
643288e7d53SBarry Smith    PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
64485317021SBarry Smith    this causes later ones to use the fill ratio computed in the initial factorization.
64585317021SBarry Smith 
646ad4df100SBarry Smith    Logically Collective on PC
64785317021SBarry Smith 
64885317021SBarry Smith    Input Parameters:
64985317021SBarry Smith +  pc - the preconditioner context
65085317021SBarry Smith -  flag - PETSC_TRUE to reuse else PETSC_FALSE
65185317021SBarry Smith 
65285317021SBarry Smith    Options Database Key:
65385317021SBarry Smith .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
65485317021SBarry Smith 
65585317021SBarry Smith    Level: intermediate
65685317021SBarry Smith 
657db781477SPatrick Sanan .seealso: `PCFactorSetReuseOrdering()`
65885317021SBarry Smith @*/
6597087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool flag)
66085317021SBarry Smith {
66185317021SBarry Smith   PetscFunctionBegin;
662064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
663acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,flag,2);
664cac4c232SBarry Smith   PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
66585317021SBarry Smith   PetscFunctionReturn(0);
66685317021SBarry Smith }
6673d1c1ea0SBarry Smith 
6684ac6704cSBarry Smith PetscErrorCode PCFactorInitialize(PC pc,MatFactorType ftype)
6693d1c1ea0SBarry Smith {
6703d1c1ea0SBarry Smith   PC_Factor      *fact = (PC_Factor*)pc->data;
6713d1c1ea0SBarry Smith 
6723d1c1ea0SBarry Smith   PetscFunctionBegin;
6739566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&fact->info));
6744ac6704cSBarry Smith   fact->factortype           = ftype;
6753d1c1ea0SBarry Smith   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
6763d1c1ea0SBarry Smith   fact->info.shiftamount     = 100.0*PETSC_MACHINE_EPSILON;
6773d1c1ea0SBarry Smith   fact->info.zeropivot       = 100.0*PETSC_MACHINE_EPSILON;
6783d1c1ea0SBarry Smith   fact->info.pivotinblocks   = 1.0;
6793d1c1ea0SBarry Smith   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
6803d1c1ea0SBarry Smith 
6819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor));
6829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor));
6839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor));
6849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor));
6859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor));
6869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor));
6879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",PCFactorGetMatSolverType_Factor));
6889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",PCFactorSetMatSolverType_Factor));
6899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",PCFactorSetUpMatSolverType_Factor));
6909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor));
6919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor));
6929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor));
6939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor));
6949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor));
6959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor));
6969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor));
6979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor));
6989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor));
6999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor));
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor));
7013d1c1ea0SBarry Smith   PetscFunctionReturn(0);
7023d1c1ea0SBarry Smith }
703*2e956fe4SStefano Zampini 
704*2e956fe4SStefano Zampini PetscErrorCode PCFactorClearComposedFunctions(PC pc)
705*2e956fe4SStefano Zampini {
706*2e956fe4SStefano Zampini   PetscFunctionBegin;
707*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",NULL));
708*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",NULL));
709*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",NULL));
710*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",NULL));
711*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",NULL));
712*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",NULL));
713*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",NULL));
714*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",NULL));
715*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",NULL));
716*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",NULL));
717*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",NULL));
718*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",NULL));
719*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",NULL));
720*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",NULL));
721*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",NULL));
722*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",NULL));
723*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",NULL));
724*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",NULL));
725*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",NULL));
726*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",NULL));
727*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",NULL));
728*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",NULL));
729*2e956fe4SStefano Zampini   PetscFunctionReturn(0);
730*2e956fe4SStefano Zampini }
731