Actual source code: ilu.c

  1: /*$Id: ilu.c,v 1.172 2001/08/07 21:30:27 bsmith Exp $*/
  2: /*
  3:    Defines a ILU factorization preconditioner for any Mat implementation
  4: */
 5:  #include src/sles/pc/pcimpl.h
 6:  #include src/sles/pc/impls/ilu/ilu.h
 7:  #include src/mat/matimpl.h

  9: /* ------------------------------------------------------------------------------------------*/
 10: EXTERN_C_BEGIN
 11: int PCILUSetSinglePrecisionSolves_ILU(PC pc,PetscTruth flag)
 12: {
 13:   PC_ILU *dir;

 16:   dir = (PC_ILU*)pc->data;
 17:   dir->single_precision_solve = flag;
 18:   return(0);
 19: }
 20: EXTERN_C_END

 22: EXTERN_C_BEGIN
 23: int PCILUSetDamping_ILU(PC pc,PetscReal damping)
 24: {
 25:   PC_ILU *dir;

 28:   dir = (PC_ILU*)pc->data;
 29:   dir->info.damping = damping;
 30:   dir->info.damp    = 1.0;
 31:   return(0);
 32: }
 33: EXTERN_C_END

 35: EXTERN_C_BEGIN
 36: int PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,int dtcount)
 37: {
 38:   PC_ILU *ilu;

 41:   ilu = (PC_ILU*)pc->data;
 42:   ilu->usedt         = PETSC_TRUE;
 43:   ilu->info.dt       = dt;
 44:   ilu->info.dtcol    = dtcol;
 45:   ilu->info.dtcount  = dtcount;
 46:   ilu->info.fill     = PETSC_DEFAULT;
 47:   return(0);
 48: }
 49: EXTERN_C_END

 51: EXTERN_C_BEGIN
 52: int PCILUSetFill_ILU(PC pc,PetscReal fill)
 53: {
 54:   PC_ILU *dir;

 57:   dir            = (PC_ILU*)pc->data;
 58:   dir->info.fill = fill;
 59:   return(0);
 60: }
 61: EXTERN_C_END

 63: EXTERN_C_BEGIN
 64: int PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
 65: {
 66:   PC_ILU *dir = (PC_ILU*)pc->data;
 67:   int    ierr;
 68: 
 70:   PetscStrfree(dir->ordering);
 71:   PetscStrallocpy(ordering,&dir->ordering);
 72:   return(0);
 73: }
 74: EXTERN_C_END

 76: EXTERN_C_BEGIN
 77: int PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
 78: {
 79:   PC_ILU *ilu;

 82:   ilu                = (PC_ILU*)pc->data;
 83:   ilu->reuseordering = flag;
 84:   return(0);
 85: }
 86: EXTERN_C_END

 88: EXTERN_C_BEGIN
 89: int PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
 90: {
 91:   PC_ILU *ilu;

 94:   ilu = (PC_ILU*)pc->data;
 95:   ilu->reusefill = flag;
 96:   if (flag) SETERRQ(1,"Not yet supported");
 97:   return(0);
 98: }
 99: EXTERN_C_END

101: EXTERN_C_BEGIN
102: int PCILUSetLevels_ILU(PC pc,int levels)
103: {
104:   PC_ILU *ilu;

107:   ilu = (PC_ILU*)pc->data;
108:   ilu->info.levels = levels;
109:   return(0);
110: }
111: EXTERN_C_END

113: EXTERN_C_BEGIN
114: int PCILUSetUseInPlace_ILU(PC pc)
115: {
116:   PC_ILU *dir;

119:   dir          = (PC_ILU*)pc->data;
120:   dir->inplace = PETSC_TRUE;
121:   return(0);
122: }
123: EXTERN_C_END

125: EXTERN_C_BEGIN
126: int PCILUSetAllowDiagonalFill_ILU(PC pc)
127: {
128:   PC_ILU *dir;

131:   dir                 = (PC_ILU*)pc->data;
132:   dir->info.diagonal_fill = 1;
133:   return(0);
134: }
135: EXTERN_C_END

