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