xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a ILU factorization preconditioner for any Mat implementation
5 */
6 #include "private/pcimpl.h"                 /*I "petscpc.h"  I*/
7 #include "src/ksp/pc/impls/factor/ilu/ilu.h"
8 #include "src/mat/matimpl.h"
9 
10 /* ------------------------------------------------------------------------------------------*/
11 EXTERN_C_BEGIN
12 #undef __FUNCT__
13 #define __FUNCT__ "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__ "PCFactorSetMatOrdering_ILU"
148 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering_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 .keywords: PC, levels, reordering, factorization, incomplete, ILU
261 @*/
262 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
263 {
264   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
268   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
269   if (f) {
270     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
271   }
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "PCFactorSetLevels"
277 /*@
278    PCFactorSetLevels - Sets the number of levels of fill to use.
279 
280    Collective on PC
281 
282    Input Parameters:
283 +  pc - the preconditioner context
284 -  levels - number of levels of fill
285 
286    Options Database Key:
287 .  -pc_factor_levels <levels> - Sets fill level
288 
289    Level: intermediate
290 
291 .keywords: PC, levels, fill, factorization, incomplete, ILU
292 @*/
293 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
294 {
295   PetscErrorCode ierr,(*f)(PC,PetscInt);
296 
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
299   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
300   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
301   if (f) {
302     ierr = (*f)(pc,levels);CHKERRQ(ierr);
303   }
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
309 /*@
310    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
311    treated as level 0 fill even if there is no non-zero location.
312 
313    Collective on PC
314 
315    Input Parameters:
316 +  pc - the preconditioner context
317 
318    Options Database Key:
319 .  -pc_factor_diagonal_fill
320 
321    Notes:
322    Does not apply with 0 fill.
323 
324    Level: intermediate
325 
326 .keywords: PC, levels, fill, factorization, incomplete, ILU
327 @*/
328 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
329 {
330   PetscErrorCode ierr,(*f)(PC);
331 
332   PetscFunctionBegin;
333   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
334   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
335   if (f) {
336     ierr = (*f)(pc);CHKERRQ(ierr);
337   }
338   PetscFunctionReturn(0);
339 }
340 
341 /* ------------------------------------------------------------------------------------------*/
342 
343 EXTERN_C_BEGIN
344 #undef __FUNCT__
345 #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU"
346 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
347 {
348   PC_ILU *dir = (PC_ILU*)pc->data;
349 
350   PetscFunctionBegin;
351   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
352   PetscFunctionReturn(0);
353 }
354 EXTERN_C_END
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "PCSetFromOptions_ILU"
358 static PetscErrorCode PCSetFromOptions_ILU(PC pc)
359 {
360   PetscErrorCode ierr;
361   PetscInt       dtmax = 3,itmp;
362   PetscTruth     flg,set;
363   PetscReal      dt[3];
364   char           tname[256];
365   PC_ILU         *ilu = (PC_ILU*)pc->data;
366   PetscFList     ordlist;
367   PetscReal      tol;
368 
369   PetscFunctionBegin;
370   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
371   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
372     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr);
373     if (flg) ilu->info.levels = itmp;
374     ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
375     if (flg) ilu->inplace = PETSC_TRUE;
376     ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr);
377     ilu->info.diagonal_fill = (double) flg;
378     ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
379     if (flg) ilu->reusefill = PETSC_TRUE;
380     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
381     if (flg) ilu->reuseordering = PETSC_TRUE;
382     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
383     if (flg) {
384       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
385     }
386     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr);
387 
388     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
389     if (flg) {
390       ierr = PetscOptionsInt("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr);
391       if (flg && !itmp) {
392 	ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr);
393       } else {
394 	ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
395       }
396     }
397     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr);
398 
399     dt[0] = ilu->info.dt;
400     dt[1] = ilu->info.dtcol;
401     dt[2] = 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",ilu->info.fill,&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","PCFactorSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr);
416     if (flg) {
417       ierr = PCFactorSetMatOrdering(pc,tname);CHKERRQ(ierr);
418     }
419     flg = 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 = PetscOptionsTail();CHKERRQ(ierr);
425   PetscFunctionReturn(0);
426 }
427 
428 #undef __FUNCT__
429 #define __FUNCT__ "PCView_ILU"
430 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
431 {
432   PC_ILU         *ilu = (PC_ILU*)pc->data;
433   PetscErrorCode ierr;
434   PetscTruth     isstring,iascii;
435 
436   PetscFunctionBegin;
437   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
438   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
439   if (iascii) {
440     if (ilu->usedt) {
441         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr);
442         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr);
443         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr);
444     } else if (ilu->info.levels == 1) {
445         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
446     } else {
447         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
448     }
449     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr);
450     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr);
451     if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
452     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
453     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
454     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr);
455     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
456     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
457     if (ilu->fact) {
458       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr);
459       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
460       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
461       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
462       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
463       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
464       ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr);
465       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
466       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
467       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
468       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
469     }
470   } else if (isstring) {
471     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
472   } else {
473     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
474   }
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "PCSetUp_ILU"
480 static PetscErrorCode PCSetUp_ILU(PC pc)
481 {
482   PetscErrorCode ierr;
483   PC_ILU         *ilu = (PC_ILU*)pc->data;
484   MatInfo        info;
485 
486   PetscFunctionBegin;
487   if (ilu->inplace) {
488     if (!pc->setupcalled) {
489 
490       /* In-place factorization only makes sense with the natural ordering,
491          so we only need to get the ordering once, even if nonzero structure changes */
492       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
493       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
494       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
495     }
496 
497     /* In place ILU only makes sense with fill factor of 1.0 because
498        cannot have levels of fill */
499     ilu->info.fill          = 1.0;
500     ilu->info.diagonal_fill = 0;
501     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
502     ilu->fact = pc->pmat;
503   } else if (ilu->usedt) {
504     if (!pc->setupcalled) {
505       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
506       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
507       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
508       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
509       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
510     } else if (pc->flag != SAME_NONZERO_PATTERN) {
511       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
512       if (!ilu->reuseordering) {
513         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
514         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
515         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
516         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
517         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
518       }
519       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
520       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
521     } else if (!ilu->reusefill) {
522       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
523       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
524       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
525     } else {
526       ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
527     }
528   } else {
529     if (!pc->setupcalled) {
530       /* first time in so compute reordering and symbolic factorization */
531       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
532       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
533       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
534       /*  Remove zeros along diagonal?     */
535       if (ilu->nonzerosalongdiagonal) {
536         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
537       }
538       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
539       ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
540       ilu->actualfill = info.fill_ratio_needed;
541       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
542     } else if (pc->flag != SAME_NONZERO_PATTERN) {
543       if (!ilu->reuseordering) {
544         /* compute a new ordering for the ILU */
545         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
546         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
547         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
548         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
549         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
550         /*  Remove zeros along diagonal?     */
551         if (ilu->nonzerosalongdiagonal) {
552           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
553         }
554       }
555       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
556       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
557       ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
558       ilu->actualfill = info.fill_ratio_needed;
559       ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr);
560     }
561     ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
562   }
563   PetscFunctionReturn(0);
564 }
565 
566 #undef __FUNCT__
567 #define __FUNCT__ "PCDestroy_ILU"
568 static PetscErrorCode PCDestroy_ILU(PC pc)
569 {
570   PC_ILU         *ilu = (PC_ILU*)pc->data;
571   PetscErrorCode ierr;
572 
573   PetscFunctionBegin;
574   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
575   ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr);
576   ierr = PetscFree(ilu);CHKERRQ(ierr);
577   PetscFunctionReturn(0);
578 }
579 
580 #undef __FUNCT__
581 #define __FUNCT__ "PCApply_ILU"
582 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
583 {
584   PC_ILU         *ilu = (PC_ILU*)pc->data;
585   PetscErrorCode ierr;
586 
587   PetscFunctionBegin;
588   ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 #undef __FUNCT__
593 #define __FUNCT__ "PCApplyTranspose_ILU"
594 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
595 {
596   PC_ILU         *ilu = (PC_ILU*)pc->data;
597   PetscErrorCode ierr;
598 
599   PetscFunctionBegin;
600   ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr);
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "PCGetFactoredMatrix_ILU"
606 static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
607 {
608   PC_ILU *ilu = (PC_ILU*)pc->data;
609 
610   PetscFunctionBegin;
611   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
612   *mat = ilu->fact;
613   PetscFunctionReturn(0);
614 }
615 
616 /*MC
617      PCILU - Incomplete factorization preconditioners.
618 
619    Options Database Keys:
620 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
621 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
622                       its factorization (overwrites original matrix)
623 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
624 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
625 .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
626 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
627 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
628                                    this decreases the chance of getting a zero pivot
629 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
630 .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
631                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
632                              stability of the ILU factorization
633 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
634 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
635    is optional with PETSC_TRUE being the default
636 
637    Level: beginner
638 
639   Concepts: incomplete factorization
640 
641    Notes: Only implemented for some matrix formats. Not implemented in parallel
642 
643           For BAIJ matrices this implements a point block ILU
644 
645 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
646            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(),
647            PCFactorSetFill(), PCFactorSetMatOrdering(), PCFactorSetReuseOrdering(),
648            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
649            PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
650 
651 M*/
652 
653 EXTERN_C_BEGIN
654 #undef __FUNCT__
655 #define __FUNCT__ "PCCreate_ILU"
656 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
657 {
658   PetscErrorCode ierr;
659   PC_ILU         *ilu;
660 
661   PetscFunctionBegin;
662   ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr);
663   ierr = PetscLogObjectMemory(pc,sizeof(PC_ILU));CHKERRQ(ierr);
664 
665   ilu->fact                    = 0;
666   ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr);
667   ilu->info.levels             = 0;
668   ilu->info.fill               = 1.0;
669   ilu->col                     = 0;
670   ilu->row                     = 0;
671   ilu->inplace                 = PETSC_FALSE;
672   ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr);
673   ilu->reuseordering           = PETSC_FALSE;
674   ilu->usedt                   = PETSC_FALSE;
675   ilu->info.dt                 = PETSC_DEFAULT;
676   ilu->info.dtcount            = PETSC_DEFAULT;
677   ilu->info.dtcol              = PETSC_DEFAULT;
678   ilu->info.shiftnz            = 0.0;
679   ilu->info.shiftpd            = 0.0; /* false */
680   ilu->info.shift_fraction     = 0.0;
681   ilu->info.zeropivot          = 1.e-12;
682   ilu->info.pivotinblocks      = 1.0;
683   ilu->reusefill               = PETSC_FALSE;
684   ilu->info.diagonal_fill      = 0;
685   pc->data                     = (void*)ilu;
686 
687   pc->ops->destroy             = PCDestroy_ILU;
688   pc->ops->apply               = PCApply_ILU;
689   pc->ops->applytranspose      = PCApplyTranspose_ILU;
690   pc->ops->setup               = PCSetUp_ILU;
691   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
692   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ILU;
693   pc->ops->view                = PCView_ILU;
694   pc->ops->applyrichardson     = 0;
695 
696   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU",
697                     PCFactorSetZeroPivot_ILU);CHKERRQ(ierr);
698   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU",
699                     PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr);
700   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU",
701                     PCFactorSetShiftPd_ILU);CHKERRQ(ierr);
702 
703   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
704                     PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr);
705   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU",
706                     PCFactorSetFill_ILU);CHKERRQ(ierr);
707   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrdering_C","PCFactorSetMatOrdering_ILU",
708                     PCFactorSetMatOrdering_ILU);CHKERRQ(ierr);
709   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
710                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
711   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
712                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
713   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU",
714                     PCFactorSetLevels_ILU);CHKERRQ(ierr);
715   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
716                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
717   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU",
718                     PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr);
719   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU",
720                     PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr);
721   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
722                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
723   PetscFunctionReturn(0);
724 }
725 EXTERN_C_END
726