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