xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 35e7444da1e6edf29ca585d75dc3fe3d2c63a6e4)
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     if (((PC_Factor*)ilu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
163     if (((PC_Factor*)ilu)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
164     if (((PC_Factor*)ilu)->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);}
165     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
166     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
167     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
168     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
169     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
170     if (((PC_Factor*)ilu)->fact) {
171       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr);
172       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
173       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
174       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
175       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
176       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
177       ierr = MatView(((PC_Factor*)ilu)->fact,viewer);CHKERRQ(ierr);
178       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
179       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
180       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
181       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
182     }
183   } else if (isstring) {
184     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
185   } else {
186     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
187   }
188   PetscFunctionReturn(0);
189 }
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "PCSetUp_ILU"
193 static PetscErrorCode PCSetUp_ILU(PC pc)
194 {
195   PetscErrorCode ierr;
196   PC_ILU         *ilu = (PC_ILU*)pc->data;
197   MatInfo        info;
198 
199   PetscFunctionBegin;
200   if (ilu->inplace) {
201     CHKMEMQ;
202     if (!pc->setupcalled) {
203 
204       /* In-place factorization only makes sense with the natural ordering,
205          so we only need to get the ordering once, even if nonzero structure changes */
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     }
210 
211     /* In place ILU only makes sense with fill factor of 1.0 because
212        cannot have levels of fill */
213     ((PC_Factor*)ilu)->info.fill          = 1.0;
214     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
215     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKMEMQ;
216     ((PC_Factor*)ilu)->fact = pc->pmat;
217   } else {
218     if (!pc->setupcalled) {
219       /* first time in so compute reordering and symbolic factorization */
220       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
221       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
222       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
223       /*  Remove zeros along diagonal?     */
224       if (ilu->nonzerosalongdiagonal) {
225         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
226       }
227     CHKMEMQ;
228       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
229     CHKMEMQ;
230       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
231       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
232       ilu->actualfill = info.fill_ratio_needed;
233       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
234     } else if (pc->flag != SAME_NONZERO_PATTERN) {
235       if (!ilu->reuseordering) {
236         /* compute a new ordering for the ILU */
237         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
238         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
239         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
240         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
241         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
242         /*  Remove zeros along diagonal?     */
243         if (ilu->nonzerosalongdiagonal) {
244           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
245         }
246       }
247       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
248       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
249       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
250       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
251       ilu->actualfill = info.fill_ratio_needed;
252       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
253     }
254     CHKMEMQ;
255     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
256     CHKMEMQ;
257   }
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "PCDestroy_ILU"
263 static PetscErrorCode PCDestroy_ILU(PC pc)
264 {
265   PC_ILU         *ilu = (PC_ILU*)pc->data;
266   PetscErrorCode ierr;
267 
268   PetscFunctionBegin;
269   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
270   ierr = PetscStrfree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
271   ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
272   ierr = PetscFree(ilu);CHKERRQ(ierr);
273   PetscFunctionReturn(0);
274 }
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "PCApply_ILU"
278 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
279 {
280   PC_ILU         *ilu = (PC_ILU*)pc->data;
281   PetscErrorCode ierr;
282 
283   PetscFunctionBegin;
284   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "PCApplyTranspose_ILU"
290 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
291 {
292   PC_ILU         *ilu = (PC_ILU*)pc->data;
293   PetscErrorCode ierr;
294 
295   PetscFunctionBegin;
296   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 /*MC
301      PCILU - Incomplete factorization preconditioners.
302 
303    Options Database Keys:
304 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
305 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
306                       its factorization (overwrites original matrix)
307 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
308 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
309 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - use drop tolerance factorization
310 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
311 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
312                                    this decreases the chance of getting a zero pivot
313 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
314 .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
315                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
316                              stability of the ILU factorization
317 .  -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization
318 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
319 -  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
320    is optional with true being the default
321 
322    Level: beginner
323 
324   Concepts: incomplete factorization
325 
326    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
327 
328           For BAIJ matrices this implements a point block ILU
329 
330    References:
331    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
332    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
333 
334    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
335    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
336 
337    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
338       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
339       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
340       Science and Engineering, Kluwer, pp. 167--202.
341 
342 
343 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
344            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks(),
345            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
346            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks()
347 
348 M*/
349 
350 EXTERN_C_BEGIN
351 #undef __FUNCT__
352 #define __FUNCT__ "PCCreate_ILU"
353 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
354 {
355   PetscErrorCode ierr;
356   PC_ILU         *ilu;
357 
358   PetscFunctionBegin;
359   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
360 
361   ((PC_Factor*)ilu)->fact                    = 0;
362   ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
363   ((PC_Factor*)ilu)->info.levels             = 0.;
364   ((PC_Factor*)ilu)->info.fill               = 1.0;
365   ilu->col                     = 0;
366   ilu->row                     = 0;
367   ilu->inplace                 = PETSC_FALSE;
368   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
369   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
370   ilu->reuseordering           = PETSC_FALSE;
371   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
372   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
373   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
374   /* Only one of shiftnz, shiftpd and shiftinblocks is allowed to be set as non-zero.
375      Setting shiftnz=TURE as default causes confusion if user sets shiftpd or shiftinblocks as true */
376   ((PC_Factor*)ilu)->info.shiftnz            = 0.0;    /* false */
377   ((PC_Factor*)ilu)->info.shiftpd            = 0.0;    /* false */
378   ((PC_Factor*)ilu)->info.shiftinblocks      = 0.0;    /* false */
379   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
380   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
381   ilu->reusefill               = PETSC_FALSE;
382   ((PC_Factor*)ilu)->info.diagonal_fill      = 0.;
383   pc->data                     = (void*)ilu;
384 
385   pc->ops->destroy             = PCDestroy_ILU;
386   pc->ops->apply               = PCApply_ILU;
387   pc->ops->applytranspose      = PCApplyTranspose_ILU;
388   pc->ops->setup               = PCSetUp_ILU;
389   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
390   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
391   pc->ops->view                = PCView_ILU;
392   pc->ops->applyrichardson     = 0;
393 
394   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
395                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
396   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
397                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
398   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
399                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
400   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
401                     PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr);
402 
403   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
404                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
405   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
406                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
407   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
408                     PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
409   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
410                     PCFactorSetFill_Factor);CHKERRQ(ierr);
411   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
412                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
413   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
414                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
415   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
416                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
417   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
418                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
419   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
420                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
421   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
422                     PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
423   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
424                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
425   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
426                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 EXTERN_C_END
430