Actual source code: ilu.c

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

  8: /* ------------------------------------------------------------------------------------------*/

 13: PetscErrorCode PCILUSetZeroPivot_ILU(PC pc,PetscReal z)
 14: {
 15:   PC_ILU *lu;

 18:   lu                 = (PC_ILU*)pc->data;
 19:   lu->info.zeropivot = z;
 20:   return(0);
 21: }

 26: PetscErrorCode PCDestroy_ILU_Internal(PC pc)
 27: {
 28:   PC_ILU         *ilu = (PC_ILU*)pc->data;

 32:   if (!ilu->inplace && ilu->fact) {MatDestroy(ilu->fact);}
 33:   if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(ilu->row);}
 34:   if (ilu->col) {ISDestroy(ilu->col);}
 35:   return(0);
 36: }

 41: PetscErrorCode PCILUSetDamping_ILU(PC pc,PetscReal damping)
 42: {
 43:   PC_ILU *dir;

 46:   dir = (PC_ILU*)pc->data;
 47:   if (damping == (PetscReal) PETSC_DECIDE) {
 48:     dir->info.damping = 1.e-12;
 49:   } else {
 50:     dir->info.damping = damping;
 51:   }
 52:   return(0);
 53: }

 59: PetscErrorCode PCILUSetShift_ILU(PC pc,PetscTruth shift)
 60: {
 61:   PC_ILU *dir;

 64:   dir = (PC_ILU*)pc->data;
 65:   dir->info.shift = shift;
 66:   if (shift) dir->info.shift_fraction = 0.0;
 67:   return(0);
 68: }

 74: PetscErrorCode PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 75: {
 76:   PC_ILU         *ilu;

 80:   ilu = (PC_ILU*)pc->data;
 81:   if (pc->setupcalled && (ilu->usedt == PETSC_FALSE || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) {
 82:     pc->setupcalled   = 0;
 83:     PCDestroy_ILU_Internal(pc);
 84:   }
 85:   ilu->usedt        = PETSC_TRUE;
 86:   ilu->info.dt      = dt;
 87:   ilu->info.dtcol   = dtcol;
 88:   ilu->info.dtcount = dtcount;
 89:   ilu->info.fill    = PETSC_DEFAULT;
 90:   return(0);
 91: }

 97: PetscErrorCode PCILUSetFill_ILU(PC pc,PetscReal fill)
 98: {
 99:   PC_ILU *dir;

102:   dir            = (PC_ILU*)pc->data;
103:   dir->info.fill = fill;
104:   return(0);
105: }

111: PetscErrorCode PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
112: {
113:   PC_ILU         *dir = (PC_ILU*)pc->data;
115:   PetscTruth     flg;
116: 
118:   if (!pc->setupcalled) {
119:      PetscStrfree(dir->ordering);
120:      PetscStrallocpy(ordering,&dir->ordering);
121:   } else {
122:     PetscStrcmp(dir->ordering,ordering,&flg);
123:     if (!flg) {
124:       pc->setupcalled = 0;
125:       PetscStrfree(dir->ordering);
126:       PetscStrallocpy(ordering,&dir->ordering);
127:       /* free the data structures, then create them again */
128:       PCDestroy_ILU_Internal(pc);
129:     }
130:   }
131:   return(0);
132: }

138: PetscErrorCode PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
139: {
140:   PC_ILU *ilu;

143:   ilu                = (PC_ILU*)pc->data;
144:   ilu->reuseordering = flag;
145:   return(0);
146: }

152: PetscErrorCode PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
153: {
154:   PC_ILU *ilu;

157:   ilu = (PC_ILU*)pc->data;
158:   ilu->reusefill = flag;
159:   if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported");
160:   return(0);
161: }

167: PetscErrorCode PCILUSetLevels_ILU(PC pc,PetscInt levels)
168: {
169:   PC_ILU         *ilu;

173:   ilu = (PC_ILU*)pc->data;

175:   if (!pc->setupcalled) {
176:     ilu->info.levels = levels;
177:   } else if (ilu->usedt == PETSC_TRUE || ilu->info.levels != levels) {
178:     ilu->info.levels = levels;
179:     pc->setupcalled  = 0;
180:     ilu->usedt       = PETSC_FALSE;
181:     PCDestroy_ILU_Internal(pc);
182:   }
183:   return(0);
184: }

190: PetscErrorCode PCILUSetUseInPlace_ILU(PC pc)
191: {
192:   PC_ILU *dir;

195:   dir          = (PC_ILU*)pc->data;
196:   dir->inplace = PETSC_TRUE;
197:   return(0);
198: }

204: PetscErrorCode PCILUSetAllowDiagonalFill_ILU(PC pc)
205: {
206:   PC_ILU *dir;

209:   dir                 = (PC_ILU*)pc->data;
210:   dir->info.diagonal_fill = 1;
211:   return(0);
212: }

217: /*@
218:    PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero

220:    Collective on PC
221:    
222:    Input Parameters:
223: +  pc - the preconditioner context
224: -  zero - all pivots smaller than this will be considered zero

226:    Options Database Key:
227: .  -pc_ilu_zeropivot <zero> - Sets the zero pivot size

229:    Level: intermediate

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

233: .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot()
234: @*/
235: PetscErrorCode PCILUSetZeroPivot(PC pc,PetscReal zero)
236: {
237:   PetscErrorCode ierr,(*f)(PC,PetscReal);

241:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);
242:   if (f) {
243:     (*f)(pc,zero);
244:   }
245:   return(0);
246: }

