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/ksp/pc/pcimpl.h
6: #include src/ksp/pc/impls/ilu/ilu.h
7: #include src/mat/matimpl.h
9: /* ------------------------------------------------------------------------------------------*/
11: EXTERN_C_BEGIN
14: int PCILUSetZeroPivot_ILU(PC pc,PetscReal z)
15: {
16: PC_ILU *lu;
19: lu = (PC_ILU*)pc->data;
20: lu->info.zeropivot = z;
21: return(0);
22: }
23: EXTERN_C_END
27: int PCDestroy_ILU_Internal(PC pc)
28: {
29: PC_ILU *ilu = (PC_ILU*)pc->data;
30: int ierr;
33: if (!ilu->inplace && ilu->fact) {MatDestroy(ilu->fact);}
34: if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(ilu->row);}
35: if (ilu->col) {ISDestroy(ilu->col);}
36: return(0);
37: }
39: EXTERN_C_BEGIN
42: int PCILUSetDamping_ILU(PC pc,PetscReal damping)
43: {
44: PC_ILU *dir;
47: dir = (PC_ILU*)pc->data;
48: if (damping == (PetscReal) PETSC_DECIDE) {
49: dir->info.damping = 1.e-12;
50: } else {
51: dir->info.damping = damping;
52: }
53: return(0);
54: }
55: EXTERN_C_END
57: EXTERN_C_BEGIN
60: int PCILUSetShift_ILU(PC pc,PetscTruth shift)
61: {
62: PC_ILU *dir;
65: dir = (PC_ILU*)pc->data;
66: dir->info.shift = shift;
67: if (shift) dir->info.shift_fraction = 0.0;
68: return(0);
69: }
70: EXTERN_C_END
72: EXTERN_C_BEGIN
75: int PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,int dtcount)
76: {
77: PC_ILU *ilu;
78: int ierr;
81: ilu = (PC_ILU*)pc->data;
83: if (!pc->setupcalled) {
84: ilu->usedt = PETSC_TRUE;
85: ilu->info.dt = dt;
86: ilu->info.dtcol = dtcol;
87: ilu->info.dtcount = dtcount;
88: ilu->info.fill = PETSC_DEFAULT;
89: } else if ((ilu->usedt == PETSC_FALSE)
90: || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount) {
91: ilu->usedt = PETSC_TRUE;
92: ilu->info.dt = dt;
93: ilu->info.dtcol = dtcol;
94: ilu->info.dtcount = dtcount;
95: ilu->info.fill = PETSC_DEFAULT;
96: pc->setupcalled = 0;
97: PCDestroy_ILU_Internal(pc);
98: }
99: return(0);
100: }
101: EXTERN_C_END
103: EXTERN_C_BEGIN
106: int PCILUSetFill_ILU(PC pc,PetscReal fill)
107: {
108: PC_ILU *dir;
111: dir = (PC_ILU*)pc->data;
112: dir->info.fill = fill;
113: return(0);
114: }
115: EXTERN_C_END
117: EXTERN_C_BEGIN
120: int PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
121: {
122: PC_ILU *dir = (PC_ILU*)pc->data;
123: int ierr;
124: PetscTruth flg;
125:
127: if (!pc->setupcalled) {
128: PetscStrfree(dir->ordering);
129: PetscStrallocpy(ordering,&dir->ordering);
130: } else {
131: PetscStrcmp(dir->ordering,ordering,&flg);
132: if (!flg) {
133: pc->setupcalled = 0;
134: PetscStrfree(dir->ordering);
135: PetscStrallocpy(ordering,&dir->ordering);
136: /* free the data structures, then create them again */
137: PCDestroy_ILU_Internal(pc);
138: }
139: }
140: return(0);
141: }
142: EXTERN_C_END
144: EXTERN_C_BEGIN
147: int PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
148: {
149: PC_ILU *ilu;
152: ilu = (PC_ILU*)pc->data;
153: ilu->reuseordering = flag;
154: return(0);
155: }
156: EXTERN_C_END
158: EXTERN_C_BEGIN
161: int PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
162: {
163: PC_ILU *ilu;
166: ilu = (PC_ILU*)pc->data;
167: ilu->reusefill = flag;
168: if (flag) SETERRQ(1,"Not yet supported");
169: return(0);
170: }
171: EXTERN_C_END
173: EXTERN_C_BEGIN
176: int PCILUSetLevels_ILU(PC pc,int levels)
177: {
178: PC_ILU *ilu;
179: int ierr;
182: ilu = (PC_ILU*)pc->data;
184: if (!pc->setupcalled) {
185: ilu->info.levels = levels;
186: } else if (ilu->usedt == PETSC_TRUE || ilu->info.levels != levels) {
187: ilu->info.levels = levels;
188: pc->setupcalled = 0;
189: ilu->usedt = PETSC_FALSE;
190: PCDestroy_ILU_Internal(pc);
191: }
192: return(0);
193: }
194: EXTERN_C_END
196: EXTERN_C_BEGIN
199: int PCILUSetUseInPlace_ILU(PC pc)
200: {
201: PC_ILU *dir;
204: dir = (PC_ILU*)pc->data;
205: dir->inplace = PETSC_TRUE;
206: return(0);
207: }
208: EXTERN_C_END
210: EXTERN_C_BEGIN
213: int PCILUSetAllowDiagonalFill_ILU(PC pc)
214: {
215: PC_ILU *dir;
218: dir = (PC_ILU*)pc->data;
219: dir->info.diagonal_fill = 1;
220: return(0);
221: }
222: EXTERN_C_END
226: /*@
227: PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
229: Collective on PC
230:
231: Input Parameters:
232: + pc - the preconditioner context
233: - zero - all pivots smaller than this will be considered zero
235: Options Database Key:
236: . -pc_ilu_zeropivot <zero> - Sets the zero pivot size
238: Level: intermediate
240: .keywords: PC, set, factorization, direct, fill
242: .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot()
243: @*/
244: int PCILUSetZeroPivot(PC pc,PetscReal zero)
245: {
246: int ierr,(*f)(PC,PetscReal);
250: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);
251: if (f) {
252: (*f)(pc,zero);
253: }
254: return(0);
255: }
259: /*@
260: PCILUSetDamping - adds this quantity to the diagonal of the matrix during the
261: ILU numerical factorization
263: Collective on PC
264:
265: Input Parameters:
266: + pc - the preconditioner context
267: - damping - amount of damping
269: Options Database Key:
270: . -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default
272: Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
273: pivot, then the damping is doubled until this is alleviated.
275: Level: intermediate
277: .keywords: PC, set, factorization, direct, fill
279: .seealso: PCILUSetFill(), PCLUSetDamping(), PCILUSetShift()
280: @*/
281: int PCILUSetDamping(PC pc,PetscReal damping)
282: {
283: int ierr,(*f)(PC,PetscReal);
287: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);
288: if (f) {
289: (*f)(pc,damping);
290: }
291: return(0);
292: }
296: /*@
297: PCILUSetShift - specify whether to use Manteuffel shifting of ILU.
298: If an ILU factorisation breaks down because of nonpositive pivots,
299: adding sufficient identity to the diagonal will remedy this.
300: Setting this causes a bisection method to find the minimum shift that
301: will lead to a well-defined ILU.
303: Input parameters:
304: + pc - the preconditioner context
305: - shifting - PETSC_TRUE to set shift else PETSC_FALSE
307: Options Database Key:
308: . -pc_ilu_shift - Activate PCILUSetShift()
310: Level: intermediate
312: .keywords: PC, indefinite, factorization, incomplete, ILU
314: .seealso: PCILUSetDamping()
315: @*/
316: int PCILUSetShift(PC pc,PetscTruth shifting)
317: {
318: int ierr,(*f)(PC,PetscTruth);
322: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);
323: if (f) {
324: (*f)(pc,shifting);
325: }
326: return(0);
327: }
332: /*@
333: PCILUSetUseDropTolerance - The preconditioner will use an ILU
334: based on a drop tolerance.
336: Collective on PC
338: Input Parameters:
339: + pc - the preconditioner context
340: . dt - the drop tolerance, try from 1.e-10 to .1
341: . dtcol - tolerance for column pivot, good values [0.1 to 0.01]
342: - maxrowcount - the max number of nonzeros allowed in a row, best value
343: depends on the number of nonzeros in row of original matrix
345: Options Database Key:
346: . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
348: Level: intermediate
350: Notes:
351: This uses the iludt() code of Saad's SPARSKIT package
353: .keywords: PC, levels, reordering, factorization, incomplete, ILU
354: @*/
355: int PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,int maxrowcount)
356: {
357: int ierr,(*f)(PC,PetscReal,PetscReal,int);
361: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);
362: if (f) {
363: (*f)(pc,dt,dtcol,maxrowcount);
364: }
365: return(0);
366: }
370: /*@
371: PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
372: where fill = number nonzeros in factor/number nonzeros in original matrix.
374: Collective on PC
376: Input Parameters:
377: + pc - the preconditioner context
378: - fill - amount of expected fill
380: Options Database Key:
381: $ -pc_ilu_fill <fill>
383: Note:
384: For sparse matrix factorizations it is difficult to predict how much
385: fill to expect. By running with the option -log_info PETSc will print the
386: actual amount of fill used; allowing you to set the value accurately for
387: future runs. But default PETSc uses a value of 1.0
389: Level: intermediate
391: .keywords: PC, set, factorization, direct, fill
393: .seealso: PCLUSetFill()
394: @*/
395: int PCILUSetFill(PC pc,PetscReal fill)
396: {
397: int ierr,(*f)(PC,PetscReal);
401: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
402: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);
403: if (f) {
404: (*f)(pc,fill);
405: }
406: return(0);
407: }
411: /*@C
412: PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to
413: be used in the ILU factorization.
415: Collective on PC
417: Input Parameters:
418: + pc - the preconditioner context
419: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
421: Options Database Key:
422: . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
424: Level: intermediate
426: Notes: natural ordering is used by default
428: .seealso: PCLUSetMatOrdering()
430: .keywords: PC, ILU, set, matrix, reordering
432: @*/
433: int PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
434: {
435: int ierr,(*f)(PC,MatOrderingType);
439: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);
440: if (f) {
441: (*f)(pc,ordering);
442: }
443: return(0);
444: }
448: /*@
449: PCILUSetReuseOrdering - When similar matrices are factored, this
450: causes the ordering computed in the first factor to be used for all
451: following factors; applies to both fill and drop tolerance ILUs.
453: Collective on PC
455: Input Parameters:
456: + pc - the preconditioner context
457: - flag - PETSC_TRUE to reuse else PETSC_FALSE
459: Options Database Key:
460: . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()
462: Level: intermediate
464: .keywords: PC, levels, reordering, factorization, incomplete, ILU
466: .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
467: @*/
468: int PCILUSetReuseOrdering(PC pc,PetscTruth flag)
469: {
470: int ierr,(*f)(PC,PetscTruth);
474: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);
475: if (f) {
476: (*f)(pc,flag);
477: }
478: return(0);
479: }
483: /*@
484: PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
485: this causes later ones to use the fill computed in the initial factorization.
487: Collective on PC
489: Input Parameters:
490: + pc - the preconditioner context
491: - flag - PETSC_TRUE to reuse else PETSC_FALSE
493: Options Database Key:
494: . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()
496: Level: intermediate
498: .keywords: PC, levels, reordering, factorization, incomplete, ILU
500: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
501: @*/
502: int PCILUDTSetReuseFill(PC pc,PetscTruth flag)
503: {
504: int ierr,(*f)(PC,PetscTruth);
508: PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);
509: if (f) {
510: (*f)(pc,flag);
511: }
512: return(0);
513: }
517: /*@
518: PCILUSetLevels - Sets the number of levels of fill to use.
520: Collective on PC
522: Input Parameters:
523: + pc - the preconditioner context
524: - levels - number of levels of fill
526: Options Database Key:
527: . -pc_ilu_levels <levels> - Sets fill level
529: Level: intermediate
531: .keywords: PC, levels, fill, factorization, incomplete, ILU
532: @*/
533: int PCILUSetLevels(PC pc,int levels)
534: {
535: int ierr,(*f)(PC,int);
539: if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
540: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);
541: if (f) {
542: (*f)(pc,levels);
543: }
544: return(0);
545: }
549: /*@
550: PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be
551: treated as level 0 fill even if there is no non-zero location.
553: Collective on PC
555: Input Parameters:
556: + pc - the preconditioner context
558: Options Database Key:
559: . -pc_ilu_diagonal_fill
561: Notes:
562: Does not apply with 0 fill.
564: Level: intermediate
566: .keywords: PC, levels, fill, factorization, incomplete, ILU
567: @*/
568: int PCILUSetAllowDiagonalFill(PC pc)
569: {
570: int ierr,(*f)(PC);
574: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);
575: if (f) {
576: (*f)(pc);
577: }
578: return(0);
579: }
583: /*@
584: PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
585: Collective on PC
587: Input Parameters:
588: . pc - the preconditioner context
590: Options Database Key:
591: . -pc_ilu_in_place - Activates in-place factorization
593: Notes:
594: PCILUSetUseInPlace() is intended for use with matrix-free variants of
595: Krylov methods, or when a different matrices are employed for the linear
596: system and preconditioner, or with ASM preconditioning. Do NOT use
597: this option if the linear system
598: matrix also serves as the preconditioning matrix, since the factored
599: matrix would then overwrite the original matrix.
601: Only works well with ILU(0).
603: Level: intermediate
605: .keywords: PC, set, factorization, inplace, in-place, ILU
607: .seealso: PCLUSetUseInPlace()
608: @*/
609: int PCILUSetUseInPlace(PC pc)
610: {
611: int ierr,(*f)(PC);
615: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);
616: if (f) {
617: (*f)(pc);
618: }
619: return(0);
620: }
624: /*@
625: PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
626: with BAIJ or SBAIJ matrices
628: Collective on PC
630: Input Parameters:
631: + pc - the preconditioner context
632: - pivot - PETSC_TRUE or PETSC_FALSE
634: Options Database Key:
635: . -pc_ilu_pivot_in_blocks <true,false>
637: Level: intermediate
639: .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
640: @*/
641: int PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
642: {
643: int ierr,(*f)(PC,PetscTruth);
646: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);
647: if (f) {
648: (*f)(pc,pivot);
649: }
650: return(0);
651: }
653: /* ------------------------------------------------------------------------------------------*/
655: EXTERN_C_BEGIN
658: int PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
659: {
660: PC_ILU *dir = (PC_ILU*)pc->data;
663: dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
664: return(0);
665: }
666: EXTERN_C_END
670: static int PCSetFromOptions_ILU(PC pc)
671: {
672: int ierr,dtmax = 3,itmp;
673: PetscTruth flg,set;
674: PetscReal dt[3];
675: char tname[256];
676: PC_ILU *ilu = (PC_ILU*)pc->data;
677: PetscFList ordlist;
678: PetscReal tol;
681: MatOrderingRegisterAll(PETSC_NULL);
682: PetscOptionsHead("ILU Options");
683: PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(int)ilu->info.levels,&itmp,&flg);
684: if (flg) ilu->info.levels = itmp;
685: PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);
686: PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);
687: ilu->info.diagonal_fill = (double) flg;
688: PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);
689: PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);
690: PetscOptionsName("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",&flg);
691: if (flg) {
692: PCILUSetDamping(pc,(PetscReal) PETSC_DECIDE);
693: }
694: PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);
695: PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);
696: if (flg) {
697: PCILUSetShift(pc,PETSC_TRUE);
698: }
699: PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);
701: dt[0] = ilu->info.dt;
702: dt[1] = ilu->info.dtcol;
703: dt[2] = ilu->info.dtcount;
704: PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);
705: if (flg) {
706: PCILUSetUseDropTolerance(pc,dt[0],dt[1],(int)dt[2]);
707: }
708: PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);
709: PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);
711: MatGetOrderingList(&ordlist);
712: PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);
713: if (flg) {
714: PCILUSetMatOrdering(pc,tname);
715: }
716: flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
717: PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);
718: if (set) {
719: PCILUSetPivotInBlocks(pc,flg);
720: }
721: PetscOptionsTail();
722: return(0);
723: }
727: static int PCView_ILU(PC pc,PetscViewer viewer)
728: {
729: PC_ILU *ilu = (PC_ILU*)pc->data;
730: int ierr;
731: PetscTruth isstring,isascii;
734: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
735: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
736: if (isascii) {
737: if (ilu->usedt) {
738: PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %g\n",ilu->info.dt);
739: PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %d\n",(int)ilu->info.dtcount);
740: PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %g\n",ilu->info.dtcol);
741: } else if (ilu->info.levels == 1) {
742: PetscViewerASCIIPrintf(viewer," ILU: %d level of fill\n",(int)ilu->info.levels);
743: } else {
744: PetscViewerASCIIPrintf(viewer," ILU: %d levels of fill\n",(int)ilu->info.levels);
745: }
746: PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %g\n",ilu->info.fill);
747: PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);
748: if (ilu->info.shift) {PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");}
749: if (ilu->inplace) {PetscViewerASCIIPrintf(viewer," in-place factorization\n");}
750: else {PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");}
751: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);
752: if (ilu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
753: if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
754: if (ilu->fact) {
755: PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");
756: PetscViewerASCIIPushTab(viewer);
757: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
758: MatView(ilu->fact,viewer);
759: PetscViewerPopFormat(viewer);
760: PetscViewerASCIIPopTab(viewer);
761: }
762: } else if (isstring) {
763: PetscViewerStringSPrintf(viewer," lvls=%d,order=%s",(int)ilu->info.levels,ilu->ordering);
764: } else {
765: SETERRQ1(1,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
766: }
767: return(0);
768: }
772: static int PCSetUp_ILU(PC pc)
773: {
774: int ierr;
775: PetscTruth flg;
776: PC_ILU *ilu = (PC_ILU*)pc->data;
779: if (ilu->inplace) {
780: if (!pc->setupcalled) {
782: /* In-place factorization only makes sense with the natural ordering,
783: so we only need to get the ordering once, even if nonzero structure changes */
784: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
785: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
786: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
787: }
789: /* In place ILU only makes sense with fill factor of 1.0 because
790: cannot have levels of fill */
791: ilu->info.fill = 1.0;
792: ilu->info.diagonal_fill = 0;
793: MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);
794: ilu->fact = pc->pmat;
795: } else if (ilu->usedt) {
796: if (!pc->setupcalled) {
797: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
798: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
799: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
800: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
801: PetscLogObjectParent(pc,ilu->fact);
802: } else if (pc->flag != SAME_NONZERO_PATTERN) {
803: MatDestroy(ilu->fact);
804: if (!ilu->reuseordering) {
805: if (ilu->row) {ISDestroy(ilu->row);}
806: if (ilu->col) {ISDestroy(ilu->col);}
807: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
808: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
809: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
810: }
811: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
812: PetscLogObjectParent(pc,ilu->fact);
813: } else if (!ilu->reusefill) {
814: MatDestroy(ilu->fact);
815: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
816: PetscLogObjectParent(pc,ilu->fact);
817: } else {
818: MatLUFactorNumeric(pc->pmat,&ilu->fact);
819: }
820: } else {
821: if (!pc->setupcalled) {
822: /* first time in so compute reordering and symbolic factorization */
823: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
824: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
825: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
826: /* Remove zeros along diagonal? */
827: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
828: if (flg) {
829: PetscReal ntol = 1.e-10;
830: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
831: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
832: }
834: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
835: PetscLogObjectParent(pc,ilu->fact);
836: } else if (pc->flag != SAME_NONZERO_PATTERN) {
837: if (!ilu->reuseordering) {
838: /* compute a new ordering for the ILU */
839: ISDestroy(ilu->row);
840: ISDestroy(ilu->col);
841: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
842: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
843: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
844: /* Remove zeros along diagonal? */
845: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
846: if (flg) {
847: PetscReal ntol = 1.e-10;
848: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
849: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
850: }
851: }
852: MatDestroy(ilu->fact);
853: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
854: PetscLogObjectParent(pc,ilu->fact);
855: }
856: MatLUFactorNumeric(pc->pmat,&ilu->fact);
857: }
858: return(0);
859: }
863: static int PCDestroy_ILU(PC pc)
864: {
865: PC_ILU *ilu = (PC_ILU*)pc->data;
866: int ierr;
869: PCDestroy_ILU_Internal(pc);
870: PetscStrfree(ilu->ordering);
871: PetscFree(ilu);
872: return(0);
873: }
877: static int PCApply_ILU(PC pc,Vec x,Vec y)
878: {
879: PC_ILU *ilu = (PC_ILU*)pc->data;
880: int ierr;
883: MatSolve(ilu->fact,x,y);
884: return(0);
885: }
889: static int PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
890: {
891: PC_ILU *ilu = (PC_ILU*)pc->data;
892: int ierr;
895: MatSolveTranspose(ilu->fact,x,y);
896: return(0);
897: }
901: static int PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
902: {
903: PC_ILU *ilu = (PC_ILU*)pc->data;
906: if (!ilu->fact) SETERRQ(1,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
907: *mat = ilu->fact;
908: return(0);
909: }
911: /*MC
912: PCILU - Incomplete factorization preconditioners.
914: Options Database Keys:
915: + -pc_ilu_levels <k> - number of levels of fill for ILU(k)
916: . -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
917: its factorization (overwrites original matrix)
918: . -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
919: . -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization
920: . -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots
921: . -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner
922: . -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot
923: . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
924: . -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
925: . -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
926: this decreases the chance of getting a zero pivot
927: . -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
928: - -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
929: than 1 the diagonal blocks are factored with partial pivoting (this increases the
930: stability of the ILU factorization
932: Level: beginner
934: Concepts: incomplete factorization
936: Notes: Only implemented for some matrix formats. Not implemented in parallel
938: For BAIJ matrices this implements a point block ILU
940: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
941: PCILUSetSetZeroPivot(), PCILUSetDamping(), PCILUSetShift(), PCILUSetUseDropTolerance(),
942: PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(),
943: PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(),
945: M*/
947: EXTERN_C_BEGIN
950: int PCCreate_ILU(PC pc)
951: {
952: int ierr;
953: PC_ILU *ilu;
956: PetscNew(PC_ILU,&ilu);
957: PetscLogObjectMemory(pc,sizeof(PC_ILU));
959: ilu->fact = 0;
960: ilu->info.levels = 0;
961: ilu->info.fill = 1.0;
962: ilu->col = 0;
963: ilu->row = 0;
964: ilu->inplace = PETSC_FALSE;
965: PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);
966: ilu->reuseordering = PETSC_FALSE;
967: ilu->usedt = PETSC_FALSE;
968: ilu->info.dt = PETSC_DEFAULT;
969: ilu->info.dtcount = PETSC_DEFAULT;
970: ilu->info.dtcol = PETSC_DEFAULT;
971: ilu->info.damping = 0.0;
972: ilu->info.shift = PETSC_FALSE;
973: ilu->info.shift_fraction = 0.0;
974: ilu->info.zeropivot = 1.e-12;
975: ilu->info.pivotinblocks = 1.0;
976: ilu->reusefill = PETSC_FALSE;
977: ilu->info.diagonal_fill = 0;
978: pc->data = (void*)ilu;
980: pc->ops->destroy = PCDestroy_ILU;
981: pc->ops->apply = PCApply_ILU;
982: pc->ops->applytranspose = PCApplyTranspose_ILU;
983: pc->ops->setup = PCSetUp_ILU;
984: pc->ops->setfromoptions = PCSetFromOptions_ILU;
985: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU;
986: pc->ops->view = PCView_ILU;
987: pc->ops->applyrichardson = 0;
989: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
990: PCILUSetUseDropTolerance_ILU);
991: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
992: PCILUSetFill_ILU);
993: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
994: PCILUSetDamping_ILU);
995: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU",
996: PCILUSetShift_ILU);
997: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
998: PCILUSetMatOrdering_ILU);
999: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
1000: PCILUSetReuseOrdering_ILU);
1001: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
1002: PCILUDTSetReuseFill_ILUDT);
1003: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
1004: PCILUSetLevels_ILU);
1005: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
1006: PCILUSetUseInPlace_ILU);
1007: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
1008: PCILUSetAllowDiagonalFill_ILU);
1009: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
1010: PCILUSetPivotInBlocks_ILU);
1011: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU",
1012: PCILUSetZeroPivot_ILU);
1013: return(0);
1014: }
1015: EXTERN_C_END