Actual source code: ilu.c
petsc-3.6.4 2016-04-12
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*/
7: /* ------------------------------------------------------------------------------------------*/
10: PetscErrorCode PCFactorSetReuseFill_ILU(PC pc,PetscBool flag)
11: {
12: PC_ILU *lu = (PC_ILU*)pc->data;
15: lu->reusefill = flag;
16: return(0);
17: }
21: PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
22: {
23: PC_ILU *ilu = (PC_ILU*)pc->data;
26: ilu->nonzerosalongdiagonal = PETSC_TRUE;
27: if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
28: else ilu->nonzerosalongdiagonaltol = z;
29: return(0);
30: }
34: PetscErrorCode PCReset_ILU(PC pc)
35: {
36: PC_ILU *ilu = (PC_ILU*)pc->data;
40: if (!ilu->inplace) {MatDestroy(&((PC_Factor*)ilu)->fact);}
41: if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(&ilu->row);}
42: ISDestroy(&ilu->col);
43: return(0);
44: }
48: PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
49: {
50: PC_ILU *ilu = (PC_ILU*)pc->data;
53: if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
54: SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
55: }
56: ((PC_Factor*)ilu)->info.dt = dt;
57: ((PC_Factor*)ilu)->info.dtcol = dtcol;
58: ((PC_Factor*)ilu)->info.dtcount = dtcount;
59: ((PC_Factor*)ilu)->info.usedt = 1.0;
60: return(0);
61: }
65: PetscErrorCode PCFactorSetReuseOrdering_ILU(PC pc,PetscBool flag)
66: {
67: PC_ILU *ilu = (PC_ILU*)pc->data;
70: ilu->reuseordering = flag;
71: return(0);
72: }
76: PetscErrorCode PCFactorSetUseInPlace_ILU(PC pc,PetscBool flg)
77: {
78: PC_ILU *dir = (PC_ILU*)pc->data;
81: dir->inplace = flg;
82: return(0);
83: }
87: PetscErrorCode PCFactorGetUseInPlace_ILU(PC pc,PetscBool *flg)
88: {
89: PC_ILU *dir = (PC_ILU*)pc->data;
92: *flg = dir->inplace;
93: return(0);
94: }
98: static PetscErrorCode PCSetFromOptions_ILU(PetscOptions *PetscOptionsObject,PC pc)
99: {
101: PetscInt itmp;
102: PetscBool flg,set;
103: PC_ILU *ilu = (PC_ILU*)pc->data;
104: PetscReal tol;
107: PetscOptionsHead(PetscOptionsObject,"ILU Options");
108: PCSetFromOptions_Factor(PetscOptionsObject,pc);
110: PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);
111: if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
113: PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
114: if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg;
115: PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
116: if (flg) {
117: tol = PETSC_DECIDE;
118: PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);
119: PCFactorReorderForNonzeroDiagonal(pc,tol);
120: }
122: PetscOptionsTail();
123: return(0);
124: }
128: static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
129: {
130: PC_ILU *ilu = (PC_ILU*)pc->data;
132: PetscBool iascii;
135: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
136: if (iascii) {
137: if (ilu->inplace) {
138: PetscViewerASCIIPrintf(viewer," ILU: in-place factorization\n");
139: } else {
140: PetscViewerASCIIPrintf(viewer," ILU: out-of-place factorization\n");
141: }
143: if (ilu->reusefill) {PetscViewerASCIIPrintf(viewer," ILU: Reusing fill from past factorization\n");}
144: if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer," ILU: Reusing reordering from past factorization\n");}
145: }
146: PCView_Factor(pc,viewer);
147: return(0);
148: }
152: static PetscErrorCode PCSetUp_ILU(PC pc)
153: {
155: PC_ILU *ilu = (PC_ILU*)pc->data;
156: MatInfo info;
157: PetscBool flg;
160: /* ugly hack to change default, since it is not support by some matrix types */
161: if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
162: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);
163: if (!flg) {
164: PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);
165: if (!flg) {
166: ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
167: PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");
168: }
169: }
170: }
172: if (ilu->inplace) {
173: if (!pc->setupcalled) {
175: /* In-place factorization only makes sense with the natural ordering,
176: so we only need to get the ordering once, even if nonzero structure changes */
177: MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
178: if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
179: if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
180: }
182: /* In place ILU only makes sense with fill factor of 1.0 because
183: cannot have levels of fill */
184: ((PC_Factor*)ilu)->info.fill = 1.0;
185: ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
187: MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
188: ((PC_Factor*)ilu)->fact = pc->pmat;
189: /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
190: PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);
191: } else {
192: if (!pc->setupcalled) {
193: /* first time in so compute reordering and symbolic factorization */
194: MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
195: if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
196: if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
197: /* Remove zeros along diagonal? */
198: if (ilu->nonzerosalongdiagonal) {
199: MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
200: }
201: if (!((PC_Factor*)ilu)->fact) {
202: MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
203: }
204: MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
205: MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
207: ilu->actualfill = info.fill_ratio_needed;
209: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);
210: } else if (pc->flag != SAME_NONZERO_PATTERN) {
211: if (!ilu->reuseordering) {
212: /* compute a new ordering for the ILU */
213: ISDestroy(&ilu->row);
214: ISDestroy(&ilu->col);
215: MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
216: if (ilu->row) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);}
217: if (ilu->col) {PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);}
218: /* Remove zeros along diagonal? */
219: if (ilu->nonzerosalongdiagonal) {
220: MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
221: }
222: }
223: MatDestroy(&((PC_Factor*)ilu)->fact);
224: MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
225: MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
226: MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
228: ilu->actualfill = info.fill_ratio_needed;
230: PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);
231: }
232: MatSetErrorIfFPE(pc->pmat,pc->erroriffailure);
233: MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);
234: }
235: return(0);
236: }
240: static PetscErrorCode PCDestroy_ILU(PC pc)
241: {
242: PC_ILU *ilu = (PC_ILU*)pc->data;
246: PCReset_ILU(pc);
247: PetscFree(((PC_Factor*)ilu)->solvertype);
248: PetscFree(((PC_Factor*)ilu)->ordering);
249: PetscFree(pc->data);
250: return(0);
251: }
255: static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
256: {
257: PC_ILU *ilu = (PC_ILU*)pc->data;
261: MatSolve(((PC_Factor*)ilu)->fact,x,y);
262: return(0);
263: }
267: static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
268: {
269: PC_ILU *ilu = (PC_ILU*)pc->data;
273: MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);
274: return(0);
275: }
279: static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
280: {
282: PC_ILU *icc = (PC_ILU*)pc->data;
285: MatForwardSolve(((PC_Factor*)icc)->fact,x,y);
286: return(0);
287: }
291: static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
292: {
294: PC_ILU *icc = (PC_ILU*)pc->data;
297: MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);
298: return(0);
299: }
301: /*MC
302: PCILU - Incomplete factorization preconditioners.
304: Options Database Keys:
305: + -pc_factor_levels <k> - number of levels of fill for ILU(k)
306: . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
307: its factorization (overwrites original matrix)
308: . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
309: . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
310: . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
311: . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
312: this decreases the chance of getting a zero pivot
313: . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
314: - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
315: than 1 the diagonal blocks are factored with partial pivoting (this increases the
316: stability of the ILU factorization
318: Level: beginner
320: Concepts: incomplete factorization
322: Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
324: For BAIJ matrices this implements a point block ILU
326: The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
327: even when the matrix is not symmetric since the U stores the diagonals of the factorization.
329: If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization
330: is never done on the GPU).
332: References:
333: T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
334: self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
336: T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
337: fusion problems. Quart. Appl. Math., 19:221--229, 1961.
339: Review article:
340: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
341: http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
342: Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
343: Science and Engineering, Kluwer, pp. 167--202.
346: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
347: PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
348: PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
349: PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
350: PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace()
352: M*/
356: PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
357: {
359: PC_ILU *ilu;
362: PetscNewLog(pc,&ilu);
364: ((PC_Factor*)ilu)->fact = 0;
365: MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);
366: ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU;
367: ((PC_Factor*)ilu)->info.levels = 0.;
368: ((PC_Factor*)ilu)->info.fill = 1.0;
369: ilu->col = 0;
370: ilu->row = 0;
371: ilu->inplace = PETSC_FALSE;
372: PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);
373: PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);
374: ilu->reuseordering = PETSC_FALSE;
375: ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT;
376: ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT;
377: ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT;
378: ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
379: ((PC_Factor*)ilu)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
380: ((PC_Factor*)ilu)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON;
381: ((PC_Factor*)ilu)->info.pivotinblocks = 1.0;
382: ilu->reusefill = PETSC_FALSE;
383: ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
384: pc->data = (void*)ilu;
386: pc->ops->reset = PCReset_ILU;
387: pc->ops->destroy = PCDestroy_ILU;
388: pc->ops->apply = PCApply_ILU;
389: pc->ops->applytranspose = PCApplyTranspose_ILU;
390: pc->ops->setup = PCSetUp_ILU;
391: pc->ops->setfromoptions = PCSetFromOptions_ILU;
392: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
393: pc->ops->view = PCView_ILU;
394: pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU;
395: pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
396: pc->ops->applyrichardson = 0;
398: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
399: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
400: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
401: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);
402: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);
403: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);
404: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);
405: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);
406: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);
407: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_ILU);
408: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_ILU);
409: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);
410: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);
411: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_ILU);
412: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ILU);
413: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);
414: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);
415: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);
416: PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);
417: return(0);
418: }