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