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