xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 71fd2e92b2fee25f70818ffc1dca9fc57f8584d7)
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__ "PCFactorReorderForNonzeroDiagonal_ILU"
25 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
26 {
27   PC_ILU *ilu = (PC_ILU*)pc->data;
28 
29   PetscFunctionBegin;
30   ilu->nonzerosalongdiagonal = PETSC_TRUE;
31   if (z == PETSC_DECIDE) {
32     ilu->nonzerosalongdiagonaltol = 1.e-10;
33   } else {
34     ilu->nonzerosalongdiagonaltol = z;
35   }
36   PetscFunctionReturn(0);
37 }
38 EXTERN_C_END
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "PCDestroy_ILU_Internal"
42 PetscErrorCode PCDestroy_ILU_Internal(PC pc)
43 {
44   PC_ILU         *ilu = (PC_ILU*)pc->data;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);}
49   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
50   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
51   PetscFunctionReturn(0);
52 }
53 
54 EXTERN_C_BEGIN
55 #undef __FUNCT__
56 #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU"
57 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
58 {
59   PC_ILU         *ilu = (PC_ILU*)pc->data;
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   if (pc->setupcalled && (!ilu->usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
64     pc->setupcalled   = 0;
65     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
66   }
67   ilu->usedt                      = PETSC_TRUE;
68   ((PC_Factor*)ilu)->info.dt      = dt;
69   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
70   ((PC_Factor*)ilu)->info.dtcount = dtcount;
71   ((PC_Factor*)ilu)->info.fill    = PETSC_DEFAULT;
72   PetscFunctionReturn(0);
73 }
74 EXTERN_C_END
75 
76 EXTERN_C_BEGIN
77 #undef __FUNCT__
78 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
79 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
80 {
81   PC_ILU *ilu = (PC_ILU*)pc->data;
82 
83   PetscFunctionBegin;
84   ilu->reuseordering = flag;
85   PetscFunctionReturn(0);
86 }
87 EXTERN_C_END
88 
89 EXTERN_C_BEGIN
90 #undef __FUNCT__
91 #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
92 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc)
93 {
94   PC_ILU *dir = (PC_ILU*)pc->data;
95 
96   PetscFunctionBegin;
97   dir->inplace = PETSC_TRUE;
98   PetscFunctionReturn(0);
99 }
100 EXTERN_C_END
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "PCSetFromOptions_ILU"
104 static PetscErrorCode PCSetFromOptions_ILU(PC pc)
105 {
106   PetscErrorCode ierr;
107   PetscInt       dtmax = 3,itmp;
108   PetscTruth     flg;
109   PetscReal      dt[3];
110   PC_ILU         *ilu = (PC_ILU*)pc->data;
111   PetscReal      tol;
112 
113   PetscFunctionBegin;
114   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
115     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
116 
117     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
118     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
119     flg  = PETSC_FALSE;
120     ierr = PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
121     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
122 
123     dt[0] = ((PC_Factor*)ilu)->info.dt;
124     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
125     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
126     ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
127     if (flg) {
128       ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
129     }
130 
131     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
132     if (flg) {
133       tol = PETSC_DECIDE;
134       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
135       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
136     }
137 
138   ierr = PetscOptionsTail();CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "PCView_ILU"
144 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
145 {
146   PC_ILU         *ilu = (PC_ILU*)pc->data;
147   PetscErrorCode ierr;
148   PetscTruth     isstring,iascii;
149 
150   PetscFunctionBegin;
151   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
152   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
153   if (iascii) {
154     if (ilu->usedt) {
155         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %G\n",((PC_Factor*)ilu)->info.dt);CHKERRQ(ierr);
156         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)((PC_Factor*)ilu)->info.dtcount);CHKERRQ(ierr);
157         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %G\n",((PC_Factor*)ilu)->info.dtcol);CHKERRQ(ierr);
158     } else if (((PC_Factor*)ilu)->info.levels == 1) {
159         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr);
160     } else {
161         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr);
162     }
163     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio allocated %G\n",((PC_Factor*)ilu)->info.fill);CHKERRQ(ierr);
164     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %G\n",((PC_Factor*)ilu)->info.zeropivot);CHKERRQ(ierr);
165     if (((PC_Factor*)ilu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
166     if (((PC_Factor*)ilu)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
167     if (((PC_Factor*)ilu)->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);}
168     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
169     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
170     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
171     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
172     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
173     if (((PC_Factor*)ilu)->fact) {
174       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr);
175       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
176       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
177       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
178       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
179       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
180       ierr = MatView(((PC_Factor*)ilu)->fact,viewer);CHKERRQ(ierr);
181       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
182       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
183       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
184       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
185     }
186   } else if (isstring) {
187     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
188   } else {
189     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
190   }
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "PCSetUp_ILU"
196 static PetscErrorCode PCSetUp_ILU(PC pc)
197 {
198   PetscErrorCode ierr;
199   PC_ILU         *ilu = (PC_ILU*)pc->data;
200   MatInfo        info;
201 
202   PetscFunctionBegin;
203   if (ilu->inplace) {
204     CHKMEMQ;
205     if (!pc->setupcalled) {
206 
207       /* In-place factorization only makes sense with the natural ordering,
208          so we only need to get the ordering once, even if nonzero structure changes */
209       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
210       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
211       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
212     }
213 
214     /* In place ILU only makes sense with fill factor of 1.0 because
215        cannot have levels of fill */
216     ((PC_Factor*)ilu)->info.fill          = 1.0;
217     ((PC_Factor*)ilu)->info.diagonal_fill = 0;
218     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
219     CHKMEMQ;
220     ((PC_Factor*)ilu)->fact = pc->pmat;
221   } else if (ilu->usedt) {
222     if (!pc->setupcalled) {
223       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
224     CHKMEMQ;
225       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
226       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
227       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
228       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
229       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
230     } else if (pc->flag != SAME_NONZERO_PATTERN) {
231     CHKMEMQ;
232       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
233     CHKMEMQ;
234       if (!ilu->reuseordering) {
235         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
236         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
237         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
238         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
239         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
240       }
241       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
242       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
243       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
244     } else if (!ilu->reusefill) {
245       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
246       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
247       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
248       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
249     } else {
250       ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
251     }
252   } else {
253     if (!pc->setupcalled) {
254       /* first time in so compute reordering and symbolic factorization */
255       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
256       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
257       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
258       /*  Remove zeros along diagonal?     */
259       if (ilu->nonzerosalongdiagonal) {
260         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
261       }
262     CHKMEMQ;
263       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
264     CHKMEMQ;
265       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
266       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
267       ilu->actualfill = info.fill_ratio_needed;
268       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
269     } else if (pc->flag != SAME_NONZERO_PATTERN) {
270       if (!ilu->reuseordering) {
271         /* compute a new ordering for the ILU */
272         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
273         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
274         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
275         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
276         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
277         /*  Remove zeros along diagonal?     */
278         if (ilu->nonzerosalongdiagonal) {
279           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
280         }
281       }
282       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
283       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
284       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
285       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
286       ilu->actualfill = info.fill_ratio_needed;
287       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
288     }
289     CHKMEMQ;
290     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
291     CHKMEMQ;
292   }
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "PCDestroy_ILU"
298 static PetscErrorCode PCDestroy_ILU(PC pc)
299 {
300   PC_ILU         *ilu = (PC_ILU*)pc->data;
301   PetscErrorCode ierr;
302 
303   PetscFunctionBegin;
304   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
305   ierr = PetscStrfree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
306   ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
307   ierr = PetscFree(ilu);CHKERRQ(ierr);
308   PetscFunctionReturn(0);
309 }
310 
311 #undef __FUNCT__
312 #define __FUNCT__ "PCApply_ILU"
313 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
314 {
315   PC_ILU         *ilu = (PC_ILU*)pc->data;
316   PetscErrorCode ierr;
317 
318   PetscFunctionBegin;
319   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
320   PetscFunctionReturn(0);
321 }
322 
323 #undef __FUNCT__
324 #define __FUNCT__ "PCApplyTranspose_ILU"
325 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
326 {
327   PC_ILU         *ilu = (PC_ILU*)pc->data;
328   PetscErrorCode ierr;
329 
330   PetscFunctionBegin;
331   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
332   PetscFunctionReturn(0);
333 }
334 
335 /*MC
336      PCILU - Incomplete factorization preconditioners.
337 
338    Options Database Keys:
339 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
340 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
341                       its factorization (overwrites original matrix)
342 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
343 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
344 .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
345 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
346 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
347                                    this decreases the chance of getting a zero pivot
348 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
349 .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
350                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
351                              stability of the ILU factorization
352 .  -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization
353 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
354 -  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
355    is optional with true being the default
356 
357    Level: beginner
358 
359   Concepts: incomplete factorization
360 
361    Notes: Only implemented for some matrix formats. (for parallel use you
362              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended
363              unless you really want a parallel ILU).
364 
365           For BAIJ matrices this implements a point block ILU
366 
367    References:
368    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
369    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
370 
371    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
372    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
373 
374    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
375       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
376       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
377       Science and Engineering, Kluwer, pp. 167--202.
378 
379 
380 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
381            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks(),
382            PCFactorSetUseDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
383            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks()
384 
385 M*/
386 
387 EXTERN_C_BEGIN
388 #undef __FUNCT__
389 #define __FUNCT__ "PCCreate_ILU"
390 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
391 {
392   PetscErrorCode ierr;
393   PC_ILU         *ilu;
394 
395   PetscFunctionBegin;
396   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
397 
398   ((PC_Factor*)ilu)->fact                    = 0;
399   ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
400   ((PC_Factor*)ilu)->info.levels             = 0;
401   ((PC_Factor*)ilu)->info.fill               = 1.0;
402   ilu->col                     = 0;
403   ilu->row                     = 0;
404   ilu->inplace                 = PETSC_FALSE;
405   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
406   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
407   ilu->reuseordering           = PETSC_FALSE;
408   ilu->usedt                   = PETSC_FALSE;
409   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
410   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
411   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
412   ((PC_Factor*)ilu)->info.shiftnz            = 1.e-12;
413   ((PC_Factor*)ilu)->info.shiftpd            = 0.0; /* false */
414   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
415   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
416   ((PC_Factor*)ilu)->info.shiftinblocks      = 1.e-12;
417   ilu->reusefill               = PETSC_FALSE;
418   ((PC_Factor*)ilu)->info.diagonal_fill      = 0;
419   pc->data                     = (void*)ilu;
420 
421   pc->ops->destroy             = PCDestroy_ILU;
422   pc->ops->apply               = PCApply_ILU;
423   pc->ops->applytranspose      = PCApplyTranspose_ILU;
424   pc->ops->setup               = PCSetUp_ILU;
425   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
426   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
427   pc->ops->view                = PCView_ILU;
428   pc->ops->applyrichardson     = 0;
429 
430   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
431                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
432   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
433                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
434   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
435                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
436   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
437                     PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr);
438 
439   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
440                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
441   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
442                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
443   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
444                     PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr);
445   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
446                     PCFactorSetFill_Factor);CHKERRQ(ierr);
447   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
448                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
449   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
450                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
451   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
452                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
453   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
454                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
455   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
456                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
457   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
458                     PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
459   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
460                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
461   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
462                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 EXTERN_C_END
466