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