137: /* ------------------------------------------------------------------------------------------*/
138: /*@
139:    PCILUSetSinglePrecisionSolves - Sets precision of triangular solves to single precision

141:    Collective on PC
142:    
143:    Input Parameters:
144: +  pc - the preconditioner context
145: -  flag - PETSC_TRUE to use single precision, default is PETSC_FALSE

147:    Options Database Key:
148: .  -pc_ilu_single_precision_solves - activates PCILUSetSinglePrecisionSolves
149: .  -pc_ilu_single_precision_solves <flag> - flag of TRUE/FALSE will override this call in the code.

151:    Level: intermediate

153: .keywords: PC, set, solve, precision
154: @*/
155: int PCILUSetSinglePrecisionSolves(PC pc,PetscTruth flag)
156: {
157:   int ierr,(*f)(PC,PetscTruth);

161:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetSinglePrecisionSolves_C",(void (**)(void))&f);
162:   if (f) {
163:     (*f)(pc,flag);
164:   }
165:   return(0);
166: }

168: /*@
169:    PCILUSetDamping - adds this quantity to the diagonal of the matrix during the 
170:      ILU numerical factorization

172:    Collective on PC
173:    
174:    Input Parameters:
175: +  pc - the preconditioner context
176: -  damping - amount of damping

178:    Options Database Key:
179: .  -pc_ilu_damping <damping> - Sets damping amount

181:    Level: intermediate

183: .keywords: PC, set, factorization, direct, fill

185: .seealso: PCILUSetFill(), PCLUSetDamp()
186: @*/
187: int PCILUSetDamping(PC pc,PetscReal damping)
188: {
189:   int ierr,(*f)(PC,PetscReal);

193:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);
194:   if (f) {
195:     (*f)(pc,damping);
196:   }
197:   return(0);
198: }

200: /*@
201:    PCILUSetUseDropTolerance - The preconditioner will use an ILU 
202:    based on a drop tolerance.

204:    Collective on PC

206:    Input Parameters:
207: +  pc - the preconditioner context
208: .  dt - the drop tolerance, try from 1.e-10 to .1
209: .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
210: -  maxrowcount - the max number of nonzeros allowed in a row, best value
211:                  depends on the number of nonzeros in row of original matrix

213:    Options Database Key:
214: .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance

216:    Level: intermediate

218:     Notes:
219:       This uses the iludt() code of Saad's SPARSKIT package

221: .keywords: PC, levels, reordering, factorization, incomplete, ILU
222: @*/
223: int PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,int maxrowcount)
224: {
225:   int ierr,(*f)(PC,PetscReal,PetscReal,int);

229:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);
230:   if (f) {
231:     (*f)(pc,dt,dtcol,maxrowcount);
232:   }
233:   return(0);
234: }

236: /*@
237:    PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
238:    where fill = number nonzeros in factor/number nonzeros in original matrix.

240:    Collective on PC

242:    Input Parameters:
243: +  pc - the preconditioner context
244: -  fill - amount of expected fill

246:    Options Database Key:
247: $  -pc_ilu_fill <fill>

249:    Note:
250:    For sparse matrix factorizations it is difficult to predict how much 
251:    fill to expect. By running with the option -log_info PETSc will print the 
252:    actual amount of fill used; allowing you to set the value accurately for
253:    future runs. But default PETSc uses a value of 1.0

255:    Level: intermediate

257: .keywords: PC, set, factorization, direct, fill

259: .seealso: PCLUSetFill()
260: @*/
261: int PCILUSetFill(PC pc,PetscReal fill)
262: {
263:   int ierr,(*f)(PC,PetscReal);

267:   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
268:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);
269:   if (f) {
270:     (*f)(pc,fill);
271:   }
272:   return(0);
273: }

