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: /* ------------------------------------------------------------------------------------------*/
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:
126: if (!pc->setupcalled) {
127: PetscStrfree(dir->ordering);
128: PetscStrallocpy(ordering,&dir->ordering);
129: } else if (dir->ordering != ordering) {
130: pc->setupcalled = 0;
131: PetscStrfree(dir->ordering);
132: PetscStrallocpy(ordering,&dir->ordering);
133: /* free the data structures, then create them again */
134: PCDestroy_ILU_Internal(pc);
135: }
136: return(0);
137: }
138: EXTERN_C_END
140: EXTERN_C_BEGIN
143: int PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
144: {
145: PC_ILU *ilu;
148: ilu = (PC_ILU*)pc->data;
149: ilu->reuseordering = flag;
150: return(0);
151: }
152: EXTERN_C_END
154: EXTERN_C_BEGIN
157: int PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
158: {
159: PC_ILU *ilu;
162: ilu = (PC_ILU*)pc->data;
163: ilu->reusefill = flag;
164: if (flag) SETERRQ(1,"Not yet supported");
165: return(0);
166: }
167: EXTERN_C_END
169: EXTERN_C_BEGIN
172: int PCILUSetLevels_ILU(PC pc,int levels)
173: {
174: PC_ILU *ilu;
175: int ierr;
178: ilu = (PC_ILU*)pc->data;
180: if (!pc->setupcalled) {
181: ilu->info.levels = levels;
182: } else if (ilu->usedt == PETSC_TRUE || ilu->info.levels != levels) {
183: ilu->info.levels = levels;
184: pc->setupcalled = 0;
185: ilu->usedt = PETSC_FALSE;
186: PCDestroy_ILU_Internal(pc);
187: }
188: return(0);
189: }
190: EXTERN_C_END
192: EXTERN_C_BEGIN
195: int PCILUSetUseInPlace_ILU(PC pc)
196: {
197: PC_ILU *dir;
200: dir = (PC_ILU*)pc->data;
201: dir->inplace = PETSC_TRUE;
202: return(0);
203: }
204: EXTERN_C_END
206: EXTERN_C_BEGIN
209: int PCILUSetAllowDiagonalFill_ILU(PC pc)
210: {
211: PC_ILU *dir;
214: dir = (PC_ILU*)pc->data;
215: dir->info.diagonal_fill = 1;
216: return(0);
217: }
218: EXTERN_C_END
222: /*@
223: PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
225: Collective on PC
226:
227: Input Parameters:
228: + pc - the preconditioner context
229: - zero - all pivots smaller than this will be considered zero
231: Options Database Key:
232: . -pc_ilu_zeropivot <zero> - Sets the zero pivot size
234: Level: intermediate
236: .keywords: PC, set, factorization, direct, fill
238: .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot()
239: @*/
240: int PCILUSetZeroPivot(PC pc,PetscReal zero)
241: {
242: int ierr,(*f)(PC,PetscReal);
246: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);
247: if (f) {
248: (*f)(pc,zero);
249: }
250: return(0);
251: }
255: /*@
256: PCILUSetDamping - adds this quantity to the diagonal of the matrix during the
257: ILU numerical factorization
259: Collective on PC
260:
261: Input Parameters:
262: + pc - the preconditioner context
263: - damping - amount of damping
265: Options Database Key:
266: . -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default
268: Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
269: pivot, then the damping is doubled until this is alleviated.
271: Level: intermediate
273: .keywords: PC, set, factorization, direct, fill
275: .seealso: PCILUSetFill(), PCLUSetDamp()
276: @*/
277: int PCILUSetDamping(PC pc,PetscReal damping)
278: {
279: int ierr,(*f)(PC,PetscReal);
283: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);
284: if (f) {
285: (*f)(pc,damping);
286: }
287: return(0);
288: }
292: /*@
293: PCILUSetShift - specify whether to use Manteuffel shifting of ILU.
294: If an ILU factorisation breaks down because of nonpositive pivots,
295: adding sufficient identity to the diagonal will remedy this.
296: Setting this causes a bisection method to find the minimum shift that
297: will lead to a well-defined ILU.
299: Input parameters:
300: + pc - the preconditioner context
301: - shifting - PETSC_TRUE to set shift else PETSC_FALSE
303: Options Database Key:
304: . -pc_ilu_shift - Activate PCILUSetShift()
306: Level: intermediate
308: .keywords: PC, indefinite, factorization, incomplete, ILU
310: .seealso: PCILUSetDamping()
311: @*/
312: int PCILUSetShift(PC pc,PetscTruth shifting)
313: {
314: int ierr,(*f)(PC,PetscTruth);
318: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);
319: if (f) {
320: (*f)(pc,shifting);
321: }
322: return(0);
323: }
328: /*@
329: PCILUSetUseDropTolerance - The preconditioner will use an ILU
330: based on a drop tolerance.
332: Collective on PC
334: Input Parameters:
335: + pc - the preconditioner context
336: . dt - the drop tolerance, try from 1.e-10 to .1
337: . dtcol - tolerance for column pivot, good values [0.1 to 0.01]
338: - maxrowcount - the max number of nonzeros allowed in a row, best value
339: depends on the number of nonzeros in row of original matrix
341: Options Database Key:
342: . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
344: Level: intermediate
346: Notes:
347: This uses the iludt() code of Saad's SPARSKIT package
349: .keywords: PC, levels, reordering, factorization, incomplete, ILU
350: @*/
351: int PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,int maxrowcount)
352: {
353: int ierr,(*f)(PC,PetscReal,PetscReal,int);
357: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);
358: if (f) {
359: (*f)(pc,dt,dtcol,maxrowcount);
360: }
361: return(0);
362: }
366: /*@
367: PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
368: where fill = number nonzeros in factor/number nonzeros in original matrix.
370: Collective on PC
372: Input Parameters:
373: + pc - the preconditioner context
374: - fill - amount of expected fill
376: Options Database Key:
377: $ -pc_ilu_fill <fill>
379: Note:
380: For sparse matrix factorizations it is difficult to predict how much
381: fill to expect. By running with the option -log_info PETSc will print the
382: actual amount of fill used; allowing you to set the value accurately for
383: future runs. But default PETSc uses a value of 1.0
385: Level: intermediate
387: .keywords: PC, set, factorization, direct, fill
389: .seealso: PCLUSetFill()
390: @*/
391: int PCILUSetFill(PC pc,PetscReal fill)
392: {
393: int ierr,(*f)(PC,PetscReal);
397: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
398: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);
399: if (f) {
400: (*f)(pc,fill);
401: }
402: return(0);
403: }
407: /*@C
408: PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to
409: be used in the ILU factorization.
411: Collective on PC
413: Input Parameters:
414: + pc - the preconditioner context
415: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
417: Options Database Key:
418: . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
420: Level: intermediate
422: Notes: natural ordering is used by default
424: .seealso: PCLUSetMatOrdering()
426: .keywords: PC, ILU, set, matrix, reordering
428: @*/
429: int PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
430: {
431: int ierr,(*f)(PC,MatOrderingType);
435: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);
436: if (f) {
437: (*f)(pc,ordering);
438: }
439: return(0);
440: }
444: /*@
445: PCILUSetReuseOrdering - When similar matrices are factored, this
446: causes the ordering computed in the first factor to be used for all
447: following factors; applies to both fill and drop tolerance ILUs.
449: Collective on PC
451: Input Parameters:
452: + pc - the preconditioner context
453: - flag - PETSC_TRUE to reuse else PETSC_FALSE
455: Options Database Key:
456: . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()
458: Level: intermediate
460: .keywords: PC, levels, reordering, factorization, incomplete, ILU
462: .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
463: @*/
464: int PCILUSetReuseOrdering(PC pc,PetscTruth flag)
465: {
466: int ierr,(*f)(PC,PetscTruth);
470: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);
471: if (f) {
472: (*f)(pc,flag);
473: }
474: return(0);
475: }
479: /*@
480: PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
481: this causes later ones to use the fill computed in the initial factorization.
483: Collective on PC
485: Input Parameters:
486: + pc - the preconditioner context
487: - flag - PETSC_TRUE to reuse else PETSC_FALSE
489: Options Database Key:
490: . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()
492: Level: intermediate
494: .keywords: PC, levels, reordering, factorization, incomplete, ILU
496: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
497: @*/
498: int PCILUDTSetReuseFill(PC pc,PetscTruth flag)
499: {
500: int ierr,(*f)(PC,PetscTruth);
504: PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);
505: if (f) {
506: (*f)(pc,flag);
507: }
508: return(0);
509: }
513: /*@
514: PCILUSetLevels - Sets the number of levels of fill to use.
516: Collective on PC
518: Input Parameters:
519: + pc - the preconditioner context
520: - levels - number of levels of fill
522: Options Database Key:
523: . -pc_ilu_levels <levels> - Sets fill level
525: Level: intermediate
527: .keywords: PC, levels, fill, factorization, incomplete, ILU
528: @*/
529: int PCILUSetLevels(PC pc,int levels)
530: {
531: int ierr,(*f)(PC,int);
535: if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
536: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);
537: if (f) {
538: (*f)(pc,levels);
539: }
540: return(0);
541: }
545: /*@
546: PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be
547: treated as level 0 fill even if there is no non-zero location.
549: Collective on PC
551: Input Parameters:
552: + pc - the preconditioner context
554: Options Database Key:
555: . -pc_ilu_diagonal_fill
557: Notes:
558: Does not apply with 0 fill.
560: Level: intermediate
562: .keywords: PC, levels, fill, factorization, incomplete, ILU
563: @*/
564: int PCILUSetAllowDiagonalFill(PC pc)
565: {
566: int ierr,(*f)(PC);
570: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);
571: if (f) {
572: (*f)(pc);
573: }
574: return(0);
575: }
579: /*@
580: PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
581: Collective on PC
583: Input Parameters:
584: . pc - the preconditioner context
586: Options Database Key:
587: . -pc_ilu_in_place - Activates in-place factorization
589: Notes:
590: PCILUSetUseInPlace() is intended for use with matrix-free variants of
591: Krylov methods, or when a different matrices are employed for the linear
592: system and preconditioner, or with ASM preconditioning. Do NOT use
593: this option if the linear system
594: matrix also serves as the preconditioning matrix, since the factored
595: matrix would then overwrite the original matrix.
597: Only works well with ILU(0).
599: Level: intermediate
601: .keywords: PC, set, factorization, inplace, in-place, ILU
603: .seealso: PCLUSetUseInPlace()
604: @*/
605: int PCILUSetUseInPlace(PC pc)
606: {
607: int ierr,(*f)(PC);
611: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);
612: if (f) {
613: (*f)(pc);
614: }
615: return(0);
616: }
620: /*@
621: PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
622: with BAIJ or SBAIJ matrices
624: Collective on PC
626: Input Parameters:
627: + pc - the preconditioner context
628: - pivot - PETSC_TRUE or PETSC_FALSE
630: Options Database Key:
631: . -pc_ilu_pivot_in_blocks <true,false>
633: Level: intermediate
635: .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
636: @*/
637: int PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
638: {
639: int ierr,(*f)(PC,PetscTruth);
642: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);
643: if (f) {
644: (*f)(pc,pivot);
645: }
646: return(0);
647: }
649: /* ------------------------------------------------------------------------------------------*/
651: EXTERN_C_BEGIN
654: int PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
655: {
656: PC_ILU *dir = (PC_ILU*)pc->data;
659: dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
660: return(0);
661: }
662: EXTERN_C_END
666: static int PCSetFromOptions_ILU(PC pc)
667: {
668: int ierr,dtmax = 3,itmp;
669: PetscTruth flg,set;
670: PetscReal dt[3];
671: char tname[256];
672: PC_ILU *ilu = (PC_ILU*)pc->data;
673: PetscFList ordlist;
674: PetscReal tol;
677: MatOrderingRegisterAll(PETSC_NULL);
678: PetscOptionsHead("ILU Options");
679: PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(int)ilu->info.levels,&itmp,&flg);
680: if (flg) ilu->info.levels = itmp;
681: PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);
682: PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);
683: ilu->info.diagonal_fill = (double) flg;
684: PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);
685: PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);
686: PetscOptionsName("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",&flg);
687: if (flg) {
688: PCILUSetDamping(pc,(PetscReal) PETSC_DECIDE);
689: }
690: PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);
691: PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);
692: if (flg) {
693: PCILUSetShift(pc,PETSC_TRUE);
694: }
695: PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);
697: dt[0] = ilu->info.dt;
698: dt[1] = ilu->info.dtcol;
699: dt[2] = ilu->info.dtcount;
700: PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);
701: if (flg) {
702: PCILUSetUseDropTolerance(pc,dt[0],dt[1],(int)dt[2]);
703: }
704: PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);
705: PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);
707: MatGetOrderingList(&ordlist);
708: PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);
709: if (flg) {
710: PCILUSetMatOrdering(pc,tname);
711: }
712: flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
713: PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);
714: if (set) {
715: PCILUSetPivotInBlocks(pc,flg);
716: }
717: PetscOptionsTail();
718: return(0);
719: }
723: static int PCView_ILU(PC pc,PetscViewer viewer)
724: {
725: PC_ILU *ilu = (PC_ILU*)pc->data;
726: int ierr;
727: PetscTruth isstring,isascii;
730: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
731: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
732: if (isascii) {
733: if (ilu->usedt) {
734: PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %g\n",ilu->info.dt);
735: PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %d\n",(int)ilu->info.dtcount);
736: PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %g\n",ilu->info.dtcol);
737: } else if (ilu->info.levels == 1) {
738: PetscViewerASCIIPrintf(viewer," ILU: %d level of fill\n",(int)ilu->info.levels);
739: } else {
740: PetscViewerASCIIPrintf(viewer," ILU: %d levels of fill\n",(int)ilu->info.levels);
741: }
742: PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %g\n",ilu->info.fill);
743: PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);
744: if (ilu->info.shift) {PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");}
745: if (ilu->inplace) {PetscViewerASCIIPrintf(viewer," in-place factorization\n");}
746: else {PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");}
747: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);
748: if (ilu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
749: if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
750: if (ilu->fact) {
751: PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");
752: PetscViewerASCIIPushTab(viewer);
753: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
754: MatView(ilu->fact,viewer);
755: PetscViewerPopFormat(viewer);
756: PetscViewerASCIIPopTab(viewer);
757: }
758: } else if (isstring) {
759: PetscViewerStringSPrintf(viewer," lvls=%d,order=%s",(int)ilu->info.levels,ilu->ordering);
760: } else {
761: SETERRQ1(1,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
762: }
763: return(0);
764: }
768: static int PCSetUp_ILU(PC pc)
769: {
770: int ierr;
771: PetscTruth flg;
772: PC_ILU *ilu = (PC_ILU*)pc->data;
775: if (ilu->inplace) {
776: if (!pc->setupcalled) {
778: /* In-place factorization only makes sense with the natural ordering,
779: so we only need to get the ordering once, even if nonzero structure changes */
780: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
781: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
782: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
783: }
785: /* In place ILU only makes sense with fill factor of 1.0 because
786: cannot have levels of fill */
787: ilu->info.fill = 1.0;
788: ilu->info.diagonal_fill = 0;
789: MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);
790: ilu->fact = pc->pmat;
791: } else if (ilu->usedt) {
792: if (!pc->setupcalled) {
793: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
794: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
795: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
796: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
797: PetscLogObjectParent(pc,ilu->fact);
798: } else if (pc->flag != SAME_NONZERO_PATTERN) {
799: MatDestroy(ilu->fact);
800: if (!ilu->reuseordering) {
801: if (ilu->row) {ISDestroy(ilu->row);}
802: if (ilu->col) {ISDestroy(ilu->col);}
803: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
804: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
805: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
806: }
807: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
808: PetscLogObjectParent(pc,ilu->fact);
809: } else if (!ilu->reusefill) {
810: MatDestroy(ilu->fact);
811: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
812: PetscLogObjectParent(pc,ilu->fact);
813: } else {
814: MatLUFactorNumeric(pc->pmat,&ilu->fact);
815: }
816: } else {
817: if (!pc->setupcalled) {
818: /* first time in so compute reordering and symbolic factorization */
819: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
820: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
821: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
822: /* Remove zeros along diagonal? */
823: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
824: if (flg) {
825: PetscReal ntol = 1.e-10;
826: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
827: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
828: }
830: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
831: PetscLogObjectParent(pc,ilu->fact);
832: } else if (pc->flag != SAME_NONZERO_PATTERN) {
833: if (!ilu->reuseordering) {
834: /* compute a new ordering for the ILU */
835: ISDestroy(ilu->row);
836: ISDestroy(ilu->col);
837: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
838: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
839: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
840: /* Remove zeros along diagonal? */
841: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
842: if (flg) {
843: PetscReal ntol = 1.e-10;
844: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
845: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
846: }
847: }
848: MatDestroy(ilu->fact);
849: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
850: PetscLogObjectParent(pc,ilu->fact);
851: }
852: MatLUFactorNumeric(pc->pmat,&ilu->fact);
853: }
854: return(0);
855: }
859: static int PCDestroy_ILU(PC pc)
860: {
861: PC_ILU *ilu = (PC_ILU*)pc->data;
862: int ierr;
865: PCDestroy_ILU_Internal(pc);
866: PetscStrfree(ilu->ordering);
867: PetscFree(ilu);
868: return(0);
869: }
873: static int PCApply_ILU(PC pc,Vec x,Vec y)
874: {
875: PC_ILU *ilu = (PC_ILU*)pc->data;
876: int ierr;
879: MatSolve(ilu->fact,x,y);
880: return(0);
881: }
885: static int PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
886: {
887: PC_ILU *ilu = (PC_ILU*)pc->data;
888: int ierr;
891: MatSolveTranspose(ilu->fact,x,y);
892: return(0);
893: }
897: static int PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
898: {
899: PC_ILU *ilu = (PC_ILU*)pc->data;
902: if (!ilu->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
903: *mat = ilu->fact;
904: return(0);
905: }
907: EXTERN_C_BEGIN
910: int PCCreate_ILU(PC pc)
911: {
912: int ierr;
913: PC_ILU *ilu;
916: PetscNew(PC_ILU,&ilu);
917: PetscLogObjectMemory(pc,sizeof(PC_ILU));
919: ilu->fact = 0;
920: ilu->info.levels = 0;
921: ilu->info.fill = 1.0;
922: ilu->col = 0;
923: ilu->row = 0;
924: ilu->inplace = PETSC_FALSE;
925: PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);
926: ilu->reuseordering = PETSC_FALSE;
927: ilu->usedt = PETSC_FALSE;
928: ilu->info.dt = PETSC_DEFAULT;
929: ilu->info.dtcount = PETSC_DEFAULT;
930: ilu->info.dtcol = PETSC_DEFAULT;
931: ilu->info.damping = 0.0;
932: ilu->info.shift = PETSC_FALSE;
933: ilu->info.shift_fraction = 0.0;
934: ilu->info.zeropivot = 1.e-12;
935: ilu->info.pivotinblocks = 1.0;
936: ilu->reusefill = PETSC_FALSE;
937: ilu->info.diagonal_fill = 0;
938: pc->data = (void*)ilu;
940: pc->ops->destroy = PCDestroy_ILU;
941: pc->ops->apply = PCApply_ILU;
942: pc->ops->applytranspose = PCApplyTranspose_ILU;
943: pc->ops->setup = PCSetUp_ILU;
944: pc->ops->setfromoptions = PCSetFromOptions_ILU;
945: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU;
946: pc->ops->view = PCView_ILU;
947: pc->ops->applyrichardson = 0;
949: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
950: PCILUSetUseDropTolerance_ILU);
951: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
952: PCILUSetFill_ILU);
953: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
954: PCILUSetDamping_ILU);
955: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU",
956: PCILUSetShift_ILU);
957: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
958: PCILUSetMatOrdering_ILU);
959: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
960: PCILUSetReuseOrdering_ILU);
961: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
962: PCILUDTSetReuseFill_ILUDT);
963: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
964: PCILUSetLevels_ILU);
965: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
966: PCILUSetUseInPlace_ILU);
967: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
968: PCILUSetAllowDiagonalFill_ILU);
969: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
970: PCILUSetPivotInBlocks_ILU);
971: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU",
972: PCILUSetZeroPivot_ILU);
973: return(0);
974: }
975: EXTERN_C_END