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