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