275: /*@C
276:     PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to 
277:     be used in the ILU factorization.

279:     Collective on PC

281:     Input Parameters:
282: +   pc - the preconditioner context
283: -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM

285:     Options Database Key:
286: .   -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine

288:     Level: intermediate

290:     Notes: natural ordering is used by default

292: .seealso: PCLUSetMatOrdering()

294: .keywords: PC, ILU, set, matrix, reordering

296: @*/
297: int PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
298: {
299:   int ierr,(*f)(PC,MatOrderingType);

303:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);
304:   if (f) {
305:     (*f)(pc,ordering);
306:   }
307:   return(0);
308: }

310: /*@
311:    PCILUSetReuseOrdering - When similar matrices are factored, this
312:    causes the ordering computed in the first factor to be used for all
313:    following factors; applies to both fill and drop tolerance ILUs.

315:    Collective on PC

317:    Input Parameters:
318: +  pc - the preconditioner context
319: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

321:    Options Database Key:
322: .  -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()

324:    Level: intermediate

326: .keywords: PC, levels, reordering, factorization, incomplete, ILU

328: .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
329: @*/
330: int PCILUSetReuseOrdering(PC pc,PetscTruth flag)
331: {
332:   int ierr,(*f)(PC,PetscTruth);

336:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);
337:   if (f) {
338:     (*f)(pc,flag);
339:   }
340:   return(0);
341: }

343: /*@
344:    PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
345:    this causes later ones to use the fill computed in the initial factorization.

347:    Collective on PC

349:    Input Parameters:
350: +  pc - the preconditioner context
351: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

353:    Options Database Key:
354: .  -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()

356:    Level: intermediate

358: .keywords: PC, levels, reordering, factorization, incomplete, ILU

360: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
361: @*/
362: int PCILUDTSetReuseFill(PC pc,PetscTruth flag)
363: {
364:   int ierr,(*f)(PC,PetscTruth);

368:   PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);
369:   if (f) {
370:     (*f)(pc,flag);
371:   }
372:   return(0);
373: }

375: /*@
376:    PCILUSetLevels - Sets the number of levels of fill to use.

378:    Collective on PC

380:    Input Parameters:
381: +  pc - the preconditioner context
382: -  levels - number of levels of fill

384:    Options Database Key:
385: .  -pc_ilu_levels <levels> - Sets fill level

387:    Level: intermediate

389: .keywords: PC, levels, fill, factorization, incomplete, ILU
390: @*/
391: int PCILUSetLevels(PC pc,int levels)
392: {
393:   int ierr,(*f)(PC,int);

397:   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
398:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);
399:   if (f) {
400:     (*f)(pc,levels);
401:   }
402:   return(0);
403: }

405: /*@
406:    PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be 
407:    treated as level 0 fill even if there is no non-zero location.

409:    Collective on PC

411:    Input Parameters:
412: +  pc - the preconditioner context

414:    Options Database Key:
415: .  -pc_ilu_diagonal_fill

417:    Notes:
418:    Does not apply with 0 fill.

420:    Level: intermediate

422: .keywords: PC, levels, fill, factorization, incomplete, ILU
423: @*/
424: int PCILUSetAllowDiagonalFill(PC pc)
425: {
426:   int ierr,(*f)(PC);

430:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);
431:   if (f) {
432:     (*f)(pc);
433:   }
434:   return(0);
435: }

437: /*@
438:    PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
439:    Collective on PC

441:    Input Parameters:
442: .  pc - the preconditioner context

444:    Options Database Key:
445: .  -pc_ilu_in_place - Activates in-place factorization

447:    Notes:
448:    PCILUSetUseInPlace() is intended for use with matrix-free variants of
449:    Krylov methods, or when a different matrices are employed for the linear
450:    system and preconditioner, or with ASM preconditioning.  Do NOT use 
451:    this option if the linear system
452:    matrix also serves as the preconditioning matrix, since the factored
453:    matrix would then overwrite the original matrix. 

455:    Only works well with ILU(0).

457:    Level: intermediate

459: .keywords: PC, set, factorization, inplace, in-place, ILU

461: .seealso:  PCLUSetUseInPlace()
462: @*/
463: int PCILUSetUseInPlace(PC pc)
464: {
465:   int ierr,(*f)(PC);

469:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);
470:   if (f) {
471:     (*f)(pc);
472:   }
473:   return(0);
474: }

