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