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