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       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:     /*
124:     PetscInt       dtmax = 3;
125:     PetscOptionsRealArray("-pc_factor_drop_tolerance,","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);
126:     if (flg) {
127:       PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);
128:     }
129:     */
130:     PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
131:     if (flg) {
132:       tol = PETSC_DECIDE;
133:       PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);
134:       PCFactorReorderForNonzeroDiagonal(pc,tol);
135:     }

137:   PetscOptionsTail();
138:   return(0);
139: }

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

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

167: static PetscErrorCode PCSetUp_ILU(PC pc)
168: {
170:   PC_ILU         *ilu = (PC_ILU*)pc->data;
171:   MatInfo        info;

174:   if (ilu->inplace) {
175:     CHKMEMQ;
176:     if (!pc->setupcalled) {

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

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

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

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

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

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

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

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

274: /*MC
275:      PCILU - Incomplete factorization preconditioners.

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

291:    Level: beginner

293:   Concepts: incomplete factorization

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

297:           For BAIJ matrices this implements a point block ILU

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

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

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


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

317: M*/

322: PetscErrorCode  PCCreate_ILU(PC pc)
323: {
325:   PC_ILU         *ilu;

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

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

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

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