Actual source code: ilu.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    Defines a ILU factorization preconditioner for any Mat implementation
  5: */
 6:  #include ../src/ksp/pc/impls/factor/ilu/ilu.h

  8: /* ------------------------------------------------------------------------------------------*/
 12: PetscErrorCode  PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag)
 13: {
 14:   PC_ILU *lu = (PC_ILU*)pc->data;
 15: 
 17:   lu->reusefill = flag;
 18:   return(0);
 19: }

 25: PetscErrorCode  PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
 26: {
 27:   PC_ILU *ilu = (PC_ILU*)pc->data;

 30:   ilu->nonzerosalongdiagonal = PETSC_TRUE;
 31:   if (z == PETSC_DECIDE) {
 32:     ilu->nonzerosalongdiagonaltol = 1.e-10;
 33:   } else {
 34:     ilu->nonzerosalongdiagonaltol = z;
 35:   }
 36:   return(0);
 37: }

 42: PetscErrorCode PCDestroy_ILU_Internal(PC pc)
 43: {
 44:   PC_ILU         *ilu = (PC_ILU*)pc->data;

 48:   if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {MatDestroy(((PC_Factor*)ilu)->fact);}
 49:   if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(ilu->row);}
 50:   if (ilu->col) {ISDestroy(ilu->col);}
 51:   return(0);
 52: }

 57: PetscErrorCode  PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 58: {
 59:   PC_ILU         *ilu = (PC_ILU*)pc->data;

 62:   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 63:     SETERRQ(PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
 64:   }
 65:   ((PC_Factor*)ilu)->info.dt      = dt;
 66:   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
 67:   ((PC_Factor*)ilu)->info.dtcount = dtcount;
 68:   ((PC_Factor*)ilu)->info.usedt   = 1.0;
 69:   return(0);
 70: }

 76: PetscErrorCode  PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
 77: {
 78:   PC_ILU *ilu = (PC_ILU*)pc->data;

 81:   ilu->reuseordering = flag;
 82:   return(0);
 83: }

 89: PetscErrorCode  PCFactorSetUseInPlace_ILU(PC pc)
 90: {
 91:   PC_ILU *dir = (PC_ILU*)pc->data;

 94:   dir->inplace = PETSC_TRUE;
 95:   return(0);
 96: }

101: static PetscErrorCode PCSetFromOptions_ILU(PC pc)
102: {
104:   PetscInt       dtmax = 3,itmp;
105:   PetscTruth     flg;
106:   PetscReal      dt[3];
107:   PC_ILU         *ilu = (PC_ILU*)pc->data;
108:   PetscReal      tol;

111:   PetscOptionsHead("ILU Options");
112:     PCSetFromOptions_Factor(pc);

114:     PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);
115:     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
116:     flg  = PETSC_FALSE;
117:     PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);
118:     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
119: 
120:     dt[0] = ((PC_Factor*)ilu)->info.dt;
121:     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
122:     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
123:     PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);
124:     if (flg) {
125:       PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);
126:     }

128:     PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
129:     if (flg) {
130:       tol = PETSC_DECIDE;
131:       PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);
132:       PCFactorReorderForNonzeroDiagonal(pc,tol);
133:     }

135:   PetscOptionsTail();
136:   return(0);
137: }

141: static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
142: {
143:   PC_ILU         *ilu = (PC_ILU*)pc->data;
145:   PetscTruth     iascii;

148:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
149:   if (iascii) {
150:     if (ilu->inplace) {
151:       PetscViewerASCIIPrintf(viewer,"  ILU: in-place factorization\n");
152:     } else {
153:       PetscViewerASCIIPrintf(viewer,"  ILU: out-of-place factorization\n");
154:     }
155: 
156:     if (ilu->reusefill)     {PetscViewerASCIIPrintf(viewer,"  ILU: Reusing fill from past factorization\n");}
157:     if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer,"  ILU: Reusing reordering from past factorization\n");}
158:   }
159:   PCView_Factor(pc,viewer);
160:   return(0);
161: }

165: static PetscErrorCode PCSetUp_ILU(PC pc)
166: {
168:   PC_ILU         *ilu = (PC_ILU*)pc->data;
169:   MatInfo        info;

172:   if (ilu->inplace) {
173:     CHKMEMQ;
174:     if (!pc->setupcalled) {

176:       /* In-place factorization only makes sense with the natural ordering,
177:          so we only need to get the ordering once, even if nonzero structure changes */
178:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
179:       if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
180:       if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
181:     }

183:     /* In place ILU only makes sense with fill factor of 1.0 because 
184:        cannot have levels of fill */
185:     ((PC_Factor*)ilu)->info.fill          = 1.0;
186:     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
187:     MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKMEMQ;
188:     ((PC_Factor*)ilu)->fact = pc->pmat;
189:   } else {
190:     if (!pc->setupcalled) {
191:       /* first time in so compute reordering and symbolic factorization */
192:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
193:       if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
194:       if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
195:       /*  Remove zeros along diagonal?     */
196:       if (ilu->nonzerosalongdiagonal) {
197:         MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
198:       }
199:     CHKMEMQ;
200:       MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
201:     CHKMEMQ;
202:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
203:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
204:       ilu->actualfill = info.fill_ratio_needed;
205:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
206:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
207:       if (!ilu->reuseordering) {
208:         /* compute a new ordering for the ILU */
209:         ISDestroy(ilu->row);
210:         ISDestroy(ilu->col);
211:         MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
212:         if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
213:         if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
214:         /*  Remove zeros along diagonal?     */
215:         if (ilu->nonzerosalongdiagonal) {
216:           MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
217:         }
218:       }
219:       MatDestroy(((PC_Factor*)ilu)->fact);
220:       MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
221:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
222:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
223:       ilu->actualfill = info.fill_ratio_needed;
224:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
225:     }
226:     CHKMEMQ;
227:     MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);
228:     CHKMEMQ;
229:   }
230:   return(0);
231: }