250: /*@
251:    PCILUSetDamping - adds this quantity to the diagonal of the matrix during the 
252:      ILU numerical factorization

254:    Collective on PC
255:    
256:    Input Parameters:
257: +  pc - the preconditioner context
258: -  damping - amount of damping

260:    Options Database Key:
261: .  -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default

263:    Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
264:          pivot, then the damping is doubled until this is alleviated.

266:    Level: intermediate

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

270: .seealso: PCILUSetFill(), PCLUSetDamping(), PCILUSetShift()
271: @*/
272: PetscErrorCode PCILUSetDamping(PC pc,PetscReal damping)
273: {
274:   PetscErrorCode ierr,(*f)(PC,PetscReal);

278:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);
279:   if (f) {
280:     (*f)(pc,damping);
281:   }
282:   return(0);
283: }

287: /*@
288:    PCILUSetShift - specify whether to use Manteuffel shifting of ILU.
289:    If an ILU factorisation breaks down because of nonpositive pivots,
290:    adding sufficient identity to the diagonal will remedy this.
291:    Setting this causes a bisection method to find the minimum shift that
292:    will lead to a well-defined ILU.

294:    Input parameters:
295: +  pc - the preconditioner context
296: -  shifting - PETSC_TRUE to set shift else PETSC_FALSE

298:    Options Database Key:
299: .  -pc_ilu_shift [1/0] - Activate/Deactivate PCILUSetShift(); the value
300:    is optional with 1 being the default

302:    Level: intermediate

304: .keywords: PC, indefinite, factorization, incomplete, ILU

306: .seealso: PCILUSetDamping()
307: @*/
308: PetscErrorCode PCILUSetShift(PC pc,PetscTruth shifting)
309: {
310:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

314:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);
315:   if (f) {
316:     (*f)(pc,shifting);
317:   }
318:   return(0);
319: }


324: /*@
325:    PCILUSetUseDropTolerance - The preconditioner will use an ILU 
326:    based on a drop tolerance.

328:    Collective on PC

330:    Input Parameters:
331: +  pc - the preconditioner context
332: .  dt - the drop tolerance, try from 1.e-10 to .1
333: .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
334: -  maxrowcount - the max number of nonzeros allowed in a row, best value
335:                  depends on the number of nonzeros in row of original matrix

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

340:    Level: intermediate

342:     Notes:
343:       This uses the iludt() code of Saad's SPARSKIT package

345: .keywords: PC, levels, reordering, factorization, incomplete, ILU
346: @*/
347: PetscErrorCode PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
348: {
349:   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);

353:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);
354:   if (f) {
355:     (*f)(pc,dt,dtcol,maxrowcount);
356:   }
357:   return(0);
358: }

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

366:    Collective on PC

368:    Input Parameters:
369: +  pc - the preconditioner context
370: -  fill - amount of expected fill

372:    Options Database Key:
373: $  -pc_ilu_fill <fill>

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

381:    Level: intermediate

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

385: .seealso: PCLUSetFill()
386: @*/
387: PetscErrorCode PCILUSetFill(PC pc,PetscReal fill)
388: {
389:   PetscErrorCode ierr,(*f)(PC,PetscReal);

393:   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
394:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);
395:   if (f) {
396:     (*f)(pc,fill);
397:   }
398:   return(0);
399: }

403: /*@C
404:     PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to 
405:     be used in the ILU factorization.

407:     Collective on PC

409:     Input Parameters:
410: +   pc - the preconditioner context
411: -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM

413:     Options Database Key:
414: .   -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine

416:     Level: intermediate

418:     Notes: natural ordering is used by default

420: .seealso: PCLUSetMatOrdering()

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

424: @*/
425: PetscErrorCode PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
426: {
427:   PetscErrorCode ierr,(*f)(PC,MatOrderingType);

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

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

445:    Collective on PC

447:    Input Parameters:
448: +  pc - the preconditioner context
449: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

451:    Options Database Key:
452: .  -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()

454:    Level: intermediate

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

458: .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
459: @*/
460: PetscErrorCode PCILUSetReuseOrdering(PC pc,PetscTruth flag)
461: {
462:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

466:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);
467:   if (f) {
468:     (*f)(pc,flag);
469:   }
470:   return(0);
471: }

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

479:    Collective on PC

481:    Input Parameters:
482: +  pc - the preconditioner context
483: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

485:    Options Database Key:
486: .  -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()

488:    Level: intermediate

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

492: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
493: @*/
494: PetscErrorCode PCILUDTSetReuseFill(PC pc,PetscTruth flag)
495: {
496:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

500:   PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);
501:   if (f) {
502:     (*f)(pc,flag);
503:   }
504:   return(0);
505: }

509: /*@
510:    PCILUSetLevels - Sets the number of levels of fill to use.

512:    Collective on PC

514:    Input Parameters:
515: +  pc - the preconditioner context
516: -  levels - number of levels of fill

518:    Options Database Key:
519: .  -pc_ilu_levels <levels> - Sets fill level

521:    Level: intermediate

523: .keywords: PC, levels, fill, factorization, incomplete, ILU
524: @*/
525: PetscErrorCode PCILUSetLevels(PC pc,PetscInt levels)
526: {
527:   PetscErrorCode ierr,(*f)(PC,PetscInt);

531:   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
532:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);
533:   if (f) {
534:     (*f)(pc,levels);
535:   }
536:   return(0);
537: }

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

545:    Collective on PC

547:    Input Parameters:
548: +  pc - the preconditioner context

550:    Options Database Key:
551: .  -pc_ilu_diagonal_fill

553:    Notes:
554:    Does not apply with 0 fill.

556:    Level: intermediate

558: .keywords: PC, levels, fill, factorization, incomplete, ILU
559: @*/
560: PetscErrorCode PCILUSetAllowDiagonalFill(PC pc)
561: {
562:   PetscErrorCode ierr,(*f)(PC);

566:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);
567:   if (f) {
568:     (*f)(pc);
569:   }
570:   return(0);
571: }

