xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision ff00295018bcf3807910866b6079563e62bacbd4)
1 /*
2    Defines a ILU factorization preconditioner for any Mat implementation
3 */
4 #include "src/ksp/pc/pcimpl.h"                 /*I "petscpc.h"  I*/
5 #include "src/ksp/pc/impls/factor/ilu/ilu.h"
6 #include "src/mat/matimpl.h"
7 
8 /* ------------------------------------------------------------------------------------------*/
9 
10 EXTERN_C_BEGIN
11 #undef __FUNCT__
12 #define __FUNCT__ "PCILUReorderForNonzeroDiagonal_ILU"
13 PetscErrorCode PCILUReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
14 {
15   PC_ILU *ilu = (PC_ILU*)pc->data;
16 
17   PetscFunctionBegin;
18   ilu->nonzerosalongdiagonal = PETSC_TRUE;
19   if (z == PETSC_DECIDE) {
20     ilu->nonzerosalongdiagonaltol = 1.e-10;
21   } else {
22     ilu->nonzerosalongdiagonaltol = z;
23   }
24   PetscFunctionReturn(0);
25 }
26 EXTERN_C_END
27 
28 EXTERN_C_BEGIN
29 #undef __FUNCT__
30 #define __FUNCT__ "PCILUSetSetZeroPivot_ILU"
31 PetscErrorCode PCILUSetZeroPivot_ILU(PC pc,PetscReal z)
32 {
33   PC_ILU *lu;
34 
35   PetscFunctionBegin;
36   lu                 = (PC_ILU*)pc->data;
37   lu->info.zeropivot = z;
38   PetscFunctionReturn(0);
39 }
40 EXTERN_C_END
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "PCDestroy_ILU_Internal"
44 PetscErrorCode PCDestroy_ILU_Internal(PC pc)
45 {
46   PC_ILU         *ilu = (PC_ILU*)pc->data;
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);}
51   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
52   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
53   PetscFunctionReturn(0);
54 }
55 
56 EXTERN_C_BEGIN
57 #undef __FUNCT__
58 #define __FUNCT__ "PCILUSetDamping_ILU"
59 PetscErrorCode PCILUSetDamping_ILU(PC pc,PetscReal damping)
60 {
61   PC_ILU *dir;
62 
63   PetscFunctionBegin;
64   dir = (PC_ILU*)pc->data;
65   if (damping == (PetscReal) PETSC_DECIDE) {
66     dir->info.damping = 1.e-12;
67   } else {
68     dir->info.damping = damping;
69   }
70   PetscFunctionReturn(0);
71 }
72 EXTERN_C_END
73 
74 EXTERN_C_BEGIN
75 #undef __FUNCT__
76 #define __FUNCT__ "PCILUSetShift_ILU"
77 PetscErrorCode PCILUSetShift_ILU(PC pc,PetscTruth shift)
78 {
79   PC_ILU *dir;
80 
81   PetscFunctionBegin;
82   dir = (PC_ILU*)pc->data;
83   dir->info.shift = shift;
84   if (shift) dir->info.shift_fraction = 0.0;
85   PetscFunctionReturn(0);
86 }
87 EXTERN_C_END
88 
89 EXTERN_C_BEGIN
90 #undef __FUNCT__
91 #define __FUNCT__ "PCILUSetUseDropTolerance_ILU"
92 PetscErrorCode PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
93 {
94   PC_ILU         *ilu;
95   PetscErrorCode ierr;
96 
97   PetscFunctionBegin;
98   ilu = (PC_ILU*)pc->data;
99   if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) {
100     pc->setupcalled   = 0;
101     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
102   }
103   ilu->usedt        = PETSC_TRUE;
104   ilu->info.dt      = dt;
105   ilu->info.dtcol   = dtcol;
106   ilu->info.dtcount = dtcount;
107   ilu->info.fill    = PETSC_DEFAULT;
108   PetscFunctionReturn(0);
109 }
110 EXTERN_C_END
111 
112 EXTERN_C_BEGIN
113 #undef __FUNCT__
114 #define __FUNCT__ "PCILUSetFill_ILU"
115 PetscErrorCode PCILUSetFill_ILU(PC pc,PetscReal fill)
116 {
117   PC_ILU *dir;
118 
119   PetscFunctionBegin;
120   dir            = (PC_ILU*)pc->data;
121   dir->info.fill = fill;
122   PetscFunctionReturn(0);
123 }
124 EXTERN_C_END
125 
126 EXTERN_C_BEGIN
127 #undef __FUNCT__
128 #define __FUNCT__ "PCILUSetMatOrdering_ILU"
129 PetscErrorCode PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
130 {
131   PC_ILU         *dir = (PC_ILU*)pc->data;
132   PetscErrorCode ierr;
133   PetscTruth     flg;
134 
135   PetscFunctionBegin;
136   if (!pc->setupcalled) {
137      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
138      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
139   } else {
140     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
141     if (!flg) {
142       pc->setupcalled = 0;
143       ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
144       ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
145       /* free the data structures, then create them again */
146       ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
147     }
148   }
149   PetscFunctionReturn(0);
150 }
151 EXTERN_C_END
152 
153 EXTERN_C_BEGIN
154 #undef __FUNCT__
155 #define __FUNCT__ "PCILUSetReuseOrdering_ILU"
156 PetscErrorCode PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
157 {
158   PC_ILU *ilu;
159 
160   PetscFunctionBegin;
161   ilu                = (PC_ILU*)pc->data;
162   ilu->reuseordering = flag;
163   PetscFunctionReturn(0);
164 }
165 EXTERN_C_END
166 
167 EXTERN_C_BEGIN
168 #undef __FUNCT__
169 #define __FUNCT__ "PCILUDTSetReuseFill_ILUDT"
170 PetscErrorCode PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
171 {
172   PC_ILU *ilu;
173 
174   PetscFunctionBegin;
175   ilu = (PC_ILU*)pc->data;
176   ilu->reusefill = flag;
177   if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported");
178   PetscFunctionReturn(0);
179 }
180 EXTERN_C_END
181 
182 EXTERN_C_BEGIN
183 #undef __FUNCT__
184 #define __FUNCT__ "PCILUSetLevels_ILU"
185 PetscErrorCode PCILUSetLevels_ILU(PC pc,PetscInt levels)
186 {
187   PC_ILU         *ilu;
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   ilu = (PC_ILU*)pc->data;
192 
193   if (!pc->setupcalled) {
194     ilu->info.levels = levels;
195   } else if (ilu->usedt || ilu->info.levels != levels) {
196     ilu->info.levels = levels;
197     pc->setupcalled  = 0;
198     ilu->usedt       = PETSC_FALSE;
199     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
200   }
201   PetscFunctionReturn(0);
202 }
203 EXTERN_C_END
204 
205 EXTERN_C_BEGIN
206 #undef __FUNCT__
207 #define __FUNCT__ "PCILUSetUseInPlace_ILU"
208 PetscErrorCode PCILUSetUseInPlace_ILU(PC pc)
209 {
210   PC_ILU *dir;
211 
212   PetscFunctionBegin;
213   dir          = (PC_ILU*)pc->data;
214   dir->inplace = PETSC_TRUE;
215   PetscFunctionReturn(0);
216 }
217 EXTERN_C_END
218 
219 EXTERN_C_BEGIN
220 #undef __FUNCT__
221 #define __FUNCT__ "PCILUSetAllowDiagonalFill"
222 PetscErrorCode PCILUSetAllowDiagonalFill_ILU(PC pc)
223 {
224   PC_ILU *dir;
225 
226   PetscFunctionBegin;
227   dir                 = (PC_ILU*)pc->data;
228   dir->info.diagonal_fill = 1;
229   PetscFunctionReturn(0);
230 }
231 EXTERN_C_END
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "PCILUSetZeroPivot"
235 /*@
236    PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
237 
238    Collective on PC
239 
240    Input Parameters:
241 +  pc - the preconditioner context
242 -  zero - all pivots smaller than this will be considered zero
243 
244    Options Database Key:
245 .  -pc_ilu_zeropivot <zero> - Sets the zero pivot size
246 
247    Level: intermediate
248 
249 .keywords: PC, set, factorization, direct, fill
250 
251 .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot()
252 @*/
253 PetscErrorCode PCILUSetZeroPivot(PC pc,PetscReal zero)
254 {
255   PetscErrorCode ierr,(*f)(PC,PetscReal);
256 
257   PetscFunctionBegin;
258   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
259   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
260   if (f) {
261     ierr = (*f)(pc,zero);CHKERRQ(ierr);
262   }
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "PCILUSetDamping"
268 /*@
269    PCILUSetDamping - adds this quantity to the diagonal of the matrix during the
270      ILU numerical factorization
271 
272    Collective on PC
273 
274    Input Parameters:
275 +  pc - the preconditioner context
276 -  damping - amount of damping
277 
278    Options Database Key:
279 .  -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default
280 
281    Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
282          pivot, then the damping is doubled until this is alleviated.
283 
284    Level: intermediate
285 
286 .keywords: PC, set, factorization, direct, fill
287 
288 .seealso: PCILUSetFill(), PCLUSetDamping(), PCILUSetShift()
289 @*/
290 PetscErrorCode PCILUSetDamping(PC pc,PetscReal damping)
291 {
292   PetscErrorCode ierr,(*f)(PC,PetscReal);
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
296   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
297   if (f) {
298     ierr = (*f)(pc,damping);CHKERRQ(ierr);
299   }
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "PCILUSetShift"
305 /*@
306    PCILUSetShift - specify whether to use Manteuffel shifting of ILU.
307    If an ILU factorisation breaks down because of nonpositive pivots,
308    adding sufficient identity to the diagonal will remedy this.
309    Setting this causes a bisection method to find the minimum shift that
310    will lead to a well-defined ILU.
311 
312    Input parameters:
313 +  pc - the preconditioner context
314 -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
315 
316    Options Database Key:
317 .  -pc_ilu_shift [1/0] - Activate/Deactivate PCILUSetShift(); the value
318    is optional with 1 being the default
319 
320    Level: intermediate
321 
322 .keywords: PC, indefinite, factorization, incomplete, ILU
323 
324 .seealso: PCILUSetDamping()
325 @*/
326 PetscErrorCode PCILUSetShift(PC pc,PetscTruth shifting)
327 {
328   PetscErrorCode ierr,(*f)(PC,PetscTruth);
329 
330   PetscFunctionBegin;
331   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
332   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr);
333   if (f) {
334     ierr = (*f)(pc,shifting);CHKERRQ(ierr);
335   }
336   PetscFunctionReturn(0);
337 }
338 
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "PCILUSetUseDropTolerance"
342 /*@
343    PCILUSetUseDropTolerance - The preconditioner will use an ILU
344    based on a drop tolerance.
345 
346    Collective on PC
347 
348    Input Parameters:
349 +  pc - the preconditioner context
350 .  dt - the drop tolerance, try from 1.e-10 to .1
351 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
352 -  maxrowcount - the max number of nonzeros allowed in a row, best value
353                  depends on the number of nonzeros in row of original matrix
354 
355    Options Database Key:
356 .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
357 
358    Level: intermediate
359 
360     Notes:
361       This uses the iludt() code of Saad's SPARSKIT package
362 
363 .keywords: PC, levels, reordering, factorization, incomplete, ILU
364 @*/
365 PetscErrorCode PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
366 {
367   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
368 
369   PetscFunctionBegin;
370   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
371   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
372   if (f) {
373     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
374   }
375   PetscFunctionReturn(0);
376 }
377 
378 #undef __FUNCT__
379 #define __FUNCT__ "PCILUSetFill"
380 /*@
381    PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
382    where fill = number nonzeros in factor/number nonzeros in original matrix.
383 
384    Collective on PC
385 
386    Input Parameters:
387 +  pc - the preconditioner context
388 -  fill - amount of expected fill
389 
390    Options Database Key:
391 $  -pc_ilu_fill <fill>
392 
393    Note:
394    For sparse matrix factorizations it is difficult to predict how much
395    fill to expect. By running with the option -log_info PETSc will print the
396    actual amount of fill used; allowing you to set the value accurately for
397    future runs. But default PETSc uses a value of 1.0
398 
399    Level: intermediate
400 
401 .keywords: PC, set, factorization, direct, fill
402 
403 .seealso: PCLUSetFill()
404 @*/
405 PetscErrorCode PCILUSetFill(PC pc,PetscReal fill)
406 {
407   PetscErrorCode ierr,(*f)(PC,PetscReal);
408 
409   PetscFunctionBegin;
410   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
411   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
412   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
413   if (f) {
414     ierr = (*f)(pc,fill);CHKERRQ(ierr);
415   }
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "PCILUSetMatOrdering"
421 /*@C
422     PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to
423     be used in the ILU factorization.
424 
425     Collective on PC
426 
427     Input Parameters:
428 +   pc - the preconditioner context
429 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
430 
431     Options Database Key:
432 .   -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
433 
434     Level: intermediate
435 
436     Notes: natural ordering is used by default
437 
438 .seealso: PCLUSetMatOrdering()
439 
440 .keywords: PC, ILU, set, matrix, reordering
441 
442 @*/
443 PetscErrorCode PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
444 {
445   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
446 
447   PetscFunctionBegin;
448   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
449   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
450   if (f) {
451     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
452   }
453   PetscFunctionReturn(0);
454 }
455 
456 #undef __FUNCT__
457 #define __FUNCT__ "PCILUSetReuseOrdering"
458 /*@
459    PCILUSetReuseOrdering - When similar matrices are factored, this
460    causes the ordering computed in the first factor to be used for all
461    following factors; applies to both fill and drop tolerance ILUs.
462 
463    Collective on PC
464 
465    Input Parameters:
466 +  pc - the preconditioner context
467 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
468 
469    Options Database Key:
470 .  -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()
471 
472    Level: intermediate
473 
474 .keywords: PC, levels, reordering, factorization, incomplete, ILU
475 
476 .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
477 @*/
478 PetscErrorCode PCILUSetReuseOrdering(PC pc,PetscTruth flag)
479 {
480   PetscErrorCode ierr,(*f)(PC,PetscTruth);
481 
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
484   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
485   if (f) {
486     ierr = (*f)(pc,flag);CHKERRQ(ierr);
487   }
488   PetscFunctionReturn(0);
489 }
490 
491 #undef __FUNCT__
492 #define __FUNCT__ "PCILUDTSetReuseFill"
493 /*@
494    PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
495    this causes later ones to use the fill computed in the initial factorization.
496 
497    Collective on PC
498 
499    Input Parameters:
500 +  pc - the preconditioner context
501 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
502 
503    Options Database Key:
504 .  -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()
505 
506    Level: intermediate
507 
508 .keywords: PC, levels, reordering, factorization, incomplete, ILU
509 
510 .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
511 @*/
512 PetscErrorCode PCILUDTSetReuseFill(PC pc,PetscTruth flag)
513 {
514   PetscErrorCode ierr,(*f)(PC,PetscTruth);
515 
516   PetscFunctionBegin;
517   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
518   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
519   if (f) {
520     ierr = (*f)(pc,flag);CHKERRQ(ierr);
521   }
522   PetscFunctionReturn(0);
523 }
524 
525 #undef __FUNCT__
526 #define __FUNCT__ "PCILUSetLevels"
527 /*@
528    PCILUSetLevels - Sets the number of levels of fill to use.
529 
530    Collective on PC
531 
532    Input Parameters:
533 +  pc - the preconditioner context
534 -  levels - number of levels of fill
535 
536    Options Database Key:
537 .  -pc_ilu_levels <levels> - Sets fill level
538 
539    Level: intermediate
540 
541 .keywords: PC, levels, fill, factorization, incomplete, ILU
542 @*/
543 PetscErrorCode PCILUSetLevels(PC pc,PetscInt levels)
544 {
545   PetscErrorCode ierr,(*f)(PC,PetscInt);
546 
547   PetscFunctionBegin;
548   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
549   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
550   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
551   if (f) {
552     ierr = (*f)(pc,levels);CHKERRQ(ierr);
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "PCILUSetAllowDiagonalFill"
559 /*@
560    PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be
561    treated as level 0 fill even if there is no non-zero location.
562 
563    Collective on PC
564 
565    Input Parameters:
566 +  pc - the preconditioner context
567 
568    Options Database Key:
569 .  -pc_ilu_diagonal_fill
570 
571    Notes:
572    Does not apply with 0 fill.
573 
574    Level: intermediate
575 
576 .keywords: PC, levels, fill, factorization, incomplete, ILU
577 @*/
578 PetscErrorCode PCILUSetAllowDiagonalFill(PC pc)
579 {
580   PetscErrorCode ierr,(*f)(PC);
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
584   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
585   if (f) {
586     ierr = (*f)(pc);CHKERRQ(ierr);
587   }
588   PetscFunctionReturn(0);
589 }
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "PCILUSetUseInPlace"
593 /*@
594    PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
595    Collective on PC
596 
597    Input Parameters:
598 .  pc - the preconditioner context
599 
600    Options Database Key:
601 .  -pc_ilu_in_place - Activates in-place factorization
602 
603    Notes:
604    PCILUSetUseInPlace() is intended for use with matrix-free variants of
605    Krylov methods, or when a different matrices are employed for the linear
606    system and preconditioner, or with ASM preconditioning.  Do NOT use
607    this option if the linear system
608    matrix also serves as the preconditioning matrix, since the factored
609    matrix would then overwrite the original matrix.
610 
611    Only works well with ILU(0).
612 
613    Level: intermediate
614 
615 .keywords: PC, set, factorization, inplace, in-place, ILU
616 
617 .seealso:  PCLUSetUseInPlace()
618 @*/
619 PetscErrorCode PCILUSetUseInPlace(PC pc)
620 {
621   PetscErrorCode ierr,(*f)(PC);
622 
623   PetscFunctionBegin;
624   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
625   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
626   if (f) {
627     ierr = (*f)(pc);CHKERRQ(ierr);
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "PCILUReorderForNonzeroDiagonal"
634 /*@
635    PCILUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
636 
637    Collective on PC
638 
639    Input Parameters:
640 +  pc - the preconditioner context
641 -  tol - diagonal entries smaller than this in absolute value are considered zero
642 
643    Options Database Key:
644 .  -pc_lu_nonzeros_along_diagonal
645 
646    Level: intermediate
647 
648 .keywords: PC, set, factorization, direct, fill
649 
650 .seealso: PCILUSetFill(), PCILUSetDamp(), PCIILUSetZeroPivot(), MatReorderForNonzeroDiagonal()
651 @*/
652 PetscErrorCode PCILUReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
653 {
654   PetscErrorCode ierr,(*f)(PC,PetscReal);
655 
656   PetscFunctionBegin;
657   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
658   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
659   if (f) {
660     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
661   }
662   PetscFunctionReturn(0);
663 }
664 
665 #undef __FUNCT__
666 #define __FUNCT__ "PCILUSetPivotInBlocks"
667 /*@
668     PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
669       with BAIJ or SBAIJ matrices
670 
671     Collective on PC
672 
673     Input Parameters:
674 +   pc - the preconditioner context
675 -   pivot - PETSC_TRUE or PETSC_FALSE
676 
677     Options Database Key:
678 .   -pc_ilu_pivot_in_blocks <true,false>
679 
680     Level: intermediate
681 
682 .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
683 @*/
684 PetscErrorCode PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
685 {
686   PetscErrorCode ierr,(*f)(PC,PetscTruth);
687 
688   PetscFunctionBegin;
689   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
690   if (f) {
691     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
692   }
693   PetscFunctionReturn(0);
694 }
695 
696 /* ------------------------------------------------------------------------------------------*/
697 
698 EXTERN_C_BEGIN
699 #undef __FUNCT__
700 #define __FUNCT__ "PCILUSetPivotInBlocks_ILU"
701 PetscErrorCode PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
702 {
703   PC_ILU *dir = (PC_ILU*)pc->data;
704 
705   PetscFunctionBegin;
706   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
707   PetscFunctionReturn(0);
708 }
709 EXTERN_C_END
710 
711 #undef __FUNCT__
712 #define __FUNCT__ "PCSetFromOptions_ILU"
713 static PetscErrorCode PCSetFromOptions_ILU(PC pc)
714 {
715   PetscErrorCode ierr;
716   PetscInt       dtmax = 3,itmp;
717   PetscTruth     flg,set;
718   PetscReal      dt[3];
719   char           tname[256];
720   PC_ILU         *ilu = (PC_ILU*)pc->data;
721   PetscFList     ordlist;
722   PetscReal      tol;
723 
724   PetscFunctionBegin;
725   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
726   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
727     ierr = PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr);
728     if (flg) ilu->info.levels = itmp;
729     ierr = PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);CHKERRQ(ierr);
730     ierr = PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);CHKERRQ(ierr);
731     ilu->info.diagonal_fill = (double) flg;
732     ierr = PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);CHKERRQ(ierr);
733     ierr = PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr);
734 
735     ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
736     if (flg) {
737       ierr = PCFactorSetShiftNonzero((PetscReal) PETSC_DECIDE,&ilu->info);CHKERRQ(ierr);
738     }
739     ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.damping,&ilu->info.damping,0);CHKERRQ(ierr);
740 
741     ierr = PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);CHKERRQ(ierr);
742     if (flg) {
743       ierr = PetscOptionsInt("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",(PetscInt)ilu->info.shift,&itmp,&flg); CHKERRQ(ierr);
744       if (flg && !itmp) {
745 	ierr = PCILUSetShift(pc,PETSC_FALSE);CHKERRQ(ierr);
746       } else {
747 	ierr = PCILUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr);
748       }
749     }
750     ierr = PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr);
751 
752     dt[0] = ilu->info.dt;
753     dt[1] = ilu->info.dtcol;
754     dt[2] = ilu->info.dtcount;
755     ierr = PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
756     if (flg) {
757       ierr = PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
758     }
759     ierr = PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr);
760     ierr = PetscOptionsName("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
761     if (flg) {
762       tol = PETSC_DECIDE;
763       ierr = PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
764       ierr = PCILUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
765     }
766 
767     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
768     ierr = PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr);
769     if (flg) {
770       ierr = PCILUSetMatOrdering(pc,tname);CHKERRQ(ierr);
771     }
772     flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
773     ierr = PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
774     if (set) {
775       ierr = PCILUSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
776     }
777   ierr = PetscOptionsTail();CHKERRQ(ierr);
778   PetscFunctionReturn(0);
779 }
780 
781 #undef __FUNCT__
782 #define __FUNCT__ "PCView_ILU"
783 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
784 {
785   PC_ILU         *ilu = (PC_ILU*)pc->data;
786   PetscErrorCode ierr;
787   PetscTruth     isstring,iascii;
788 
789   PetscFunctionBegin;
790   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
791   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
792   if (iascii) {
793     if (ilu->usedt) {
794         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %g\n",ilu->info.dt);CHKERRQ(ierr);
795         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr);
796         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %g\n",ilu->info.dtcol);CHKERRQ(ierr);
797     } else if (ilu->info.levels == 1) {
798         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
799     } else {
800         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
801     }
802     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max fill ratio allocated %g\n",ilu->info.fill);CHKERRQ(ierr);
803     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);CHKERRQ(ierr);
804     if (ilu->info.shift) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
805     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
806     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
807     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr);
808     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
809     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
810     if (ilu->fact) {
811       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
812       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
813       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
814       ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr);
815       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
816       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
817     }
818   } else if (isstring) {
819     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
820   } else {
821     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
822   }
823   PetscFunctionReturn(0);
824 }
825 
826 #undef __FUNCT__
827 #define __FUNCT__ "PCSetUp_ILU"
828 static PetscErrorCode PCSetUp_ILU(PC pc)
829 {
830   PetscErrorCode ierr;
831   PC_ILU         *ilu = (PC_ILU*)pc->data;
832 
833   PetscFunctionBegin;
834   if (ilu->inplace) {
835     if (!pc->setupcalled) {
836 
837       /* In-place factorization only makes sense with the natural ordering,
838          so we only need to get the ordering once, even if nonzero structure changes */
839       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
840       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
841       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
842     }
843 
844     /* In place ILU only makes sense with fill factor of 1.0 because
845        cannot have levels of fill */
846     ilu->info.fill          = 1.0;
847     ilu->info.diagonal_fill = 0;
848     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
849     ilu->fact = pc->pmat;
850   } else if (ilu->usedt) {
851     if (!pc->setupcalled) {
852       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
853       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
854       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
855       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
856       PetscLogObjectParent(pc,ilu->fact);
857     } else if (pc->flag != SAME_NONZERO_PATTERN) {
858       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
859       if (!ilu->reuseordering) {
860         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
861         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
862         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
863         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
864         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
865       }
866       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
867       PetscLogObjectParent(pc,ilu->fact);
868     } else if (!ilu->reusefill) {
869       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
870       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
871       PetscLogObjectParent(pc,ilu->fact);
872     } else {
873       ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
874     }
875    } else {
876     if (!pc->setupcalled) {
877       /* first time in so compute reordering and symbolic factorization */
878       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
879       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
880       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
881       /*  Remove zeros along diagonal?     */
882       if (ilu->nonzerosalongdiagonal) {
883         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
884       }
885       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
886       PetscLogObjectParent(pc,ilu->fact);
887     } else if (pc->flag != SAME_NONZERO_PATTERN) {
888       if (!ilu->reuseordering) {
889         /* compute a new ordering for the ILU */
890         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
891         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
892         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
893         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
894         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
895         /*  Remove zeros along diagonal?     */
896         if (ilu->nonzerosalongdiagonal) {
897           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
898         }
899       }
900       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
901       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
902       PetscLogObjectParent(pc,ilu->fact);
903     }
904     ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
905   }
906   PetscFunctionReturn(0);
907 }
908 
909 #undef __FUNCT__
910 #define __FUNCT__ "PCDestroy_ILU"
911 static PetscErrorCode PCDestroy_ILU(PC pc)
912 {
913   PC_ILU         *ilu = (PC_ILU*)pc->data;
914   PetscErrorCode ierr;
915 
916   PetscFunctionBegin;
917   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
918   ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr);
919   ierr = PetscFree(ilu);CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "PCApply_ILU"
925 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
926 {
927   PC_ILU         *ilu = (PC_ILU*)pc->data;
928   PetscErrorCode ierr;
929 
930   PetscFunctionBegin;
931   ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr);
932   PetscFunctionReturn(0);
933 }
934 
935 #undef __FUNCT__
936 #define __FUNCT__ "PCApplyTranspose_ILU"
937 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
938 {
939   PC_ILU         *ilu = (PC_ILU*)pc->data;
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 #undef __FUNCT__
948 #define __FUNCT__ "PCGetFactoredMatrix_ILU"
949 static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
950 {
951   PC_ILU *ilu = (PC_ILU*)pc->data;
952 
953   PetscFunctionBegin;
954   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
955   *mat = ilu->fact;
956   PetscFunctionReturn(0);
957 }
958 
959 /*MC
960      PCILU - Incomplete factorization preconditioners.
961 
962    Options Database Keys:
963 +  -pc_ilu_levels <k> - number of levels of fill for ILU(k)
964 .  -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
965                       its factorization (overwrites original matrix)
966 .  -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
967 .  -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization
968 .  -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots
969 .  -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner
970 .  -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot
971 .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
972 .  -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
973 .  -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
974                                    this decreases the chance of getting a zero pivot
975 .  -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
976 -  -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
977                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
978                              stability of the ILU factorization
979 
980    Level: beginner
981 
982   Concepts: incomplete factorization
983 
984    Notes: Only implemented for some matrix formats. Not implemented in parallel
985 
986           For BAIJ matrices this implements a point block ILU
987 
988 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
989            PCILUSetSetZeroPivot(), PCILUSetDamping(), PCILUSetShift(), PCILUSetUseDropTolerance(),
990            PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(),
991            PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(),
992 
993 M*/
994 
995 EXTERN_C_BEGIN
996 #undef __FUNCT__
997 #define __FUNCT__ "PCCreate_ILU"
998 PetscErrorCode PCCreate_ILU(PC pc)
999 {
1000   PetscErrorCode ierr;
1001   PC_ILU         *ilu;
1002 
1003   PetscFunctionBegin;
1004   ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr);
1005   PetscLogObjectMemory(pc,sizeof(PC_ILU));
1006 
1007   ilu->fact                    = 0;
1008   ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr);
1009   ilu->info.levels             = 0;
1010   ilu->info.fill               = 1.0;
1011   ilu->col                     = 0;
1012   ilu->row                     = 0;
1013   ilu->inplace                 = PETSC_FALSE;
1014   ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr);
1015   ilu->reuseordering           = PETSC_FALSE;
1016   ilu->usedt                   = PETSC_FALSE;
1017   ilu->info.dt                 = PETSC_DEFAULT;
1018   ilu->info.dtcount            = PETSC_DEFAULT;
1019   ilu->info.dtcol              = PETSC_DEFAULT;
1020   ilu->info.damping            = 0.0;
1021   ilu->info.shift              = PETSC_FALSE;
1022   ilu->info.shift_fraction     = 0.0;
1023   ilu->info.zeropivot          = 1.e-12;
1024   ilu->info.pivotinblocks      = 1.0;
1025   ilu->reusefill               = PETSC_FALSE;
1026   ilu->info.diagonal_fill      = 0;
1027   pc->data                     = (void*)ilu;
1028 
1029   pc->ops->destroy             = PCDestroy_ILU;
1030   pc->ops->apply               = PCApply_ILU;
1031   pc->ops->applytranspose      = PCApplyTranspose_ILU;
1032   pc->ops->setup               = PCSetUp_ILU;
1033   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
1034   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ILU;
1035   pc->ops->view                = PCView_ILU;
1036   pc->ops->applyrichardson     = 0;
1037 
1038   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
1039                     PCILUSetUseDropTolerance_ILU);CHKERRQ(ierr);
1040   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
1041                     PCILUSetFill_ILU);CHKERRQ(ierr);
1042   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
1043                     PCILUSetDamping_ILU);CHKERRQ(ierr);
1044   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU",
1045  		    PCILUSetShift_ILU);CHKERRQ(ierr);
1046   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
1047                     PCILUSetMatOrdering_ILU);CHKERRQ(ierr);
1048   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
1049                     PCILUSetReuseOrdering_ILU);CHKERRQ(ierr);
1050   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
1051                     PCILUDTSetReuseFill_ILUDT);CHKERRQ(ierr);
1052   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
1053                     PCILUSetLevels_ILU);CHKERRQ(ierr);
1054   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
1055                     PCILUSetUseInPlace_ILU);CHKERRQ(ierr);
1056   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
1057                     PCILUSetAllowDiagonalFill_ILU);CHKERRQ(ierr);
1058   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
1059                     PCILUSetPivotInBlocks_ILU);CHKERRQ(ierr);
1060   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU",
1061                     PCILUSetZeroPivot_ILU);CHKERRQ(ierr);
1062   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C","PCILUReorderForNonzeroDiagonal_ILU",
1063                     PCILUReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
1064   PetscFunctionReturn(0);
1065 }
1066 EXTERN_C_END
1067