xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 3463a7ba26fa12b272fa42291bbeadb1398fb231)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a ILU factorization preconditioner for any Mat implementation
5 */
6 #include "../src/ksp/pc/impls/factor/ilu/ilu.h"     /*I "petscpc.h"  I*/
7 
8 /* ------------------------------------------------------------------------------------------*/
9 EXTERN_C_BEGIN
10 #undef __FUNCT__
11 #define __FUNCT__ "PCFactorSetReuseFill_ILU"
12 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag)
13 {
14   PC_ILU *lu = (PC_ILU*)pc->data;
15 
16   PetscFunctionBegin;
17   lu->reusefill = flag;
18   PetscFunctionReturn(0);
19 }
20 EXTERN_C_END
21 
22 EXTERN_C_BEGIN
23 #undef __FUNCT__
24 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU"
25 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
26 {
27   PC_ILU *ilu = (PC_ILU*)pc->data;
28 
29   PetscFunctionBegin;
30   ilu->nonzerosalongdiagonal = PETSC_TRUE;
31   if (z == PETSC_DECIDE) {
32     ilu->nonzerosalongdiagonaltol = 1.e-10;
33   } else {
34     ilu->nonzerosalongdiagonaltol = z;
35   }
36   PetscFunctionReturn(0);
37 }
38 EXTERN_C_END
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "PCDestroy_ILU_Internal"
42 PetscErrorCode PCDestroy_ILU_Internal(PC pc)
43 {
44   PC_ILU         *ilu = (PC_ILU*)pc->data;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);}
49   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
50   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
51   PetscFunctionReturn(0);
52 }
53 
54 EXTERN_C_BEGIN
55 #undef __FUNCT__
56 #define __FUNCT__ "PCFactorSetDropTolerance_ILU"
57 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
58 {
59   PC_ILU         *ilu = (PC_ILU*)pc->data;
60 
61   PetscFunctionBegin;
62   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
63     SETERRQ(PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
64   }
65   ((PC_Factor*)ilu)->info.dt      = dt;
66   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
67   ((PC_Factor*)ilu)->info.dtcount = dtcount;
68   ((PC_Factor*)ilu)->info.usedt   = 1.0;
69   PetscFunctionReturn(0);
70 }
71 EXTERN_C_END
72 
73 EXTERN_C_BEGIN
74 #undef __FUNCT__
75 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
76 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
77 {
78   PC_ILU *ilu = (PC_ILU*)pc->data;
79 
80   PetscFunctionBegin;
81   ilu->reuseordering = flag;
82   PetscFunctionReturn(0);
83 }
84 EXTERN_C_END
85 
86 EXTERN_C_BEGIN
87 #undef __FUNCT__
88 #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
89 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc)
90 {
91   PC_ILU *dir = (PC_ILU*)pc->data;
92 
93   PetscFunctionBegin;
94   dir->inplace = PETSC_TRUE;
95   PetscFunctionReturn(0);
96 }
97 EXTERN_C_END
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "PCSetFromOptions_ILU"
101 static PetscErrorCode PCSetFromOptions_ILU(PC pc)
102 {
103   PetscErrorCode ierr;
104   PetscInt       dtmax = 3,itmp;
105   PetscTruth     flg;
106   PetscReal      dt[3];
107   PC_ILU         *ilu = (PC_ILU*)pc->data;
108   PetscReal      tol;
109 
110   PetscFunctionBegin;
111   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
112     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
113 
114     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
115     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
116     flg  = PETSC_FALSE;
117     ierr = PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
118     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
119 
120     dt[0] = ((PC_Factor*)ilu)->info.dt;
121     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
122     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
123     ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
124     if (flg) {
125       ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
126     }
127 
128     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
129     if (flg) {
130       tol = PETSC_DECIDE;
131       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
132       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
133     }
134 
135   ierr = PetscOptionsTail();CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "PCView_ILU"
141 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
142 {
143   PC_ILU         *ilu = (PC_ILU*)pc->data;
144   PetscErrorCode ierr;
145   PetscTruth     iascii;
146 
147   PetscFunctionBegin;
148   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
149   if (iascii) {
150     if (ilu->inplace) {
151       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: in-place factorization\n");CHKERRQ(ierr);
152     } else {
153       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: out-of-place factorization\n");CHKERRQ(ierr);
154     }
155 
156     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing fill from past factorization\n");CHKERRQ(ierr);}
157     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing reordering from past factorization\n");CHKERRQ(ierr);}
158   }
159   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "PCSetUp_ILU"
165 static PetscErrorCode PCSetUp_ILU(PC pc)
166 {
167   PetscErrorCode ierr;
168   PC_ILU         *ilu = (PC_ILU*)pc->data;
169   MatInfo        info;
170 
171   PetscFunctionBegin;
172   if (ilu->inplace) {
173     CHKMEMQ;
174     if (!pc->setupcalled) {
175 
176       /* In-place factorization only makes sense with the natural ordering,
177          so we only need to get the ordering once, even if nonzero structure changes */
178       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
179       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
180       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
181     }
182 
183     /* In place ILU only makes sense with fill factor of 1.0 because
184        cannot have levels of fill */
185     ((PC_Factor*)ilu)->info.fill          = 1.0;
186     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
187     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKMEMQ;
188     ((PC_Factor*)ilu)->fact = pc->pmat;
189   } else {
190     if (!pc->setupcalled) {
191       /* first time in so compute reordering and symbolic factorization */
192       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
193       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
194       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
195       /*  Remove zeros along diagonal?     */
196       if (ilu->nonzerosalongdiagonal) {
197         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
198       }
199     CHKMEMQ;
200       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
201     CHKMEMQ;
202       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
203       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
204       ilu->actualfill = info.fill_ratio_needed;
205       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
206     } else if (pc->flag != SAME_NONZERO_PATTERN) {
207       if (!ilu->reuseordering) {
208         /* compute a new ordering for the ILU */
209         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
210         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
211         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
212         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
213         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
214         /*  Remove zeros along diagonal?     */
215         if (ilu->nonzerosalongdiagonal) {
216           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
217         }
218       }
219       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
220       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
221       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
222       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
223       ilu->actualfill = info.fill_ratio_needed;
224       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
225     }
226     CHKMEMQ;
227     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
228     CHKMEMQ;
229   }
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "PCDestroy_ILU"
235 static PetscErrorCode PCDestroy_ILU(PC pc)
236 {
237   PC_ILU         *ilu = (PC_ILU*)pc->data;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
242   ierr = PetscStrfree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
243   ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
244   ierr = PetscFree(ilu);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "PCApply_ILU"
250 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
251 {
252   PC_ILU         *ilu = (PC_ILU*)pc->data;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "PCApplyTranspose_ILU"
262 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
263 {
264   PC_ILU         *ilu = (PC_ILU*)pc->data;
265   PetscErrorCode ierr;
266 
267   PetscFunctionBegin;
268   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 
272 /*MC
273      PCILU - Incomplete factorization preconditioners.
274 
275    Options Database Keys:
276 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
277 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
278                       its factorization (overwrites original matrix)
279 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
280 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
281 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - use drop tolerance factorization
282 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
283 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
284                                    this decreases the chance of getting a zero pivot
285 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
286 .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
287                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
288                              stability of the ILU factorization
289 
290    Level: beginner
291 
292   Concepts: incomplete factorization
293 
294    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
295 
296           For BAIJ matrices this implements a point block ILU
297 
298    References:
299    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
300    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
301 
302    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
303    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
304 
305    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
306       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
307       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
308       Science and Engineering, Kluwer, pp. 167--202.
309 
310 
311 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
312            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
313            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
314            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks()
315 
316 M*/
317 
318 EXTERN_C_BEGIN
319 #undef __FUNCT__
320 #define __FUNCT__ "PCCreate_ILU"
321 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
322 {
323   PetscErrorCode ierr;
324   PC_ILU         *ilu;
325 
326   PetscFunctionBegin;
327   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
328 
329   ((PC_Factor*)ilu)->fact                    = 0;
330   ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
331   ((PC_Factor*)ilu)->info.levels             = 0.;
332   ((PC_Factor*)ilu)->info.fill               = 1.0;
333   ilu->col                     = 0;
334   ilu->row                     = 0;
335   ilu->inplace                 = PETSC_FALSE;
336   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
337   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
338   ilu->reuseordering           = PETSC_FALSE;
339   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
340   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
341   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
342   ((PC_Factor*)ilu)->info.shifttype          = MAT_SHIFT_NONZERO;
343   ((PC_Factor*)ilu)->info.shiftamount        = 1.e-12;
344   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
345   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
346   ilu->reusefill               = PETSC_FALSE;
347   ((PC_Factor*)ilu)->info.diagonal_fill      = 0.;
348   pc->data                     = (void*)ilu;
349 
350   pc->ops->destroy             = PCDestroy_ILU;
351   pc->ops->apply               = PCApply_ILU;
352   pc->ops->applytranspose      = PCApplyTranspose_ILU;
353   pc->ops->setup               = PCSetUp_ILU;
354   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
355   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
356   pc->ops->view                = PCView_ILU;
357   pc->ops->applyrichardson     = 0;
358 
359   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
360                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
361   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
362                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
363   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
364                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
365   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
366                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
367   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
368                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
369   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
370                     PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
371   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
372                     PCFactorSetFill_Factor);CHKERRQ(ierr);
373   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
374                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
375   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
376                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
377   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
378                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
379   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
380                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
381   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
382                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
383   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
384                     PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
385   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
386                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
387   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
388                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 EXTERN_C_END
392