xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision b271bb04eaec0d477b629c8be871cf42cb4980f5)
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 "include/private/matimpl.h"
9 
10 /* ------------------------------------------------------------------------------------------*/
11 EXTERN_C_BEGIN
12 #undef __FUNCT__
13 #define __FUNCT__ "PCFactorSetReuseFill_ILU"
14 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag)
15 {
16   PC_ILU *lu;
17 
18   PetscFunctionBegin;
19   lu = (PC_ILU*)pc->data;
20   lu->reusefill = flag;
21   PetscFunctionReturn(0);
22 }
23 EXTERN_C_END
24 
25 EXTERN_C_BEGIN
26 #undef __FUNCT__
27 #define __FUNCT__ "PCFactorSetZeroPivot_ILU"
28 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ILU(PC pc,PetscReal z)
29 {
30   PC_ILU *ilu;
31 
32   PetscFunctionBegin;
33   ilu                 = (PC_ILU*)pc->data;
34   ilu->info.zeropivot = z;
35   PetscFunctionReturn(0);
36 }
37 EXTERN_C_END
38 
39 EXTERN_C_BEGIN
40 #undef __FUNCT__
41 #define __FUNCT__ "PCFactorSetShiftNonzero_ILU"
42 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ILU(PC pc,PetscReal shift)
43 {
44   PC_ILU *dir;
45 
46   PetscFunctionBegin;
47   dir = (PC_ILU*)pc->data;
48   if (shift == (PetscReal) PETSC_DECIDE) {
49     dir->info.shiftnz = 1.e-12;
50   } else {
51     dir->info.shiftnz = shift;
52   }
53   PetscFunctionReturn(0);
54 }
55 EXTERN_C_END
56 
57 EXTERN_C_BEGIN
58 #undef __FUNCT__
59 #define __FUNCT__ "PCFactorSetShiftPd_ILU"
60 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ILU(PC pc,PetscTruth shift)
61 {
62   PC_ILU *dir;
63 
64   PetscFunctionBegin;
65   dir = (PC_ILU*)pc->data;
66   if (shift) {
67     dir->info.shiftpd = 1.0;
68   } else {
69     dir->info.shiftpd = 0.0;
70   }
71   PetscFunctionReturn(0);
72 }
73 EXTERN_C_END
74 
75 EXTERN_C_BEGIN
76 #undef __FUNCT__
77 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU"
78 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
79 {
80   PC_ILU *ilu = (PC_ILU*)pc->data;
81 
82   PetscFunctionBegin;
83   ilu->nonzerosalongdiagonal = PETSC_TRUE;
84   if (z == PETSC_DECIDE) {
85     ilu->nonzerosalongdiagonaltol = 1.e-10;
86   } else {
87     ilu->nonzerosalongdiagonaltol = z;
88   }
89   PetscFunctionReturn(0);
90 }
91 EXTERN_C_END
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "PCDestroy_ILU_Internal"
95 PetscErrorCode PCDestroy_ILU_Internal(PC pc)
96 {
97   PC_ILU         *ilu = (PC_ILU*)pc->data;
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);}
102   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
103   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
104   PetscFunctionReturn(0);
105 }
106 
107 EXTERN_C_BEGIN
108 #undef __FUNCT__
109 #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU"
110 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
111 {
112   PC_ILU         *ilu;
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin;
116   ilu = (PC_ILU*)pc->data;
117   if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) {
118     pc->setupcalled   = 0;
119     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
120   }
121   ilu->usedt        = PETSC_TRUE;
122   ilu->info.dt      = dt;
123   ilu->info.dtcol   = dtcol;
124   ilu->info.dtcount = dtcount;
125   ilu->info.fill    = PETSC_DEFAULT;
126   PetscFunctionReturn(0);
127 }
128 EXTERN_C_END
129 
130 EXTERN_C_BEGIN
131 #undef __FUNCT__
132 #define __FUNCT__ "PCFactorSetFill_ILU"
133 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ILU(PC pc,PetscReal fill)
134 {
135   PC_ILU *dir;
136 
137   PetscFunctionBegin;
138   dir            = (PC_ILU*)pc->data;
139   dir->info.fill = fill;
140   PetscFunctionReturn(0);
141 }
142 EXTERN_C_END
143 
144 EXTERN_C_BEGIN
145 #undef __FUNCT__
146 #define __FUNCT__ "PCFactorSetMatOrderingType_ILU"
147 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ILU(PC pc,const MatOrderingType ordering)
148 {
149   PC_ILU         *dir = (PC_ILU*)pc->data;
150   PetscErrorCode ierr;
151   PetscTruth     flg;
152 
153   PetscFunctionBegin;
154   if (!pc->setupcalled) {
155      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
156      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
157   } else {
158     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
159     if (!flg) {
160       pc->setupcalled = 0;
161       ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
162       ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
163       /* free the data structures, then create them again */
164       ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
165     }
166   }
167   PetscFunctionReturn(0);
168 }
169 EXTERN_C_END
170 
171 EXTERN_C_BEGIN
172 #undef __FUNCT__
173 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
174 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
175 {
176   PC_ILU *ilu;
177 
178   PetscFunctionBegin;
179   ilu                = (PC_ILU*)pc->data;
180   ilu->reuseordering = flag;
181   PetscFunctionReturn(0);
182 }
183 EXTERN_C_END
184 
185 EXTERN_C_BEGIN
186 #undef __FUNCT__
187 #define __FUNCT__ "PCFactorSetLevels_ILU"
188 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ILU(PC pc,PetscInt levels)
189 {
190   PC_ILU         *ilu;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   ilu = (PC_ILU*)pc->data;
195 
196   if (!pc->setupcalled) {
197     ilu->info.levels = levels;
198   } else if (ilu->usedt || ilu->info.levels != levels) {
199     ilu->info.levels = levels;
200     pc->setupcalled  = 0;
201     ilu->usedt       = PETSC_FALSE;
202     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
203   }
204   PetscFunctionReturn(0);
205 }
206 EXTERN_C_END
207 
208 EXTERN_C_BEGIN
209 #undef __FUNCT__
210 #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
211 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc)
212 {
213   PC_ILU *dir;
214 
215   PetscFunctionBegin;
216   dir          = (PC_ILU*)pc->data;
217   dir->inplace = PETSC_TRUE;
218   PetscFunctionReturn(0);
219 }
220 EXTERN_C_END
221 
222 EXTERN_C_BEGIN
223 #undef __FUNCT__
224 #define __FUNCT__ "PCFactorSetAllowDiagonalFill_ILU"
225 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_ILU(PC pc)
226 {
227   PC_ILU *dir;
228 
229   PetscFunctionBegin;
230   dir                 = (PC_ILU*)pc->data;
231   dir->info.diagonal_fill = 1;
232   PetscFunctionReturn(0);
233 }
234 EXTERN_C_END
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "PCFactorSetUseDropTolerance"
238 /*@
239    PCFactorSetUseDropTolerance - The preconditioner will use an ILU
240    based on a drop tolerance.
241 
242    Collective on PC
243 
244    Input Parameters:
245 +  pc - the preconditioner context
246 .  dt - the drop tolerance, try from 1.e-10 to .1
247 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
248 -  maxrowcount - the max number of nonzeros allowed in a row, best value
249                  depends on the number of nonzeros in row of original matrix
250 
251    Options Database Key:
252 .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
253 
254    Level: intermediate
255 
256     Notes:
257       This uses the iludt() code of Saad's SPARSKIT package
258 
259       There are NO default values for the 3 parameters, you must set them with reasonable values for your
260       matrix. We don't know how to compute reasonable values.
261 
262 .keywords: PC, levels, reordering, factorization, incomplete, ILU
263 @*/
264 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
265 {
266   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
267 
268   PetscFunctionBegin;
269   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
270   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
271   if (f) {
272     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
273   }
274   PetscFunctionReturn(0);
275 }
276 
277 #undef __FUNCT__
278 #define __FUNCT__ "PCFactorSetLevels"
279 /*@
280    PCFactorSetLevels - Sets the number of levels of fill to use.
281 
282    Collective on PC
283 
284    Input Parameters:
285 +  pc - the preconditioner context
286 -  levels - number of levels of fill
287 
288    Options Database Key:
289 .  -pc_factor_levels <levels> - Sets fill level
290 
291    Level: intermediate
292 
293 .keywords: PC, levels, fill, factorization, incomplete, ILU
294 @*/
295 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
296 {
297   PetscErrorCode ierr,(*f)(PC,PetscInt);
298 
299   PetscFunctionBegin;
300   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
301   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
302   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
303   if (f) {
304     ierr = (*f)(pc,levels);CHKERRQ(ierr);
305   }
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
311 /*@
312    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
313    treated as level 0 fill even if there is no non-zero location.
314 
315    Collective on PC
316 
317    Input Parameters:
318 +  pc - the preconditioner context
319 
320    Options Database Key:
321 .  -pc_factor_diagonal_fill
322 
323    Notes:
324    Does not apply with 0 fill.
325 
326    Level: intermediate
327 
328 .keywords: PC, levels, fill, factorization, incomplete, ILU
329 @*/
330 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
331 {
332   PetscErrorCode ierr,(*f)(PC);
333 
334   PetscFunctionBegin;
335   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
336   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
337   if (f) {
338     ierr = (*f)(pc);CHKERRQ(ierr);
339   }
340   PetscFunctionReturn(0);
341 }
342 
343 /* ------------------------------------------------------------------------------------------*/
344 
345 EXTERN_C_BEGIN
346 #undef __FUNCT__
347 #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU"
348 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
349 {
350   PC_ILU *dir = (PC_ILU*)pc->data;
351 
352   PetscFunctionBegin;
353   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
354   PetscFunctionReturn(0);
355 }
356 EXTERN_C_END
357 
358 EXTERN_C_BEGIN
359 #undef __FUNCT__
360 #define __FUNCT__ "PCFactorSetShiftInBlocks_ILU"
361 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_ILU(PC pc,PetscReal shift)
362 {
363   PC_ILU *dir = (PC_ILU*)pc->data;
364 
365   PetscFunctionBegin;
366   if (shift == PETSC_DEFAULT) {
367     dir->info.shiftinblocks = 1.e-12;
368   } else {
369     dir->info.shiftinblocks = shift;
370   }
371   PetscFunctionReturn(0);
372 }
373 EXTERN_C_END
374 
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "PCSetFromOptions_ILU"
378 static PetscErrorCode PCSetFromOptions_ILU(PC pc)
379 {
380   PetscErrorCode ierr;
381   PetscInt       dtmax = 3,itmp;
382   PetscTruth     flg,set;
383   PetscReal      dt[3];
384   char           tname[256];
385   PC_ILU         *ilu = (PC_ILU*)pc->data;
386   PetscFList     ordlist;
387   PetscReal      tol;
388 
389   PetscFunctionBegin;
390   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
391   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
392     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr);
393     if (flg) ilu->info.levels = itmp;
394     ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
395     if (flg) ilu->inplace = PETSC_TRUE;
396     ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr);
397     ilu->info.diagonal_fill = (double) flg;
398     ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
399     if (flg) ilu->reusefill = PETSC_TRUE;
400     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
401     if (flg) ilu->reuseordering = PETSC_TRUE;
402     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
403     if (flg) {
404       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
405     }
406     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr);
407     flg = (ilu->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE;
408     ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
409     ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr);
410     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr);
411 
412     dt[0] = ilu->info.dt;
413     dt[1] = ilu->info.dtcol;
414     dt[2] = ilu->info.dtcount;
415     ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
416     if (flg) {
417       ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
418     }
419     ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr);
420     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
421     if (flg) {
422       tol = PETSC_DECIDE;
423       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
424       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
425     }
426 
427     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
428     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr);
429     if (flg) {
430       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
431     }
432     flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
433     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
434     if (set) {
435       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
436     }
437     ierr = PetscOptionsName("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",&flg);CHKERRQ(ierr);
438     if (flg) {
439       ierr = PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
440     }
441     ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",ilu->info.shiftinblocks,&ilu->info.shiftinblocks,0);CHKERRQ(ierr);
442 
443   ierr = PetscOptionsTail();CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 #undef __FUNCT__
448 #define __FUNCT__ "PCView_ILU"
449 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
450 {
451   PC_ILU         *ilu = (PC_ILU*)pc->data;
452   PetscErrorCode ierr;
453   PetscTruth     isstring,iascii;
454 
455   PetscFunctionBegin;
456   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
457   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
458   if (iascii) {
459     if (ilu->usedt) {
460         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr);
461         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr);
462         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr);
463     } else if (ilu->info.levels == 1) {
464         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
465     } else {
466         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
467     }
468     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr);
469     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr);
470     if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
471     if (ilu->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
472     if (ilu->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);}
473     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
474     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
475     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr);
476     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
477     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
478     if (ilu->fact) {
479       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr);
480       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
481       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
482       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
483       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
484       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
485       ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr);
486       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
487       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
488       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
489       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
490     }
491   } else if (isstring) {
492     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
493   } else {
494     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
495   }
496   PetscFunctionReturn(0);
497 }
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "PCSetUp_ILU"
501 static PetscErrorCode PCSetUp_ILU(PC pc)
502 {
503   PetscErrorCode ierr;
504   PC_ILU         *ilu = (PC_ILU*)pc->data;
505   MatInfo        info;
506 
507   PetscFunctionBegin;
508   if (ilu->inplace) {
509     if (!pc->setupcalled) {
510 
511       /* In-place factorization only makes sense with the natural ordering,
512          so we only need to get the ordering once, even if nonzero structure changes */
513       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
514       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
515       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
516     }
517 
518     /* In place ILU only makes sense with fill factor of 1.0 because
519        cannot have levels of fill */
520     ilu->info.fill          = 1.0;
521     ilu->info.diagonal_fill = 0;
522     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
523     ilu->fact = pc->pmat;
524   } else if (ilu->usedt) {
525     if (!pc->setupcalled) {
526       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
527       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
528       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
529       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
530       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
531     } else if (pc->flag != SAME_NONZERO_PATTERN) {
532       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
533       if (!ilu->reuseordering) {
534         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
535         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
536         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
537         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
538         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
539       }
540       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
541       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
542     } else if (!ilu->reusefill) {
543       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
544       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
545       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
546     } else {
547       ierr = MatLUFactorNumeric(ilu->fact,pc->pmat,&ilu->info);CHKERRQ(ierr);
548     }
549   } else {
550     if (!pc->setupcalled) {
551       /* first time in so compute reordering and symbolic factorization */
552       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
553       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
554       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
555       /*  Remove zeros along diagonal?     */
556       if (ilu->nonzerosalongdiagonal) {
557         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
558       }
559       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&ilu->fact);CHKERRQ(ierr);
560       ierr = MatILUFactorSymbolic(ilu->fact,pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
561       ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
562       ilu->actualfill = info.fill_ratio_needed;
563       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
564     } else if (pc->flag != SAME_NONZERO_PATTERN) {
565       if (!ilu->reuseordering) {
566         /* compute a new ordering for the ILU */
567         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
568         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
569         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
570         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
571         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
572         /*  Remove zeros along diagonal?     */
573         if (ilu->nonzerosalongdiagonal) {
574           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
575         }
576       }
577       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
578       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&ilu->fact);CHKERRQ(ierr);
579       ierr = MatILUFactorSymbolic(ilu->fact,pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
580       ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
581       ilu->actualfill = info.fill_ratio_needed;
582       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
583     }
584     ierr = MatLUFactorNumeric(ilu->fact,pc->pmat,&ilu->info);CHKERRQ(ierr);
585   }
586   PetscFunctionReturn(0);
587 }
588 
589 #undef __FUNCT__
590 #define __FUNCT__ "PCDestroy_ILU"
591 static PetscErrorCode PCDestroy_ILU(PC pc)
592 {
593   PC_ILU         *ilu = (PC_ILU*)pc->data;
594   PetscErrorCode ierr;
595 
596   PetscFunctionBegin;
597   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
598   ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr);
599   ierr = PetscFree(ilu);CHKERRQ(ierr);
600   PetscFunctionReturn(0);
601 }
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "PCApply_ILU"
605 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
606 {
607   PC_ILU         *ilu = (PC_ILU*)pc->data;
608   PetscErrorCode ierr;
609 
610   PetscFunctionBegin;
611   ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr);
612   PetscFunctionReturn(0);
613 }
614 
615 #undef __FUNCT__
616 #define __FUNCT__ "PCApplyTranspose_ILU"
617 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
618 {
619   PC_ILU         *ilu = (PC_ILU*)pc->data;
620   PetscErrorCode ierr;
621 
622   PetscFunctionBegin;
623   ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr);
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "PCFactorGetMatrix_ILU"
629 static PetscErrorCode PCFactorGetMatrix_ILU(PC pc,Mat *mat)
630 {
631   PC_ILU *ilu = (PC_ILU*)pc->data;
632 
633   PetscFunctionBegin;
634   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
635   *mat = ilu->fact;
636   PetscFunctionReturn(0);
637 }
638 
639 #undef __FUNCT__
640 #define __FUNCT__ "PCFactorSetShiftInBlocks"
641 /*@
642     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot.
643 
644     Collective on PC
645 
646     Input Parameters:
647 +   pc - the preconditioner context
648 -   shift - amount to shift or PETSC_DECIDE
649 
650     Options Database Key:
651 .   -pc_factor_shift_in_blocks <shift>
652 
653     Level: intermediate
654 
655 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks()
656 @*/
657 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
658 {
659   PetscErrorCode ierr,(*f)(PC,PetscReal);
660 
661   PetscFunctionBegin;
662   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
663   if (f) {
664     ierr = (*f)(pc,shift);CHKERRQ(ierr);
665   }
666   PetscFunctionReturn(0);
667 }
668 
669 /*MC
670      PCILU - Incomplete factorization preconditioners.
671 
672    Options Database Keys:
673 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
674 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
675                       its factorization (overwrites original matrix)
676 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
677 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
678 .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
679 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
680 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
681                                    this decreases the chance of getting a zero pivot
682 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
683 .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
684                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
685                              stability of the ILU factorization
686 .  -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization
687 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
688 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
689    is optional with PETSC_TRUE being the default
690 
691    Level: beginner
692 
693   Concepts: incomplete factorization
694 
695    Notes: Only implemented for some matrix formats. (for parallel use you
696              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended
697              unless you really want a parallel ILU).
698 
699           For BAIJ matrices this implements a point block ILU
700 
701    References:
702    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
703    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
704 
705    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
706    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
707 
708    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
709       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
710       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
711       Science and Engineering, Kluwer, pp. 167--202.
712 
713 
714 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
715            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(),
716            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
717            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
718            PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
719 
720 M*/
721 
722 EXTERN_C_BEGIN
723 #undef __FUNCT__
724 #define __FUNCT__ "PCCreate_ILU"
725 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
726 {
727   PetscErrorCode ierr;
728   PC_ILU         *ilu;
729 
730   PetscFunctionBegin;
731   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
732 
733   ilu->fact                    = 0;
734   ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr);
735   ilu->info.levels             = 0;
736   ilu->info.fill               = 1.0;
737   ilu->col                     = 0;
738   ilu->row                     = 0;
739   ilu->inplace                 = PETSC_FALSE;
740   ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr);
741   ilu->reuseordering           = PETSC_FALSE;
742   ilu->usedt                   = PETSC_FALSE;
743   ilu->info.dt                 = PETSC_DEFAULT;
744   ilu->info.dtcount            = PETSC_DEFAULT;
745   ilu->info.dtcol              = PETSC_DEFAULT;
746   ilu->info.shiftnz            = 1.e-12;
747   ilu->info.shiftpd            = 0.0; /* false */
748   ilu->info.zeropivot          = 1.e-12;
749   ilu->info.pivotinblocks      = 1.0;
750   ilu->info.shiftinblocks      = 1.e-12;
751   ilu->reusefill               = PETSC_FALSE;
752   ilu->info.diagonal_fill      = 0;
753   pc->data                     = (void*)ilu;
754 
755   pc->ops->destroy             = PCDestroy_ILU;
756   pc->ops->apply               = PCApply_ILU;
757   pc->ops->applytranspose      = PCApplyTranspose_ILU;
758   pc->ops->setup               = PCSetUp_ILU;
759   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
760   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_ILU;
761   pc->ops->view                = PCView_ILU;
762   pc->ops->applyrichardson     = 0;
763 
764   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU",
765                     PCFactorSetZeroPivot_ILU);CHKERRQ(ierr);
766   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU",
767                     PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr);
768   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU",
769                     PCFactorSetShiftPd_ILU);CHKERRQ(ierr);
770 
771   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
772                     PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr);
773   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU",
774                     PCFactorSetFill_ILU);CHKERRQ(ierr);
775   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ILU",
776                     PCFactorSetMatOrderingType_ILU);CHKERRQ(ierr);
777   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
778                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
779   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
780                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
781   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU",
782                     PCFactorSetLevels_ILU);CHKERRQ(ierr);
783   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
784                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
785   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU",
786                     PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr);
787   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU",
788                     PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr);
789   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_ILU",
790                     PCFactorSetShiftInBlocks_ILU);CHKERRQ(ierr);
791   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
792                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
793   PetscFunctionReturn(0);
794 }
795 EXTERN_C_END
796