476: /*@
477:     PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
478:       with BAIJ or SBAIJ matrices

480:     Collective on PC

482:     Input Parameters:
483: +   pc - the preconditioner context
484: -   pivot - PETSC_TRUE or PETSC_FALSE

486:     Options Database Key:
487: .   -pc_ilu_pivot_in_blocks <true,false>

489:     Level: intermediate

491: .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
492: @*/
493: int PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
494: {
495:   int ierr,(*f)(PC,PetscTruth);

498:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);
499:   if (f) {
500:     (*f)(pc,pivot);
501:   }
502:   return(0);
503: }

505: /* ------------------------------------------------------------------------------------------*/

507: EXTERN_C_BEGIN
508: int PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
509: {
510:   PC_ILU *dir = (PC_ILU*)pc->data;

513:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
514:   return(0);
515: }
516: EXTERN_C_END

518: static int PCSetFromOptions_ILU(PC pc)
519: {
520:   int        ierr,dtmax = 3,itmp;
521:   PetscTruth flg,single_prec,set;
522:   PetscReal  dt[3];
523:   char       tname[256];
524:   PC_ILU     *ilu = (PC_ILU*)pc->data;
525:   PetscFList ordlist;
526:   PetscReal  tol;

529:   MatOrderingRegisterAll(PETSC_NULL);
530:   PetscOptionsHead("ILU Options");
531:     PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(int)ilu->info.levels,&itmp,&flg);
532:     if (flg) ilu->info.levels = (double) itmp;
533:     PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);
534:     PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);
535:     ilu->info.diagonal_fill = (double) flg;
536:     PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);
537:     PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);
538:     PetscOptionsLogical("-pc_ilu_single_precision_solves","Precision for triangular solves","PCILUSetSinglePrecisionSolves",ilu->single_precision_solve,&single_prec,&flg);
539:     if (flg) {
540:       PCILUSetSinglePrecisionSolves(pc,single_prec);
541:     }
542:     PetscOptionsHasName(pc->prefix,"-pc_ilu_damping",&flg);
543:     if (flg) {
544:       PCILUSetDamping(pc,0.0);
545:     }
546:     PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);
547:     PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);

549:     dt[0] = ilu->info.dt;
550:     dt[1] = ilu->info.dtcol;
551:     dt[2] = ilu->info.dtcount;
552:     PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);
553:     if (flg) {
554:       PCILUSetUseDropTolerance(pc,dt[0],dt[1],(int)dt[2]);
555:     }
556:     PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);
557:     PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);

559:     MatGetOrderingList(&ordlist);
560:     PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);
561:     if (flg) {
562:       PCILUSetMatOrdering(pc,tname);
563:     }
564:     flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
565:     PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);
566:     if (set) {
567:       PCILUSetPivotInBlocks(pc,flg);
568:     }
569:   PetscOptionsTail();
570:   return(0);
571: }

