Actual source code: milu.c
1: /*$Id: milu.c,v 1.30 2001/08/07 03:03:39 balay Exp $*/
3: /*
4: Contributed by Victor Eijkhout <eijkhout@cs.utk.edu>, September 1998
5: */
7: #include src/sles/pc/pcimpl.h
9: /*
10: Manteuffel variant of ILU
11: @article{Ma:incompletefactorization,
12: author = {T.A. Manteuffel},
13: title = {An incomplete factorization technique for positive definite
14: linear systems},
15: journal = {Math. Comp.},
16: volume = {34},
17: year = {1980},
18: pages = {473--497},
19: abstract = {Extension of Meyerink/vdVorst to H-matrices;
20: shifted ICCG: if $A=D-B$ (diagonal) then
21: $A(alpha)=D-{1over 1+alpha}B$; for $alphageqalpha_n>0$
22: all pivots will be positive; find $alpha_n$ by trial and error.},
23: keywords = {incomplete factorization, positive definite matrices,
24: H-matrices}
25: }
26: */
28: /****************************************************************
29: User interface routines
30: ****************************************************************/
31: int PCmILUSetLevels(PC pc,int levels)
32: {
33: PC base_pc = (PC) pc->data;
37: PCILUSetLevels(base_pc,levels);
39: return(0);
40: }
42: int PCmILUSetBaseType(PC pc,PCType type)
43: {
44: PC base_pc = (PC) pc->data;
48: PCSetType(base_pc,type);
50: return(0);
51: }
53: /****************************************************************
54: Implementation
55: ****************************************************************/
57: static int PCSetup_mILU(PC pc)
58: {
59: PC base_pc = (PC) pc->data;
60: Mat omat = pc->pmat,pmat;
61: Vec diag;
62: PetscScalar *dia;
63: PetscReal *mprop;
64: int lsize,first,last,ierr;
67: ierr = MatGetOwnershipRange(omat,&first,&last);
68: lsize = last-first;
69: PetscMalloc((lsize+1)*sizeof(PetscReal),&mprop);
70: {
71: int irow;
72: for (irow=first; irow<last; irow++) {
73: int icol,ncols,*cols; PetscScalar *vals; PetscReal mp=0.;
74: MatGetRow(omat,irow,&ncols,&cols,&vals);
75: for (icol=0; icol<ncols; icol++) {
76: if (cols[icol]==irow) {
77: mp += PetscAbsScalar(vals[icol]);
78: } else {
79: mp -= PetscAbsScalar(vals[icol]);
80: }
81: }
82: MatRestoreRow(omat,irow,&ncols,&cols,&vals);
83: mprop[irow-first] = -PetscMin(0,mp);
84: }
85: }
86: MatConvert(omat,MATSAME,&pmat);
87: VecCreateSeq(MPI_COMM_SELF,lsize,&diag);
88: MatGetDiagonal(omat,diag);
89: VecGetArray(diag,&dia);
90: PCSetOperators(base_pc,pc->mat,pmat,SAME_NONZERO_PATTERN);
91: PCSetVector(base_pc,pc->vec);
93: #define ATTEMPTS 5
94: {
95: Mat lu; Vec piv;
96: PetscScalar *elt;
97: int bd,t,try1 = 0;
98: VecDuplicate(diag,&piv);
99: do {
100: PCSetUp(base_pc);
101: PCGetFactoredMatrix(base_pc,&lu);
102: MatGetDiagonal(lu,piv);
103: VecGetArray(piv,&elt);
104: bd = 0; for (t=0; t<lsize; t++) if (PetscRealPart(elt[t]) < 0.0) bd++;
105: VecRestoreArray(piv,&elt);
106: if (bd>0) {
107: /*printf("negative pivots %dn",bd);*/
108: try1++;
109: for (t=0; t<lsize; t++) {
110: PetscScalar v = dia[t]+(mprop[t]*try1)/ATTEMPTS;
111: int row = first+t;
112: MatSetValues(pmat,1,&row,1,&row,&v,INSERT_VALUES);
113: }
114: MatAssemblyBegin(pmat,MAT_FINAL_ASSEMBLY);
115: MatAssemblyEnd(pmat,MAT_FINAL_ASSEMBLY);
116: PCSetOperators(base_pc,pc->mat,pmat,SAME_NONZERO_PATTERN);
117: }
118: } while (bd>0);
119: VecDestroy(piv);
120: }
121:
122: VecRestoreArray(diag,&dia);
123: VecDestroy(diag);
124: PetscFree(mprop);
126: return(0);
127: }
129: static int PCApply_mILU(PC pc,Vec x,Vec y)
130: {
131: PC base_pc = (PC) pc->data;
133:
135: PCApply(base_pc,x,y);
137: return(0);
138: }
140: static int PCDestroy_mILU(PC pc)
141: {
142: PC base_pc = (PC) pc->data;
144:
146: MatDestroy(base_pc->pmat);
147: PCDestroy(base_pc);
149: return(0);
150: }
152: static int PCView_mILU(PC pc,PetscViewer viewer)
153: {
154: PC base_pc = (PC) pc->data;
155: int ierr;
156: PetscTruth isascii;
157:
159: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
160: if (isascii) {
161: PetscViewerASCIIPrintf(viewer," modified ILU preconditionern");
162: PetscViewerASCIIPrintf(viewer," see src/sles/pc/milu/milu.cn");
163: PetscViewerASCIIPrintf(viewer," base PC used by mILU nextn");
164: } else {
165: SETERRQ1(1,"Viewer type %s not supported for mILU PC",((PetscObject)viewer)->type_name);
166: }
167: PCView(base_pc,viewer);
168: return(0);
169: }
171: EXTERN_C_BEGIN
172: int PCCreate_mILU(PC pc)
173: {
174: PC base_pc;
178: pc->ops->apply = PCApply_mILU;
179: pc->ops->applyrichardson = 0;
180: pc->ops->destroy = PCDestroy_mILU;
181: pc->ops->setfromoptions = 0;
182: pc->ops->setup = PCSetup_mILU;
183: pc->ops->view = PCView_mILU;
185: PCCreate(pc->comm,&base_pc);
186: PCSetType(base_pc,PCILU);
187: pc->data = (void*)base_pc;
189: return(0);
190: }
191: EXTERN_C_END