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