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