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