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