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