xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision f42ad826bca5b91f2febf2bbee4d6ad7c168abd2)
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(PetscOptionItems *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   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
173   if (ilu->inplace) {
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((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
180       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)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 
188     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKERRQ(ierr);
189     if (((PC_Factor*)ilu)->info.errortype) { /* Factor() fails */
190       MatFactorInfo factinfo=((PC_Factor*)ilu)->info;
191       pc->failedreason = (PCFailedReason)factinfo.errortype;
192       PetscFunctionReturn(0);
193     }
194 
195     ((PC_Factor*)ilu)->fact = pc->pmat;
196     /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
197     ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr);
198   } else {
199     if (!pc->setupcalled) {
200       /* first time in so compute reordering and symbolic factorization */
201       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
202       if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
203       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
204       /*  Remove zeros along diagonal?     */
205       if (ilu->nonzerosalongdiagonal) {
206         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
207       }
208       if (!((PC_Factor*)ilu)->fact) {
209         ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
210       }
211       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
212       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
213 
214       ilu->actualfill = info.fill_ratio_needed;
215 
216       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
217     } else if (pc->flag != SAME_NONZERO_PATTERN) {
218       if (!ilu->reuseordering) {
219         /* compute a new ordering for the ILU */
220         ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);
221         ierr = ISDestroy(&ilu->col);CHKERRQ(ierr);
222         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
223         if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
224         if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
225         /*  Remove zeros along diagonal?     */
226         if (ilu->nonzerosalongdiagonal) {
227           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
228         }
229       }
230       ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
231       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
232       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
233       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
234 
235       ilu->actualfill = info.fill_ratio_needed;
236 
237       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
238     }
239     if (((PC_Factor*)ilu)->info.errortype) { /* FactorSymbolic() fails */
240       MatFactorInfo factinfo=((PC_Factor*)ilu)->info;
241       pc->failedreason = (PCFailedReason)factinfo.errortype;
242       PetscFunctionReturn(0);
243     }
244 
245     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
246     if (((PC_Factor*)ilu)->info.errortype) { /* FactorNumeric() fails */
247       MatFactorInfo factinfo=((PC_Factor*)ilu)->info;
248       pc->failedreason = (PCFailedReason)factinfo.errortype;
249     }
250   }
251   PetscFunctionReturn(0);
252 }
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "PCDestroy_ILU"
256 static PetscErrorCode PCDestroy_ILU(PC pc)
257 {
258   PC_ILU         *ilu = (PC_ILU*)pc->data;
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   ierr = PCReset_ILU(pc);CHKERRQ(ierr);
263   ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
264   ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
265   ierr = PetscFree(pc->data);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "PCApply_ILU"
271 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
272 {
273   PC_ILU         *ilu = (PC_ILU*)pc->data;
274   PetscErrorCode ierr;
275 
276   PetscFunctionBegin;
277   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "PCApplyTranspose_ILU"
283 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
284 {
285   PC_ILU         *ilu = (PC_ILU*)pc->data;
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293 #undef __FUNCT__
294 #define __FUNCT__ "PCApplySymmetricLeft_ILU"
295 static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
296 {
297   PetscErrorCode ierr;
298   PC_ILU         *icc = (PC_ILU*)pc->data;
299 
300   PetscFunctionBegin;
301   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 #undef __FUNCT__
306 #define __FUNCT__ "PCApplySymmetricRight_ILU"
307 static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
308 {
309   PetscErrorCode ierr;
310   PC_ILU         *icc = (PC_ILU*)pc->data;
311 
312   PetscFunctionBegin;
313   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 /*MC
318      PCILU - Incomplete factorization preconditioners.
319 
320    Options Database Keys:
321 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
322 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
323                       its factorization (overwrites original matrix)
324 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
325 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
326 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
327 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
328                                    this decreases the chance of getting a zero pivot
329 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
330 -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
331                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
332                              stability of the ILU factorization
333 
334    Level: beginner
335 
336   Concepts: incomplete factorization
337 
338    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
339 
340           For BAIJ matrices this implements a point block ILU
341 
342           The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
343           even when the matrix is not symmetric since the U stores the diagonals of the factorization.
344 
345           If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization
346           is never done on the GPU).
347 
348    References:
349    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
350    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
351 
352    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
353    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
354 
355    Review article:
356     APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
357       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
358       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
359       Science and Engineering, Kluwer, pp. 167--202.
360 
361 
362 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
363            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
364            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
365            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
366            PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace()
367 
368 M*/
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "PCCreate_ILU"
372 PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
373 {
374   PetscErrorCode ierr;
375   PC_ILU         *ilu;
376 
377   PetscFunctionBegin;
378   ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr);
379 
380   ((PC_Factor*)ilu)->fact               = 0;
381   ierr                                  = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
382   ((PC_Factor*)ilu)->factortype         = MAT_FACTOR_ILU;
383   ((PC_Factor*)ilu)->info.levels        = 0.;
384   ((PC_Factor*)ilu)->info.fill          = 1.0;
385   ilu->col                              = 0;
386   ilu->row                              = 0;
387   ilu->inplace                          = PETSC_FALSE;
388   ierr                                  = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
389   ierr                                  = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
390   ilu->reuseordering                    = PETSC_FALSE;
391   ((PC_Factor*)ilu)->info.dt            = PETSC_DEFAULT;
392   ((PC_Factor*)ilu)->info.dtcount       = PETSC_DEFAULT;
393   ((PC_Factor*)ilu)->info.dtcol         = PETSC_DEFAULT;
394   ((PC_Factor*)ilu)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
395   ((PC_Factor*)ilu)->info.shiftamount   = 100.0*PETSC_MACHINE_EPSILON;
396   ((PC_Factor*)ilu)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
397   ((PC_Factor*)ilu)->info.pivotinblocks = 1.0;
398   ilu->reusefill                        = PETSC_FALSE;
399   ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
400   pc->data                              = (void*)ilu;
401 
402   pc->ops->reset               = PCReset_ILU;
403   pc->ops->destroy             = PCDestroy_ILU;
404   pc->ops->apply               = PCApply_ILU;
405   pc->ops->applytranspose      = PCApplyTranspose_ILU;
406   pc->ops->setup               = PCSetUp_ILU;
407   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
408   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
409   pc->ops->view                = PCView_ILU;
410   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
411   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
412   pc->ops->applyrichardson     = 0;
413 
414   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
415   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
416   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
417   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
418   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
419   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
420   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
421   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
422   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
423   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
424   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
425   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
426   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
427   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
428   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ILU);CHKERRQ(ierr);
429   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
430   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);CHKERRQ(ierr);
431   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
432   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435