235: static PetscErrorCode PCDestroy_ILU(PC pc)
236: {
237:   PC_ILU         *ilu = (PC_ILU*)pc->data;

241:   PCDestroy_ILU_Internal(pc);
242:   PetscStrfree(((PC_Factor*)ilu)->solvertype);
243:   PetscStrfree(((PC_Factor*)ilu)->ordering);
244:   PetscFree(ilu);
245:   return(0);
246: }

250: static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
251: {
252:   PC_ILU         *ilu = (PC_ILU*)pc->data;

256:   MatSolve(((PC_Factor*)ilu)->fact,x,y);
257:   return(0);
258: }

262: static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
263: {
264:   PC_ILU         *ilu = (PC_ILU*)pc->data;

268:   MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);
269:   return(0);
270: }

272: /*MC
273:      PCILU - Incomplete factorization preconditioners.

275:    Options Database Keys:
276: +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
277: .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
278:                       its factorization (overwrites original matrix)
279: .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
280: .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
281: .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - use drop tolerance factorization
282: .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
283: .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
284:                                    this decreases the chance of getting a zero pivot
285: .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
286: .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
287:                              than 1 the diagonal blocks are factored with partial pivoting (this increases the 
288:                              stability of the ILU factorization

290:    Level: beginner

292:   Concepts: incomplete factorization

294:    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)

296:           For BAIJ matrices this implements a point block ILU

298:    References:
299:    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
300:    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.

302:    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
303:    fusion problems. Quart. Appl. Math., 19:221--229, 1961.

305:    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
306:       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
307:       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
308:       Science and Engineering, Kluwer, pp. 167--202.


311: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
312:            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
313:            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
314:            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks()

316: M*/

321: PetscErrorCode  PCCreate_ILU(PC pc)
322: {
324:   PC_ILU         *ilu;

327:   PetscNewLog(pc,PC_ILU,&ilu);

329:   ((PC_Factor*)ilu)->fact                    = 0;
330:   MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);
331:   ((PC_Factor*)ilu)->info.levels             = 0.;
332:   ((PC_Factor*)ilu)->info.fill               = 1.0;
333:   ilu->col                     = 0;
334:   ilu->row                     = 0;
335:   ilu->inplace                 = PETSC_FALSE;
336:   PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)ilu)->solvertype);
337:   PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);
338:   ilu->reuseordering           = PETSC_FALSE;
339:   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
340:   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
341:   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
342:   ((PC_Factor*)ilu)->info.shifttype          = (PetscReal)MAT_SHIFT_NONZERO;
343:   ((PC_Factor*)ilu)->info.shiftamount        = 1.e-12;
344:   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
345:   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
346:   ilu->reusefill               = PETSC_FALSE;
347:   ((PC_Factor*)ilu)->info.diagonal_fill      = 0.;
348:   pc->data                     = (void*)ilu;

350:   pc->ops->destroy             = PCDestroy_ILU;
351:   pc->ops->apply               = PCApply_ILU;
352:   pc->ops->applytranspose      = PCApplyTranspose_ILU;
353:   pc->ops->setup               = PCSetUp_ILU;
354:   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
355:   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
356:   pc->ops->view                = PCView_ILU;
357:   pc->ops->applyrichardson     = 0;

359:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
360:                     PCFactorSetZeroPivot_Factor);
361:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
362:                     PCFactorSetShiftType_Factor);
363:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
364:                     PCFactorSetShiftAmount_Factor);
365:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
366:                     PCFactorGetMatSolverPackage_Factor);
367:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
368:                     PCFactorSetMatSolverPackage_Factor);
369:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
370:                     PCFactorSetDropTolerance_ILU);
371:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
372:                     PCFactorSetFill_Factor);
373:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
374:                     PCFactorSetMatOrderingType_Factor);
375:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
376:                     PCFactorSetReuseOrdering_ILU);
377:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
378:                     PCFactorSetReuseFill_ILU);
379:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
380:                     PCFactorSetLevels_Factor);
381:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
382:                     PCFactorSetUseInPlace_ILU);
383:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
384:                     PCFactorSetAllowDiagonalFill_Factor);
385:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
386:                     PCFactorSetPivotInBlocks_Factor);
387:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
388:                     PCFactorReorderForNonzeroDiagonal_ILU);
389:   return(0);
390: }