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  PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 58: {
 59:   PC_ILU         *ilu = (PC_ILU*)pc->data;

 63:   if (pc->setupcalled && (!ilu->usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 64:     pc->setupcalled   = 0;
 65:     PCDestroy_ILU_Internal(pc);
 66:   }
 67:   ilu->usedt                      = PETSC_TRUE;
 68:   ((PC_Factor*)ilu)->info.dt      = dt;
 69:   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
 70:   ((PC_Factor*)ilu)->info.dtcount = dtcount;
 71:   ((PC_Factor*)ilu)->info.fill    = PETSC_DEFAULT;
 72:   return(0);
 73: }

 79: PetscErrorCode  PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
 80: {
 81:   PC_ILU *ilu = (PC_ILU*)pc->data;

 84:   ilu->reuseordering = flag;
 85:   return(0);
 86: }

 92: PetscErrorCode  PCFactorSetUseInPlace_ILU(PC pc)
 93: {
 94:   PC_ILU *dir = (PC_ILU*)pc->data;

 97:   dir->inplace = PETSC_TRUE;
 98:   return(0);
 99: }

104: static PetscErrorCode PCSetFromOptions_ILU(PC pc)
105: {
107:   PetscInt       dtmax = 3,itmp;
108:   PetscTruth     flg,set;
109:   PetscReal      dt[3];
110:   char           tname[256];
111:   PC_ILU         *ilu = (PC_ILU*)pc->data;
112:   PetscFList     ordlist;
113:   PetscReal      tol;

116:   if (!MatOrderingRegisterAllCalled) {MatOrderingRegisterAll(PETSC_NULL);}
117:   PetscOptionsHead("ILU Options");
118:     PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);
119:     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
120:     PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);
121:     if (flg) ilu->inplace = PETSC_TRUE;
122:     PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);
123:     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
124:     PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);
125:     if (flg) ilu->reusefill = PETSC_TRUE;
126:     PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);
127:     if (flg) ilu->reuseordering = PETSC_TRUE;
128:     PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);
129:     if (flg) {
130:       PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);
131:     }
132:     PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)ilu)->info.shiftnz,&((PC_Factor*)ilu)->info.shiftnz,0);
133:     flg = (((PC_Factor*)ilu)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE;
134:     PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);
135:     PCFactorSetShiftPd(pc,flg);
136:     PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)ilu)->info.zeropivot,&((PC_Factor*)ilu)->info.zeropivot,0);

138:     dt[0] = ((PC_Factor*)ilu)->info.dt;
139:     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
140:     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
141:     PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);
142:     if (flg) {
143:       PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);
144:     }
145:     PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)ilu)->info.fill,&((PC_Factor*)ilu)->info.fill,&flg);
146:     PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);
147:     if (flg) {
148:       tol = PETSC_DECIDE;
149:       PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);
150:       PCFactorReorderForNonzeroDiagonal(pc,tol);
151:     }

153:     MatGetOrderingList(&ordlist);
154:     PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)ilu)->ordering,tname,256,&flg);
155:     if (flg) {
156:       PCFactorSetMatOrderingType(pc,tname);
157:     }
158:     flg = ((PC_Factor*)ilu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
159:     PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);
160:     if (set) {
161:       PCFactorSetPivotInBlocks(pc,flg);
162:     }
163:     PetscOptionsName("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",&flg);
164:     if (flg) {
165:       PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);
166:     }
167:     PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",((PC_Factor*)ilu)->info.shiftinblocks,&((PC_Factor*)ilu)->info.shiftinblocks,0);

169:   PetscOptionsTail();
170:   return(0);
171: }

