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