575: /*@
576:    PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
577:    Collective on PC

579:    Input Parameters:
580: .  pc - the preconditioner context

582:    Options Database Key:
583: .  -pc_ilu_in_place - Activates in-place factorization

585:    Notes:
586:    PCILUSetUseInPlace() is intended for use with matrix-free variants of
587:    Krylov methods, or when a different matrices are employed for the linear
588:    system and preconditioner, or with ASM preconditioning.  Do NOT use 
589:    this option if the linear system
590:    matrix also serves as the preconditioning matrix, since the factored
591:    matrix would then overwrite the original matrix. 

593:    Only works well with ILU(0).

595:    Level: intermediate

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

599: .seealso:  PCLUSetUseInPlace()
600: @*/
601: PetscErrorCode PCILUSetUseInPlace(PC pc)
602: {
603:   PetscErrorCode ierr,(*f)(PC);

607:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);
608:   if (f) {
609:     (*f)(pc);
610:   }
611:   return(0);
612: }

616: /*@
617:     PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
618:       with BAIJ or SBAIJ matrices

620:     Collective on PC

622:     Input Parameters:
623: +   pc - the preconditioner context
624: -   pivot - PETSC_TRUE or PETSC_FALSE

626:     Options Database Key:
627: .   -pc_ilu_pivot_in_blocks <true,false>

629:     Level: intermediate

631: .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
632: @*/
633: PetscErrorCode PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
634: {
635:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

638:   PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);
639:   if (f) {
640:     (*f)(pc,pivot);
641:   }
642:   return(0);
643: }

645: /* ------------------------------------------------------------------------------------------*/

650: PetscErrorCode PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
651: {
652:   PC_ILU *dir = (PC_ILU*)pc->data;

655:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
656:   return(0);
657: }

662: static PetscErrorCode PCSetFromOptions_ILU(PC pc)
663: {
665:   PetscInt       dtmax = 3,itmp;
666:   PetscTruth     flg,set;
667:   PetscReal      dt[3];
668:   char           tname[256];
669:   PC_ILU         *ilu = (PC_ILU*)pc->data;
670:   PetscFList     ordlist;
671:   PetscReal      tol;

674:   MatOrderingRegisterAll(PETSC_NULL);
675:   PetscOptionsHead("ILU Options");
676:     PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);
677:     if (flg) ilu->info.levels = itmp;
678:     PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);
679:     PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);
680:     ilu->info.diagonal_fill = (double) flg;
681:     PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);
682:     PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);
683:     PetscOptionsName("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",&flg);
684:     if (flg) {
685:       PCILUSetDamping(pc,(PetscReal) PETSC_DECIDE);
686:     }
687:     PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);
688:     PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);
689:     if (flg) {
690:       PetscOptionsInt("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",(PetscInt)ilu->info.shift,&itmp,&flg);
691:       if (flg && !itmp) {
692:         PCILUSetShift(pc,PETSC_FALSE);
693:       } else {
694:         PCILUSetShift(pc,PETSC_TRUE);
695:       }
696:     }
697:     PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);

699:     dt[0] = ilu->info.dt;
700:     dt[1] = ilu->info.dtcol;
701:     dt[2] = ilu->info.dtcount;
702:     PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);
703:     if (flg) {
704:       PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);
705:     }
706:     PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);
707:     PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);

709:     MatGetOrderingList(&ordlist);
710:     PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);
711:     if (flg) {
712:       PCILUSetMatOrdering(pc,tname);
713:     }
714:     flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
715:     PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);
716:     if (set) {
717:       PCILUSetPivotInBlocks(pc,flg);
718:     }
719:   PetscOptionsTail();
720:   return(0);
721: }

