xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 6858538e710279fe46cd8279ab34c98b10293bbd)
1 
2 #include <../src/ksp/pc/impls/factor/factor.h>  /*I "petscpc.h" I*/
3 #include <petsc/private/matimpl.h>
4 
5 /*
6     If an ordering is not yet set and the matrix is available determine a default ordering
7 */
8 PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc)
9 {
10   PetscBool      foundmtype,flg;
11   const char     *prefix;
12 
13   PetscFunctionBegin;
14   if (pc->pmat) {
15     PetscCall(PCGetOptionsPrefix(pc,&prefix));
16     PetscCall(MatSetOptionsPrefixFactor(pc->pmat,prefix));
17     PC_Factor *fact = (PC_Factor*)pc->data;
18     PetscCall(MatSolverTypeGet(fact->solvertype,((PetscObject)pc->pmat)->type_name,fact->factortype,NULL,&foundmtype,NULL));
19     if (foundmtype) {
20       if (!fact->fact) {
21         PetscCall(MatGetFactor(pc->pmat,fact->solvertype,fact->factortype,&fact->fact));
22       } else if (!fact->fact->assembled) {
23         PetscCall(PetscStrcmp(fact->solvertype,fact->fact->solvertype,&flg));
24         if (!flg) {
25           Mat B;
26           PetscCall(MatGetFactor(pc->pmat,fact->solvertype,fact->factortype,&B));
27           PetscCall(MatHeaderReplace(fact->fact,&B));
28         }
29       }
30       if (!fact->ordering) {
31         PetscBool       canuseordering;
32         MatOrderingType otype;
33 
34         PetscCall(MatFactorGetCanUseOrdering(fact->fact,&canuseordering));
35         if (canuseordering) {
36           PetscCall(MatFactorGetPreferredOrdering(fact->fact,fact->factortype,&otype));
37         } else otype = MATORDERINGEXTERNAL;
38         PetscCall(PetscStrallocpy(otype,(char **)&fact->ordering));
39       }
40     }
41   }
42   PetscFunctionReturn(0);
43 }
44 
45 static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag)
46 {
47   PC_Factor *lu = (PC_Factor*)pc->data;
48 
49   PetscFunctionBegin;
50   lu->reuseordering = flag;
51   PetscFunctionReturn(0);
52 }
53 
54 static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc,PetscBool flag)
55 {
56   PC_Factor *lu = (PC_Factor*)pc->data;
57 
58   PetscFunctionBegin;
59   lu->reusefill = flag;
60   PetscFunctionReturn(0);
61 }
62 
63 static PetscErrorCode  PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg)
64 {
65   PC_Factor *dir = (PC_Factor*)pc->data;
66 
67   PetscFunctionBegin;
68   dir->inplace = flg;
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode  PCFactorGetUseInPlace_Factor(PC pc,PetscBool *flg)
73 {
74   PC_Factor *dir = (PC_Factor*)pc->data;
75 
76   PetscFunctionBegin;
77   *flg = dir->inplace;
78   PetscFunctionReturn(0);
79 }
80 
81 /*@
82     PCFactorSetUpMatSolverType - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
83        set the options for that particular factorization object.
84 
85   Input Parameter:
86 .  pc  - the preconditioner context
87 
88   Notes:
89     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.
90 
91   Level: intermediate
92 
93 .seealso: `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()`
94 @*/
95 PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
96 {
97   PetscFunctionBegin;
98   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
99   PetscTryMethod(pc,"PCFactorSetUpMatSolverType_C",(PC),(pc));
100   PetscFunctionReturn(0);
101 }
102 
103 /*@
104    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
105 
106    Logically Collective on PC
107 
108    Input Parameters:
109 +  pc - the preconditioner context
110 -  zero - all pivots smaller than this will be considered zero
111 
112    Options Database Key:
113 .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
114 
115    Level: intermediate
116 
117 .seealso: `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
118 @*/
119 PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
120 {
121   PetscFunctionBegin;
122   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123   PetscValidLogicalCollectiveReal(pc,zero,2);
124   PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));
125   PetscFunctionReturn(0);
126 }
127 
128 /*@
129    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
130      numerical factorization, thus the matrix has nonzero pivots
131 
132    Logically Collective on PC
133 
134    Input Parameters:
135 +  pc - the preconditioner context
136 -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
137 
138    Options Database Key:
139 .  -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
140 
141    Level: intermediate
142 
143 .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()`
144 @*/
145 PetscErrorCode  PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
146 {
147   PetscFunctionBegin;
148   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
149   PetscValidLogicalCollectiveEnum(pc,shifttype,2);
150   PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));
151   PetscFunctionReturn(0);
152 }
153 
154 /*@
155    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
156      numerical factorization, thus the matrix has nonzero pivots
157 
158    Logically Collective on PC
159 
160    Input Parameters:
161 +  pc - the preconditioner context
162 -  shiftamount - amount of shift
163 
164    Options Database Key:
165 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
166 
167    Level: intermediate
168 
169 .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`
170 @*/
171 PetscErrorCode  PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
172 {
173   PetscFunctionBegin;
174   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
175   PetscValidLogicalCollectiveReal(pc,shiftamount,2);
176   PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));
177   PetscFunctionReturn(0);
178 }
179 
180 /*@
181    PCFactorSetDropTolerance - The preconditioner will use an ILU
182    based on a drop tolerance. (Under development)
183 
184    Logically Collective on PC
185 
186    Input Parameters:
187 +  pc - the preconditioner context
188 .  dt - the drop tolerance, try from 1.e-10 to .1
189 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
190 -  maxrowcount - the max number of nonzeros allowed in a row, best value
191                  depends on the number of nonzeros in row of original matrix
192 
193    Options Database Key:
194 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
195 
196    Level: intermediate
197 
198       There are NO default values for the 3 parameters, you must set them with reasonable values for your
199       matrix. We don't know how to compute reasonable values.
200 
201 @*/
202 PetscErrorCode  PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
203 {
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
206   PetscValidLogicalCollectiveReal(pc,dtcol,3);
207   PetscValidLogicalCollectiveInt(pc,maxrowcount,4);
208   PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));
209   PetscFunctionReturn(0);
210 }
211 
212 /*@
213    PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
214 
215    Not Collective
216 
217    Input Parameters:
218 .  pc - the preconditioner context
219 
220    Output Parameter:
221 .  pivot - the tolerance
222 
223    Level: intermediate
224 
225 .seealso: `PCFactorSetZeroPivot()`
226 @*/
227 PetscErrorCode  PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
228 {
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
231   PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));
232   PetscFunctionReturn(0);
233 }
234 
235 /*@
236    PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
237 
238    Not Collective
239 
240    Input Parameters:
241 .  pc - the preconditioner context
242 
243    Output Parameter:
244 .  shift - how much to shift the diagonal entry
245 
246    Level: intermediate
247 
248 .seealso: `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()`
249 @*/
250 PetscErrorCode  PCFactorGetShiftAmount(PC pc,PetscReal *shift)
251 {
252   PetscFunctionBegin;
253   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
254   PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));
255   PetscFunctionReturn(0);
256 }
257 
258 /*@
259    PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
260 
261    Not Collective
262 
263    Input Parameters:
264 .  pc - the preconditioner context
265 
266    Output Parameter:
267 .  type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS
268 
269    Level: intermediate
270 
271 .seealso: `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()`
272 @*/
273 PetscErrorCode  PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
274 {
275   PetscFunctionBegin;
276   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
277   PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));
278   PetscFunctionReturn(0);
279 }
280 
281 /*@
282    PCFactorGetLevels - Gets the number of levels of fill to use.
283 
284    Logically Collective on PC
285 
286    Input Parameters:
287 .  pc - the preconditioner context
288 
289    Output Parameter:
290 .  levels - number of levels of fill
291 
292    Level: intermediate
293 
294 @*/
295 PetscErrorCode  PCFactorGetLevels(PC pc,PetscInt *levels)
296 {
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
299   PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
300   PetscFunctionReturn(0);
301 }
302 
303 /*@
304    PCFactorSetLevels - Sets the number of levels of fill to use.
305 
306    Logically Collective on PC
307 
308    Input Parameters:
309 +  pc - the preconditioner context
310 -  levels - number of levels of fill
311 
312    Options Database Key:
313 .  -pc_factor_levels <levels> - Sets fill level
314 
315    Level: intermediate
316 
317 @*/
318 PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
319 {
320   PetscFunctionBegin;
321   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
322   PetscCheck(levels >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
323   PetscValidLogicalCollectiveInt(pc,levels,2);
324   PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
325   PetscFunctionReturn(0);
326 }
327 
328 /*@
329    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
330    treated as level 0 fill even if there is no non-zero location.
331 
332    Logically Collective on PC
333 
334    Input Parameters:
335 +  pc - the preconditioner context
336 -  flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
337 
338    Options Database Key:
339 .  -pc_factor_diagonal_fill <bool> - allow the diagonal fill
340 
341    Notes:
342    Does not apply with 0 fill.
343 
344    Level: intermediate
345 
346 .seealso: `PCFactorGetAllowDiagonalFill()`
347 @*/
348 PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
349 {
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
352   PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
353   PetscFunctionReturn(0);
354 }
355 
356 /*@
357    PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
358        treated as level 0 fill even if there is no non-zero location.
359 
360    Logically Collective on PC
361 
362    Input Parameter:
363 .  pc - the preconditioner context
364 
365    Output Parameter:
366 .   flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
367 
368    Notes:
369    Does not apply with 0 fill.
370 
371    Level: intermediate
372 
373 .seealso: `PCFactorSetAllowDiagonalFill()`
374 @*/
375 PetscErrorCode  PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
376 {
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
379   PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
380   PetscFunctionReturn(0);
381 }
382 
383 /*@
384    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
385 
386    Logically Collective on PC
387 
388    Input Parameters:
389 +  pc - the preconditioner context
390 -  tol - diagonal entries smaller than this in absolute value are considered zero
391 
392    Options Database Key:
393 .  -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance
394 
395    Level: intermediate
396 
397 .seealso: `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()`
398 @*/
399 PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
400 {
401   PetscFunctionBegin;
402   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
403   PetscValidLogicalCollectiveReal(pc,rtol,2);
404   PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
405   PetscFunctionReturn(0);
406 }
407 
408 /*@C
409    PCFactorSetMatSolverType - sets the software that is used to perform the factorization
410 
411    Logically Collective on PC
412 
413    Input Parameters:
414 +  pc - the preconditioner context
415 -  stype - for example, superlu, superlu_dist
416 
417    Options Database Key:
418 .  -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
419 
420    Level: intermediate
421 
422    Note:
423      By default this will use the PETSc factorization if it exists
424 
425 .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`
426 @*/
427 PetscErrorCode  PCFactorSetMatSolverType(PC pc,MatSolverType stype)
428 {
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
431   PetscTryMethod(pc,"PCFactorSetMatSolverType_C",(PC,MatSolverType),(pc,stype));
432   PetscFunctionReturn(0);
433 }
434 
435 /*@C
436    PCFactorGetMatSolverType - gets the software that is used to perform the factorization
437 
438    Not Collective
439 
440    Input Parameter:
441 .  pc - the preconditioner context
442 
443    Output Parameter:
444 .   stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
445 
446    Level: intermediate
447 
448 .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`
449 @*/
450 PetscErrorCode  PCFactorGetMatSolverType(PC pc,MatSolverType *stype)
451 {
452   PetscErrorCode (*f)(PC,MatSolverType*);
453 
454   PetscFunctionBegin;
455   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
456   PetscValidPointer(stype,2);
457   PetscCall(PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",&f));
458   if (f) PetscCall((*f)(pc,stype));
459   else *stype = NULL;
460   PetscFunctionReturn(0);
461 }
462 
463 /*@
464    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
465    fill = number nonzeros in factor/number nonzeros in original matrix.
466 
467    Not Collective, each process can expect a different amount of fill
468 
469    Input Parameters:
470 +  pc - the preconditioner context
471 -  fill - amount of expected fill
472 
473    Options Database Key:
474 .  -pc_factor_fill <fill> - Sets fill amount
475 
476    Level: intermediate
477 
478    Note:
479    For sparse matrix factorizations it is difficult to predict how much
480    fill to expect. By running with the option -info PETSc will print the
481    actual amount of fill used; allowing you to set the value accurately for
482    future runs. Default PETSc uses a value of 5.0
483 
484    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
485 
486 @*/
487 PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
488 {
489   PetscFunctionBegin;
490   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
491   PetscCheck(fill >= 1.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
492   PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
493   PetscFunctionReturn(0);
494 }
495 
496 /*@
497    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
498    For dense matrices, this enables the solution of much larger problems.
499    For sparse matrices the factorization cannot be done truly in-place
500    so this does not save memory during the factorization, but after the matrix
501    is factored, the original unfactored matrix is freed, thus recovering that
502    space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
503 
504    Logically Collective on PC
505 
506    Input Parameters:
507 +  pc - the preconditioner context
508 -  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
509 
510    Options Database Key:
511 .  -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
512 
513    Notes:
514    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
515    a different matrix is provided for the multiply and the preconditioner in
516    a call to KSPSetOperators().
517    This is because the Krylov space methods require an application of the
518    matrix multiplication, which is not possible here because the matrix has
519    been factored in-place, replacing the original matrix.
520 
521    Level: intermediate
522 
523 .seealso: `PCFactorGetUseInPlace()`
524 @*/
525 PetscErrorCode  PCFactorSetUseInPlace(PC pc,PetscBool flg)
526 {
527   PetscFunctionBegin;
528   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
529   PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
530   PetscFunctionReturn(0);
531 }
532 
533 /*@
534    PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
535 
536    Logically Collective on PC
537 
538    Input Parameter:
539 .  pc - the preconditioner context
540 
541    Output Parameter:
542 .  flg - PETSC_TRUE to enable, PETSC_FALSE to disable
543 
544    Level: intermediate
545 
546 .seealso: `PCFactorSetUseInPlace()`
547 @*/
548 PetscErrorCode  PCFactorGetUseInPlace(PC pc,PetscBool *flg)
549 {
550   PetscFunctionBegin;
551   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
552   PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
553   PetscFunctionReturn(0);
554 }
555 
556 /*@C
557     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
558     be used in the LU, ILU, Cholesky, and ICC factorizations.
559 
560     Logically Collective on PC
561 
562     Input Parameters:
563 +   pc - the preconditioner context
564 -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
565 
566     Options Database Key:
567 .   -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine
568 
569     Level: intermediate
570 
571     Notes:
572       Nested dissection is used by default for some of PETSc's sparse matrix formats
573 
574      For Cholesky and ICC and the SBAIJ format the only reordering available is natural since only the upper half of the matrix is stored
575      and reordering this matrix is very expensive.
576 
577       You can use a SeqAIJ matrix with Cholesky and ICC and use any ordering.
578 
579       MATORDERINGEXTERNAL means PETSc will not compute an ordering and the package will use its own ordering, usable with MATSOLVERCHOLMOD, MATSOLVERUMFPACK, and others.
580 
581 .seealso: `MatOrderingType`
582 
583 @*/
584 PetscErrorCode  PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
585 {
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
588   PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));
589   PetscFunctionReturn(0);
590 }
591 
592 /*@
593     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
594       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
595       it is never done. For the MATLAB and SuperLU factorization this is used.
596 
597     Logically Collective on PC
598 
599     Input Parameters:
600 +   pc - the preconditioner context
601 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
602 
603     Options Database Key:
604 .   -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance
605 
606     Level: intermediate
607 
608 .seealso: `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()`
609 @*/
610 PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
611 {
612   PetscFunctionBegin;
613   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
614   PetscValidLogicalCollectiveReal(pc,dtcol,2);
615   PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
616   PetscFunctionReturn(0);
617 }
618 
619 /*@
620     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
621       with BAIJ or SBAIJ matrices
622 
623     Logically Collective on PC
624 
625     Input Parameters:
626 +   pc - the preconditioner context
627 -   pivot - PETSC_TRUE or PETSC_FALSE
628 
629     Options Database Key:
630 .   -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for BAIJ and SBAIJ
631 
632     Level: intermediate
633 
634 .seealso: `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()`
635 @*/
636 PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
637 {
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
640   PetscValidLogicalCollectiveBool(pc,pivot,2);
641   PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
642   PetscFunctionReturn(0);
643 }
644 
645 /*@
646    PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
647    this causes later ones to use the fill ratio computed in the initial factorization.
648 
649    Logically Collective on PC
650 
651    Input Parameters:
652 +  pc - the preconditioner context
653 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
654 
655    Options Database Key:
656 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
657 
658    Level: intermediate
659 
660 .seealso: `PCFactorSetReuseOrdering()`
661 @*/
662 PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool flag)
663 {
664   PetscFunctionBegin;
665   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
666   PetscValidLogicalCollectiveBool(pc,flag,2);
667   PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
668   PetscFunctionReturn(0);
669 }
670 
671 PetscErrorCode PCFactorInitialize(PC pc,MatFactorType ftype)
672 {
673   PC_Factor      *fact = (PC_Factor*)pc->data;
674 
675   PetscFunctionBegin;
676   PetscCall(MatFactorInfoInitialize(&fact->info));
677   fact->factortype           = ftype;
678   fact->info.shifttype       = (PetscReal)MAT_SHIFT_NONE;
679   fact->info.shiftamount     = 100.0*PETSC_MACHINE_EPSILON;
680   fact->info.zeropivot       = 100.0*PETSC_MACHINE_EPSILON;
681   fact->info.pivotinblocks   = 1.0;
682   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
683 
684   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor));
685   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor));
686   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor));
687   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor));
688   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor));
689   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor));
690   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",PCFactorGetMatSolverType_Factor));
691   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",PCFactorSetMatSolverType_Factor));
692   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",PCFactorSetUpMatSolverType_Factor));
693   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor));
694   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor));
695   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor));
696   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor));
697   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor));
698   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor));
699   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor));
700   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor));
701   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor));
702   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor));
703   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor));
704   PetscFunctionReturn(0);
705 }
706 
707 PetscErrorCode PCFactorClearComposedFunctions(PC pc)
708 {
709   PetscFunctionBegin;
710   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",NULL));
711   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",NULL));
712   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",NULL));
713   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",NULL));
714   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",NULL));
715   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",NULL));
716   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",NULL));
717   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",NULL));
718   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",NULL));
719   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",NULL));
720   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",NULL));
721   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",NULL));
722   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",NULL));
723   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",NULL));
724   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",NULL));
725   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",NULL));
726   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",NULL));
727   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",NULL));
728   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",NULL));
729   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",NULL));
730   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",NULL));
731   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",NULL));
732   PetscFunctionReturn(0);
733 }
734