175: static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
176: {
177:   PC_ILU         *ilu = (PC_ILU*)pc->data;
179:   PetscTruth     isstring,iascii;

182:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
183:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
184:   if (iascii) {
185:     if (ilu->usedt) {
186:         PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %G\n",((PC_Factor*)ilu)->info.dt);
187:         PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)((PC_Factor*)ilu)->info.dtcount);
188:         PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %G\n",((PC_Factor*)ilu)->info.dtcol);
189:     } else if (((PC_Factor*)ilu)->info.levels == 1) {
190:         PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);
191:     } else {
192:         PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);
193:     }
194:     PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio allocated %G\n",((PC_Factor*)ilu)->info.fill);
195:     PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %G\n",((PC_Factor*)ilu)->info.zeropivot);
196:     if (((PC_Factor*)ilu)->info.shiftpd) {PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");}
197:     if (((PC_Factor*)ilu)->info.shiftnz) {PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift to prevent zero pivot\n");}
198:     if (((PC_Factor*)ilu)->info.shiftinblocks) {PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift on blocks to prevent zero pivot\n");}
199:     if (ilu->inplace) {PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");}
200:     else              {PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");}
201:     PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);
202:     if (ilu->reusefill)     {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");}
203:     if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");}
204:     if (((PC_Factor*)ilu)->fact) {
205:       PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);
206:       PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");
207:       PetscViewerASCIIPushTab(viewer);
208:       PetscViewerASCIIPushTab(viewer);
209:       PetscViewerASCIIPushTab(viewer);
210:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
211:       MatView(((PC_Factor*)ilu)->fact,viewer);
212:       PetscViewerPopFormat(viewer);
213:       PetscViewerASCIIPopTab(viewer);
214:       PetscViewerASCIIPopTab(viewer);
215:       PetscViewerASCIIPopTab(viewer);
216:     }
217:   } else if (isstring) {
218:     PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);
219:   } else {
220:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
221:   }
222:   return(0);
223: }

227: static PetscErrorCode PCSetUp_ILU(PC pc)
228: {
230:   PC_ILU         *ilu = (PC_ILU*)pc->data;
231:   MatInfo        info;

234:   if (ilu->inplace) {
235:     CHKMEMQ;
236:     if (!pc->setupcalled) {

238:       /* In-place factorization only makes sense with the natural ordering,
239:          so we only need to get the ordering once, even if nonzero structure changes */
240:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
241:       if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
242:       if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
243:     }

245:     /* In place ILU only makes sense with fill factor of 1.0 because 
246:        cannot have levels of fill */
247:     ((PC_Factor*)ilu)->info.fill          = 1.0;
248:     ((PC_Factor*)ilu)->info.diagonal_fill = 0;
249:     MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
250:     CHKMEMQ;
251:     ((PC_Factor*)ilu)->fact = pc->pmat;
252:   } else if (ilu->usedt) {
253:     if (!pc->setupcalled) {
254:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
255:     CHKMEMQ;
256:       if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
257:       if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
258:       MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);
259:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
260:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
261:     CHKMEMQ;
262:       MatDestroy(((PC_Factor*)ilu)->fact);
263:     CHKMEMQ;
264:       if (!ilu->reuseordering) {
265:         if (ilu->row) {ISDestroy(ilu->row);}
266:         if (ilu->col) {ISDestroy(ilu->col);}
267:         MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
268:         if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
269:         if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
270:       }
271:       MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);
272:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
273:     } else if (!ilu->reusefill) {
274:       MatDestroy(((PC_Factor*)ilu)->fact);
275:       MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);
276:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
277:     } else {
278:       MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);
279:     }
280:   } else {
281:     if (!pc->setupcalled) {
282:       /* first time in so compute reordering and symbolic factorization */
283:       MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
284:       if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
285:       if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
286:       /*  Remove zeros along diagonal?     */
287:       if (ilu->nonzerosalongdiagonal) {
288:         MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
289:       }
290:     CHKMEMQ;
291:       MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
292:     CHKMEMQ;
293:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
294:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
295:       ilu->actualfill = info.fill_ratio_needed;
296:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
297:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
298:       if (!ilu->reuseordering) {
299:         /* compute a new ordering for the ILU */
300:         ISDestroy(ilu->row);
301:         ISDestroy(ilu->col);
302:         MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);
303:         if (ilu->row) {PetscLogObjectParent(pc,ilu->row);}
304:         if (ilu->col) {PetscLogObjectParent(pc,ilu->col);}
305:         /*  Remove zeros along diagonal?     */
306:         if (ilu->nonzerosalongdiagonal) {
307:           MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);
308:         }
309:       }
310:       MatDestroy(((PC_Factor*)ilu)->fact);
311:       MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);
312:       MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);
313:       MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);
314:       ilu->actualfill = info.fill_ratio_needed;
315:       PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);
316:     }
317:     CHKMEMQ;
318:     MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);
319:     CHKMEMQ;
320:   }
321:   return(0);
322: }

