xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 40cbb1a031ea8f2be4fe2b92dc842b003ad37be3)
1 
2 /*
3    Defines a ILU factorization preconditioner for any Mat implementation
4 */
5 #include <../src/ksp/pc/impls/factor/ilu/ilu.h>     /*I "petscpc.h"  I*/
6 
7 PetscErrorCode  PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
8 {
9   PC_ILU *ilu = (PC_ILU*)pc->data;
10 
11   PetscFunctionBegin;
12   ilu->nonzerosalongdiagonal = PETSC_TRUE;
13   if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
14   else ilu->nonzerosalongdiagonaltol = z;
15   PetscFunctionReturn(0);
16 }
17 
18 PetscErrorCode PCReset_ILU(PC pc)
19 {
20   PC_ILU         *ilu = (PC_ILU*)pc->data;
21 
22   PetscFunctionBegin;
23   if (!ilu->hdr.inplace) PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact));
24   if (ilu->row && ilu->col && ilu->row != ilu->col) PetscCall(ISDestroy(&ilu->row));
25   PetscCall(ISDestroy(&ilu->col));
26   PetscFunctionReturn(0);
27 }
28 
29 PetscErrorCode  PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
30 {
31   PC_ILU *ilu = (PC_ILU*)pc->data;
32 
33   PetscFunctionBegin;
34   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
35     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
36   }
37   ((PC_Factor*)ilu)->info.dt      = dt;
38   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
39   ((PC_Factor*)ilu)->info.dtcount = dtcount;
40   ((PC_Factor*)ilu)->info.usedt   = 1.0;
41   PetscFunctionReturn(0);
42 }
43 
44 static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc)
45 {
46   PetscInt       itmp;
47   PetscBool      flg,set;
48   PC_ILU         *ilu = (PC_ILU*)pc->data;
49   PetscReal      tol;
50 
51   PetscFunctionBegin;
52   PetscOptionsHeadBegin(PetscOptionsObject,"ILU Options");
53   PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc));
54 
55   PetscCall(PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg));
56   if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
57 
58   PetscCall(PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set));
59   if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg;
60   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg));
61   if (flg) {
62     tol  = PETSC_DECIDE;
63     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,NULL));
64     PetscCall(PCFactorReorderForNonzeroDiagonal(pc,tol));
65   }
66 
67   PetscOptionsHeadEnd();
68   PetscFunctionReturn(0);
69 }
70 
71 static PetscErrorCode PCSetUp_ILU(PC pc)
72 {
73   PC_ILU                 *ilu = (PC_ILU*)pc->data;
74   MatInfo                info;
75   PetscBool              flg;
76   MatSolverType          stype;
77   MatFactorError         err;
78   const char             *prefix;
79 
80   PetscFunctionBegin;
81   if (!((PetscObject)pc->pmat)->prefix) {
82     PetscCall(PCGetOptionsPrefix(pc,&prefix));
83     PetscCall(MatSetOptionsPrefix(pc->pmat,prefix));
84   }
85   pc->failedreason = PC_NOERROR;
86   /* ugly hack to change default, since it is not support by some matrix types */
87   if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
88     PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg));
89     if (!flg) {
90       PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg));
91       if (!flg) {
92         ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
93         PetscCall(PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n"));
94       }
95     }
96   }
97 
98   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
99   if (ilu->hdr.inplace) {
100     if (!pc->setupcalled) {
101 
102       /* In-place factorization only makes sense with the natural ordering,
103          so we only need to get the ordering once, even if nonzero structure changes */
104       /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */
105       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
106       PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact));
107       PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col));
108       if (ilu->row) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row));
109       if (ilu->col) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col));
110     }
111 
112     /* In place ILU only makes sense with fill factor of 1.0 because
113        cannot have levels of fill */
114     ((PC_Factor*)ilu)->info.fill          = 1.0;
115     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
116 
117     PetscCall(MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info));
118     PetscCall(MatFactorGetError(pc->pmat,&err));
119     if (err) { /* Factor() fails */
120       pc->failedreason = (PCFailedReason)err;
121       PetscFunctionReturn(0);
122     }
123 
124     ((PC_Factor*)ilu)->fact = pc->pmat;
125     /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
126     PetscCall(PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate));
127   } else {
128     if (!pc->setupcalled) {
129       /* first time in so compute reordering and symbolic factorization */
130       PetscBool canuseordering;
131       if (!((PC_Factor*)ilu)->fact) {
132         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact));
133         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact));
134       }
135       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering));
136       if (canuseordering) {
137         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
138         PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col));
139         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row));
140         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col));
141         /*  Remove zeros along diagonal?     */
142         if (ilu->nonzerosalongdiagonal) {
143           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col));
144         }
145       }
146       PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info));
147       PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info));
148       ilu->hdr.actualfill = info.fill_ratio_needed;
149     } else if (pc->flag != SAME_NONZERO_PATTERN) {
150       if (!ilu->hdr.reuseordering) {
151         PetscBool canuseordering;
152         PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact));
153         PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact));
154         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact));
155         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering));
156         if (canuseordering) {
157           /* compute a new ordering for the ILU */
158           PetscCall(ISDestroy(&ilu->row));
159           PetscCall(ISDestroy(&ilu->col));
160           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
161           PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col));
162           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row));
163           PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col));
164           /*  Remove zeros along diagonal?     */
165           if (ilu->nonzerosalongdiagonal) {
166             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col));
167           }
168         }
169       }
170       PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info));
171       PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info));
172       ilu->hdr.actualfill = info.fill_ratio_needed;
173     }
174     PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err));
175     if (err) { /* FactorSymbolic() fails */
176       pc->failedreason = (PCFailedReason)err;
177       PetscFunctionReturn(0);
178     }
179 
180     PetscCall(MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info));
181     PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err));
182     if (err) { /* FactorNumeric() fails */
183       pc->failedreason = (PCFailedReason)err;
184     }
185   }
186 
187   PetscCall(PCFactorGetMatSolverType(pc,&stype));
188   if (!stype) {
189     MatSolverType solverpackage;
190     PetscCall(MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage));
191     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
192   }
193   PetscFunctionReturn(0);
194 }
195 
196 static PetscErrorCode PCDestroy_ILU(PC pc)
197 {
198   PC_ILU         *ilu = (PC_ILU*)pc->data;
199 
200   PetscFunctionBegin;
201   PetscCall(PCReset_ILU(pc));
202   PetscCall(PetscFree(((PC_Factor*)ilu)->solvertype));
203   PetscCall(PetscFree(((PC_Factor*)ilu)->ordering));
204   PetscCall(PetscFree(pc->data));
205   PetscCall(PCFactorClearComposedFunctions(pc));
206   PetscFunctionReturn(0);
207 }
208 
209 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
210 {
211   PC_ILU         *ilu = (PC_ILU*)pc->data;
212 
213   PetscFunctionBegin;
214   PetscCall(MatSolve(((PC_Factor*)ilu)->fact,x,y));
215   PetscFunctionReturn(0);
216 }
217 
218 static PetscErrorCode PCMatApply_ILU(PC pc,Mat X,Mat Y)
219 {
220   PC_ILU         *ilu = (PC_ILU*)pc->data;
221 
222   PetscFunctionBegin;
223   PetscCall(MatMatSolve(((PC_Factor*)ilu)->fact,X,Y));
224   PetscFunctionReturn(0);
225 }
226 
227 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
228 {
229   PC_ILU         *ilu = (PC_ILU*)pc->data;
230 
231   PetscFunctionBegin;
232   PetscCall(MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y));
233   PetscFunctionReturn(0);
234 }
235 
236 static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
237 {
238   PC_ILU         *icc = (PC_ILU*)pc->data;
239 
240   PetscFunctionBegin;
241   PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y));
242   PetscFunctionReturn(0);
243 }
244 
245 static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
246 {
247   PC_ILU         *icc = (PC_ILU*)pc->data;
248 
249   PetscFunctionBegin;
250   PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y));
251   PetscFunctionReturn(0);
252 }
253 
254 /*MC
255      PCILU - Incomplete factorization preconditioners.
256 
257    Options Database Keys:
258 +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
259 .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
260                       its factorization (overwrites original matrix)
261 .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
262 .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
263 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
264 .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
265                                    this decreases the chance of getting a zero pivot
266 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
267 -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
268                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
269                              stability of the ILU factorization
270 
271    Level: beginner
272 
273    Notes:
274     Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
275 
276           For BAIJ matrices this implements a point block ILU
277 
278           The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
279           even when the matrix is not symmetric since the U stores the diagonals of the factorization.
280 
281           If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARSE matrices with block Jacobi), factorization
282           is never done on the GPU).
283 
284    References:
285 +  * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
286    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
287 .  * -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
288 -  * -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
289       Chapter in Parallel Numerical
290       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
291       Science and Engineering, Kluwer.
292 
293 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`,
294           `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`,
295           `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
296           `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`,
297           `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()`
298 
299 M*/
300 
301 PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
302 {
303   PC_ILU         *ilu;
304 
305   PetscFunctionBegin;
306   PetscCall(PetscNewLog(pc,&ilu));
307   pc->data = (void*)ilu;
308   PetscCall(PCFactorInitialize(pc,MAT_FACTOR_ILU));
309 
310   ((PC_Factor*)ilu)->info.levels        = 0.;
311   ((PC_Factor*)ilu)->info.fill          = 1.0;
312   ilu->col                              = NULL;
313   ilu->row                              = NULL;
314   ((PC_Factor*)ilu)->info.dt            = PETSC_DEFAULT;
315   ((PC_Factor*)ilu)->info.dtcount       = PETSC_DEFAULT;
316   ((PC_Factor*)ilu)->info.dtcol         = PETSC_DEFAULT;
317 
318   pc->ops->reset               = PCReset_ILU;
319   pc->ops->destroy             = PCDestroy_ILU;
320   pc->ops->apply               = PCApply_ILU;
321   pc->ops->matapply            = PCMatApply_ILU;
322   pc->ops->applytranspose      = PCApplyTranspose_ILU;
323   pc->ops->setup               = PCSetUp_ILU;
324   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
325   pc->ops->view                = PCView_Factor;
326   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
327   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
328   pc->ops->applyrichardson     = NULL;
329   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU));
330   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU));
331   PetscFunctionReturn(0);
332 }
333