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