725: static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
726: {
727:   PC_ILU         *ilu = (PC_ILU*)pc->data;
729:   PetscTruth     isstring,iascii;

732:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
733:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
734:   if (iascii) {
735:     if (ilu->usedt) {
736:         PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %g\n",ilu->info.dt);
737:         PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);
738:         PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %g\n",ilu->info.dtcol);
739:     } else if (ilu->info.levels == 1) {
740:         PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)ilu->info.levels);
741:     } else {
742:         PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);
743:     }
744:     PetscViewerASCIIPrintf(viewer,"  ILU: max fill ratio allocated %g\n",ilu->info.fill);
745:     PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);
746:     if (ilu->info.shift) {PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");}
747:     if (ilu->inplace) {PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");}
748:     else              {PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");}
749:     PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",ilu->ordering);
750:     if (ilu->reusefill)     {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");}
751:     if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");}
752:     if (ilu->fact) {
753:       PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");
754:       PetscViewerASCIIPushTab(viewer);
755:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
756:       MatView(ilu->fact,viewer);
757:       PetscViewerPopFormat(viewer);
758:       PetscViewerASCIIPopTab(viewer);
759:     }
760:   } else if (isstring) {
761:     PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);
762:   } else {
763:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
764:   }
765:   return(0);
766: }

770: static PetscErrorCode PCSetUp_ILU(PC pc)
771: {
773:   PetscTruth     flg;
774:   PC_ILU         *ilu = (PC_ILU*)pc->data;

777:   if (ilu->inplace) {
778:     if (!pc->setupcalled) {

780:       /* In-place factorization only makes sense with the natural ordering,
781:          so we only need to get the ordering once, even if nonzero structure changes */
782:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
783:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
784:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
785:     }

787:     /* In place ILU only makes sense with fill factor of 1.0 because 
788:        cannot have levels of fill */
789:     ilu->info.fill          = 1.0;
790:     ilu->info.diagonal_fill = 0;
791:     MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);
792:     ilu->fact = pc->pmat;
793:   } else if (ilu->usedt) {
794:     if (!pc->setupcalled) {
795:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
796:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
797:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
798:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
799:       PetscLogObjectParent(pc,ilu->fact);
800:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
801:       MatDestroy(ilu->fact);
802:       if (!ilu->reuseordering) {
803:         if (ilu->row) {ISDestroy(ilu->row);}
804:         if (ilu->col) {ISDestroy(ilu->col);}
805:         MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
806:         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
807:         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
808:       }
809:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
810:       PetscLogObjectParent(pc,ilu->fact);
811:     } else if (!ilu->reusefill) {
812:       MatDestroy(ilu->fact);
813:       MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
814:       PetscLogObjectParent(pc,ilu->fact);
815:     } else {
816:       MatLUFactorNumeric(pc->pmat,&ilu->fact);
817:     }
818:    } else {
819:     if (!pc->setupcalled) {
820:       /* first time in so compute reordering and symbolic factorization */
821:       MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
822:       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
823:       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
824:       /*  Remove zeros along diagonal?     */
825:       PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
826:       if (flg) {
827:         PetscReal ntol = 1.e-10;
828:         PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
829:         MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
830:       }

832:       MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
833:       PetscLogObjectParent(pc,ilu->fact);
834:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
835:       if (!ilu->reuseordering) {
836:         /* compute a new ordering for the ILU */
837:         ISDestroy(ilu->row);
838:         ISDestroy(ilu->col);
839:         MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
840:         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
841:         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
842:         /*  Remove zeros along diagonal?     */
843:         PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
844:         if (flg) {
845:           PetscReal ntol = 1.e-10;
846:           PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
847:           MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
848:         }
849:       }
850:       MatDestroy(ilu->fact);
851:       MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
852:       PetscLogObjectParent(pc,ilu->fact);
853:     }
854:     MatLUFactorNumeric(pc->pmat,&ilu->fact);
855:   }
856:   return(0);
857: }

861: static PetscErrorCode PCDestroy_ILU(PC pc)
862: {
863:   PC_ILU         *ilu = (PC_ILU*)pc->data;

867:   PCDestroy_ILU_Internal(pc);
868:   PetscStrfree(ilu->ordering);
869:   PetscFree(ilu);
870:   return(0);
871: }

875: static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
876: {
877:   PC_ILU         *ilu = (PC_ILU*)pc->data;

881:   MatSolve(ilu->fact,x,y);
882:   return(0);
883: }