573: static int PCView_ILU(PC pc,PetscViewer viewer)
574: {
575:   PC_ILU     *ilu = (PC_ILU*)pc->data;
576:   int        ierr;
577:   PetscTruth isstring,isascii;

580:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
581:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
582:   if (isascii) {
583:     if (ilu->usedt) {
584:         PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %gn",ilu->info.dt);
585:         PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %dn",(int)ilu->info.dtcount);
586:         PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %gn",ilu->info.dtcol);
587:     } else if (ilu->info.levels == 1) {
588:         PetscViewerASCIIPrintf(viewer,"  ILU: %d level of filln",(int)ilu->info.levels);
589:     } else {
590:         PetscViewerASCIIPrintf(viewer,"  ILU: %d levels of filln",(int)ilu->info.levels);
591:     }
592:     PetscViewerASCIIPrintf(viewer,"  ILU: max fill ratio allocated %gn",ilu->info.fill);
593:     PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %gn",ilu->info.zeropivot);
594:     if (ilu->inplace) {PetscViewerASCIIPrintf(viewer,"       in-place factorizationn");}
595:     else              {PetscViewerASCIIPrintf(viewer,"       out-of-place factorizationn");}
596:     PetscViewerASCIIPrintf(viewer,"       matrix ordering: %sn",ilu->ordering);
597:     if (ilu->reusefill)     {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorizationn");}
598:     if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorizationn");}
599:     if (ilu->fact) {
600:       PetscViewerASCIIPrintf(viewer,"       Factored matrix followsn");
601:       PetscViewerASCIIPushTab(viewer);
602:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
603:       MatView(ilu->fact,viewer);
604:       PetscViewerPopFormat(viewer);
605:       PetscViewerASCIIPopTab(viewer);
606:     }
607:   } else if (isstring) {
608:     PetscViewerStringSPrintf(viewer," lvls=%g,order=%s",ilu->info.levels,ilu->ordering);
609:   } else {
610:     SETERRQ1(1,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
611:   }
612:   return(0);
613: }

615: static int PCSetUp_ILU(PC pc)
616: {
617:   int        ierr;
618:   PetscTruth flg;
619:   PC_ILU     *ilu = (PC_ILU*)pc->data;

622:   if (ilu->inplace) {
623:     if (!pc->setupcalled) {

625:       /* In-place factorization only makes sense with the natural ordering,
626:          so we only need to get the ordering once, even if nonzero structure changes */
627:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
628:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
629:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
630:     }

632:     /* In place ILU only makes sense with fill factor of 1.0 because 
633:        cannot have levels of fill */
634:     ilu->info.fill          = 1.0;
635:     ilu->info.diagonal_fill = 0;
636:     MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);
637:     ilu->fact = pc->pmat;
638:   } else if (ilu->usedt) {
639:     if (!pc->setupcalled) {
640:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
641:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
642:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
643:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
644:       PetscLogObjectParent(pc,ilu->fact);
645:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
646:       MatDestroy(ilu->fact);
647:       if (!ilu->reuseordering) {
648:         if (ilu->row) {ISDestroy(ilu->row);}
649:         if (ilu->col) {ISDestroy(ilu->col);}
650:         MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
651:         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
652:         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
653:       }
654:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
655:       PetscLogObjectParent(pc,ilu->fact);
656:     } else if (!ilu->reusefill) {
657:       MatDestroy(ilu->fact);
658:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
659:       PetscLogObjectParent(pc,ilu->fact);
660:     } else {
661:       MatLUFactorNumeric(pc->pmat,&ilu->fact);
662:     }
663:    } else {
664:     if (!pc->setupcalled) {
665:       /* first time in so compute reordering and symbolic factorization */
666:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
667:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
668:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
669:       /*  Remove zeros along diagonal?     */
670:       PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
671:       if (flg) {
672:         PetscReal ntol = 1.e-10;
673:         PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
674:         MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
675:       }

677:       MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
678:       PetscLogObjectParent(pc,ilu->fact);
679:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
680:       if (!ilu->reuseordering) {
681:         /* compute a new ordering for the ILU */
682:         ISDestroy(ilu->row);
683:         ISDestroy(ilu->col);
684:         MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
685:         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
686:         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
687:         /*  Remove zeros along diagonal?     */
688:         PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
689:         if (flg) {
690:           PetscReal ntol = 1.e-10;
691:           PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
692:           MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
693:         }
694:       }
695:       MatDestroy(ilu->fact);
696:       MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
697:       PetscLogObjectParent(pc,ilu->fact);
698:     }
699:     MatLUFactorNumeric(pc->pmat,&ilu->fact);
700:   }
701:   if (ilu->single_precision_solve) {
702:     MatSetOption(ilu->fact,MAT_USE_SINGLE_PRECISION_SOLVES);
703:   }
704:   return(0);
705: }

707: static int PCDestroy_ILU(PC pc)
708: {
709:   PC_ILU *ilu = (PC_ILU*)pc->data;
710:   int    ierr;

713:   if (!ilu->inplace && ilu->fact) {MatDestroy(ilu->fact);}
714:   if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(ilu->row);}
715:   if (ilu->col) {ISDestroy(ilu->col);}
716:   PetscStrfree(ilu->ordering);
717:   PetscFree(ilu);
718:   return(0);
719: }

