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