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