721: static int PCApply_ILU(PC pc,Vec x,Vec y)
722: {
723:   PC_ILU *ilu = (PC_ILU*)pc->data;
724:   int    ierr;

727:   MatSolve(ilu->fact,x,y);
728:   return(0);
729: }

731: static int PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
732: {
733:   PC_ILU *ilu = (PC_ILU*)pc->data;
734:   int    ierr;

737:   MatSolveTranspose(ilu->fact,x,y);
738:   return(0);
739: }

741: static int PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
742: {
743:   PC_ILU *ilu = (PC_ILU*)pc->data;

746:   if (!ilu->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
747:   *mat = ilu->fact;
748:   return(0);
749: }

751: EXTERN_C_BEGIN
752: int PCCreate_ILU(PC pc)
753: {
754:   int    ierr;
755:   PC_ILU *ilu;

758:   PetscNew(PC_ILU,&ilu);
759:   PetscLogObjectMemory(pc,sizeof(PC_ILU));

761:   ilu->fact                    = 0;
762:   ilu->info.levels             = 0;
763:   ilu->info.fill               = 1.0;
764:   ilu->col                     = 0;
765:   ilu->row                     = 0;
766:   ilu->inplace                 = PETSC_FALSE;
767:   PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);
768:   ilu->reuseordering           = PETSC_FALSE;
769:   ilu->usedt                   = PETSC_FALSE;
770:   ilu->info.dt                 = PETSC_DEFAULT;
771:   ilu->info.dtcount            = PETSC_DEFAULT;
772:   ilu->info.dtcol              = PETSC_DEFAULT;
773:   ilu->info.damp               = 0.0;
774:   ilu->info.damping            = 0.0;
775:   ilu->info.zeropivot          = 1.e-12;
776:   ilu->info.pivotinblocks      = 1.0;
777:   ilu->reusefill               = PETSC_FALSE;
778:   ilu->info.diagonal_fill      = 0;
779:   ilu->single_precision_solve  = PETSC_FALSE;
780:   pc->data                     = (void*)ilu;

782:   pc->ops->destroy             = PCDestroy_ILU;
783:   pc->ops->apply               = PCApply_ILU;
784:   pc->ops->applytranspose      = PCApplyTranspose_ILU;
785:   pc->ops->setup               = PCSetUp_ILU;
786:   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
787:   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ILU;
788:   pc->ops->view                = PCView_ILU;
789:   pc->ops->applyrichardson     = 0;

791:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetSinglePrecisionSolves_C",
792:                     "PCILUSetSinglePrecisionSolves_ILU",
793:                      PCILUSetSinglePrecisionSolves_ILU);
794:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
795:                     PCILUSetUseDropTolerance_ILU);
796:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
797:                     PCILUSetFill_ILU);
798:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
799:                     PCILUSetDamping_ILU);
800:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
801:                     PCILUSetMatOrdering_ILU);
802:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
803:                     PCILUSetReuseOrdering_ILU);
804:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
805:                     PCILUDTSetReuseFill_ILUDT);
806:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
807:                     PCILUSetLevels_ILU);
808:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
809:                     PCILUSetUseInPlace_ILU);
810:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
811:                     PCILUSetAllowDiagonalFill_ILU);
812:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
813:                     PCILUSetPivotInBlocks_ILU);
814:   return(0);
815: }
816: EXTERN_C_END