326: static PetscErrorCode PCDestroy_ILU(PC pc)
327: {
328:   PC_ILU         *ilu = (PC_ILU*)pc->data;

332:   PCDestroy_ILU_Internal(pc);
333:   PetscStrfree(((PC_Factor*)ilu)->ordering);
334:   PetscFree(ilu);
335:   return(0);
336: }

340: static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
341: {
342:   PC_ILU         *ilu = (PC_ILU*)pc->data;

346:   MatSolve(((PC_Factor*)ilu)->fact,x,y);
347:   return(0);
348: }

352: static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
353: {
354:   PC_ILU         *ilu = (PC_ILU*)pc->data;

358:   MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);
359:   return(0);
360: }

362: /*MC
363:      PCILU - Incomplete factorization preconditioners.

365:    Options Database Keys:
366: +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
367: .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
368:                       its factorization (overwrites original matrix)
369: .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
370: .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
371: .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
372: .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
373: .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
374:                                    this decreases the chance of getting a zero pivot
375: .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
376: .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
377:                              than 1 the diagonal blocks are factored with partial pivoting (this increases the 
378:                              stability of the ILU factorization
379: .  -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization
380: .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
381: -  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
382:    is optional with true being the default

384:    Level: beginner

386:   Concepts: incomplete factorization

388:    Notes: Only implemented for some matrix formats. (for parallel use you 
389:              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended
390:              unless you really want a parallel ILU).

392:           For BAIJ matrices this implements a point block ILU

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

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

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


407: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
408:            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(),
409:            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
410:            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
411:            PCFactorSetShiftNonzero(),PCFactorSetShiftPd()

413: M*/

418: PetscErrorCode  PCCreate_ILU(PC pc)
419: {
421:   PC_ILU         *ilu;

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

426:   ((PC_Factor*)ilu)->fact                    = 0;
427:   MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);
428:   ((PC_Factor*)ilu)->info.levels             = 0;
429:   ((PC_Factor*)ilu)->info.fill               = 1.0;
430:   ilu->col                     = 0;
431:   ilu->row                     = 0;
432:   ilu->inplace                 = PETSC_FALSE;
433:   PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);
434:   ilu->reuseordering           = PETSC_FALSE;
435:   ilu->usedt                   = PETSC_FALSE;
436:   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
437:   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
438:   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
439:   ((PC_Factor*)ilu)->info.shiftnz            = 1.e-12;
440:   ((PC_Factor*)ilu)->info.shiftpd            = 0.0; /* false */
441:   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
442:   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
443:   ((PC_Factor*)ilu)->info.shiftinblocks      = 1.e-12;
444:   ilu->reusefill               = PETSC_FALSE;
445:   ((PC_Factor*)ilu)->info.diagonal_fill      = 0;
446:   pc->data                     = (void*)ilu;

448:   pc->ops->destroy             = PCDestroy_ILU;
449:   pc->ops->apply               = PCApply_ILU;
450:   pc->ops->applytranspose      = PCApplyTranspose_ILU;
451:   pc->ops->setup               = PCSetUp_ILU;
452:   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
453:   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
454:   pc->ops->view                = PCView_ILU;
455:   pc->ops->applyrichardson     = 0;

457:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
458:                     PCFactorSetZeroPivot_Factor);
459:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
460:                     PCFactorSetShiftNonzero_Factor);
461:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
462:                     PCFactorSetShiftPd_Factor);

464:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
465:                     PCFactorGetMatSolverPackage_Factor);
466:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
467:                     PCFactorSetUseDropTolerance_ILU);
468:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
469:                     PCFactorSetFill_Factor);
470:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
471:                     PCFactorSetMatOrderingType_Factor);
472:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
473:                     PCFactorSetReuseOrdering_ILU);
474:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
475:                     PCFactorSetReuseFill_ILU);
476:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
477:                     PCFactorSetLevels_Factor);
478:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
479:                     PCFactorSetUseInPlace_ILU);
480:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
481:                     PCFactorSetAllowDiagonalFill_Factor);
482:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
483:                     PCFactorSetPivotInBlocks_Factor);
484:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
485:                     PCFactorSetShiftInBlocks_Factor);
486:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
487:                     PCFactorReorderForNonzeroDiagonal_ILU);
488:   return(0);
489: }