887: static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
888: {
889:   PC_ILU         *ilu = (PC_ILU*)pc->data;

893:   MatSolveTranspose(ilu->fact,x,y);
894:   return(0);
895: }

899: static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
900: {
901:   PC_ILU *ilu = (PC_ILU*)pc->data;

904:   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
905:   *mat = ilu->fact;
906:   return(0);
907: }

909: /*MC
910:      PCILU - Incomplete factorization preconditioners.

912:    Options Database Keys:
913: +  -pc_ilu_levels <k> - number of levels of fill for ILU(k)
914: .  -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
915:                       its factorization (overwrites original matrix)
916: .  -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
917: .  -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization
918: .  -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots
919: .  -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner
920: .  -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot
921: .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
922: .  -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
923: .  -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
924:                                    this decreases the chance of getting a zero pivot
925: .  -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
926: -  -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
927:                              than 1 the diagonal blocks are factored with partial pivoting (this increases the 
928:                              stability of the ILU factorization

930:    Level: beginner

932:   Concepts: incomplete factorization

934:    Notes: Only implemented for some matrix formats. Not implemented in parallel

936:           For BAIJ matrices this implements a point block ILU

938: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
939:            PCILUSetSetZeroPivot(), PCILUSetDamping(), PCILUSetShift(), PCILUSetUseDropTolerance(),
940:            PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(),
941:            PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(),

943: M*/

948: PetscErrorCode PCCreate_ILU(PC pc)
949: {
951:   PC_ILU         *ilu;

954:   PetscNew(PC_ILU,&ilu);
955:   PetscLogObjectMemory(pc,sizeof(PC_ILU));

957:   ilu->fact                    = 0;
958:   MatFactorInfoInitialize(&ilu->info);
959:   ilu->info.levels             = 0;
960:   ilu->info.fill               = 1.0;
961:   ilu->col                     = 0;
962:   ilu->row                     = 0;
963:   ilu->inplace                 = PETSC_FALSE;
964:   PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);
965:   ilu->reuseordering           = PETSC_FALSE;
966:   ilu->usedt                   = PETSC_FALSE;
967:   ilu->info.dt                 = PETSC_DEFAULT;
968:   ilu->info.dtcount            = PETSC_DEFAULT;
969:   ilu->info.dtcol              = PETSC_DEFAULT;
970:   ilu->info.damping            = 0.0;
971:   ilu->info.shift              = PETSC_FALSE;
972:   ilu->info.shift_fraction     = 0.0;
973:   ilu->info.zeropivot          = 1.e-12;
974:   ilu->info.pivotinblocks      = 1.0;
975:   ilu->reusefill               = PETSC_FALSE;
976:   ilu->info.diagonal_fill      = 0;
977:   pc->data                     = (void*)ilu;

979:   pc->ops->destroy             = PCDestroy_ILU;
980:   pc->ops->apply               = PCApply_ILU;
981:   pc->ops->applytranspose      = PCApplyTranspose_ILU;
982:   pc->ops->setup               = PCSetUp_ILU;
983:   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
984:   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ILU;
985:   pc->ops->view                = PCView_ILU;
986:   pc->ops->applyrichardson     = 0;

988:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
989:                     PCILUSetUseDropTolerance_ILU);
990:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
991:                     PCILUSetFill_ILU);
992:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
993:                     PCILUSetDamping_ILU);
994:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU",
995:                      PCILUSetShift_ILU);
996:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
997:                     PCILUSetMatOrdering_ILU);
998:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
999:                     PCILUSetReuseOrdering_ILU);
1000:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
1001:                     PCILUDTSetReuseFill_ILUDT);
1002:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
1003:                     PCILUSetLevels_ILU);
1004:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
1005:                     PCILUSetUseInPlace_ILU);
1006:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
1007:                     PCILUSetAllowDiagonalFill_ILU);
1008:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
1009:                     PCILUSetPivotInBlocks_ILU);
1010:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU",
1011:                     PCILUSetZeroPivot_ILU);
1012:   return(0);
1013: }