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-{1\over 1+\alpha}B$; for $\alpha\geq\alpha_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: ****************************************************************/
 33: int PCmILUSetLevels(PC pc,int levels)
 34: {
 35:   PC  base_pc = (PC) pc->data;

 39:   PCILUSetLevels(base_pc,levels);

 41:   return(0);
 42: }

 46: int PCmILUSetBaseType(PC pc,PCType type)
 47: {
 48:   PC  base_pc = (PC) pc->data;

 52:   PCSetType(base_pc,type);

 54:   return(0);
 55: }

 57: /****************************************************************
 58:   Implementation
 59: ****************************************************************/

 63: static int PCSetup_mILU(PC pc)
 64: {
 65:   PC        base_pc = (PC) pc->data;
 66:   Mat       omat = pc->pmat,pmat;
 67:   Vec       diag;
 68:   PetscScalar    *dia;
 69:   PetscReal *mprop;
 70:   int       lsize,first,last,ierr;

 73:   MatGetOwnershipRange(omat,&first,&last);
 74:   lsize = last-first;
 75:   PetscMalloc((lsize+1)*sizeof(PetscReal),&mprop);
 76:   {
 77:     int irow;
 78:     for (irow=first; irow<last; irow++) {
 79:       int icol,ncols,*cols; PetscScalar *vals; PetscReal mp=0.;
 80:       MatGetRow(omat,irow,&ncols,&cols,&vals);
 81:       for (icol=0; icol<ncols; icol++) {
 82:         if (cols[icol]==irow) {
 83:           mp += PetscAbsScalar(vals[icol]);
 84:         } else {
 85:           mp -= PetscAbsScalar(vals[icol]);
 86:         }
 87:       }
 88:       MatRestoreRow(omat,irow,&ncols,&cols,&vals);
 89:       mprop[irow-first] = -PetscMin(0,mp);
 90:     }
 91:   }
 92:   MatConvert(omat,MATSAME,&pmat);
 93:   VecCreateSeq(MPI_COMM_SELF,lsize,&diag);
 94:   MatGetDiagonal(omat,diag);
 95:   VecGetArray(diag,&dia);
 96:   PCSetOperators(base_pc,pc->mat,pmat,SAME_NONZERO_PATTERN);
 97:   PCSetVector(base_pc,pc->vec);

 99: #define ATTEMPTS 5
100:   {
101:     Mat    lu; Vec piv;
102:     PetscScalar *elt;
103:     int    bd,t,try1 = 0;
104:     VecDuplicate(diag,&piv);
105:     do {
106:       PCSetUp(base_pc);
107:       PCGetFactoredMatrix(base_pc,&lu);
108:       MatGetDiagonal(lu,piv);
109:       VecGetArray(piv,&elt);
110:       bd = 0; for (t=0; t<lsize; t++) if (PetscRealPart(elt[t]) < 0.0) bd++;
111:       VecRestoreArray(piv,&elt);
112:       if (bd>0) {
113:         /*printf("negative pivots %d\n",bd);*/
114:         try1++;
115:         for (t=0; t<lsize; t++) {
116:           PetscScalar v = dia[t]+(mprop[t]*try1)/ATTEMPTS;
117:           int row  = first+t;
118:           MatSetValues(pmat,1,&row,1,&row,&v,INSERT_VALUES);
119:         }
120:         MatAssemblyBegin(pmat,MAT_FINAL_ASSEMBLY);
121:         MatAssemblyEnd(pmat,MAT_FINAL_ASSEMBLY);
122:         PCSetOperators(base_pc,pc->mat,pmat,SAME_NONZERO_PATTERN);
123:       }
124:     } while (bd>0);
125:     VecDestroy(piv);
126:   }
127: 
128:   VecRestoreArray(diag,&dia);
129:   VecDestroy(diag);
130:   PetscFree(mprop);

132:   return(0);
133: }

137: static int PCApply_mILU(PC pc,Vec x,Vec y)
138: {
139:   PC  base_pc = (PC) pc->data;
141: 
143:   PCApply(base_pc,x,y,PC_LEFT);

145:   return(0);
146: }

150: static int PCDestroy_mILU(PC pc)
151: {
152:   PC  base_pc = (PC) pc->data;
154: 
156:   MatDestroy(base_pc->pmat);
157:   PCDestroy(base_pc);

159:   return(0);
160: }

164: static int PCView_mILU(PC pc,PetscViewer viewer)
165: {
166:   PC         base_pc = (PC) pc->data;
167:   int        ierr;
168:   PetscTruth isascii;
169: 
171:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
172:   if (isascii) {
173:     PetscViewerASCIIPrintf(viewer,"  modified ILU preconditioner\n");
174:     PetscViewerASCIIPrintf(viewer,"    see src/sles/pc/milu/milu.c\n");
175:     PetscViewerASCIIPrintf(viewer,"    base PC used by mILU next\n");
176:   } else {
177:     SETERRQ1(1,"Viewer type %s not supported for mILU PC",((PetscObject)viewer)->type_name);
178:   }
179:   PCView(base_pc,viewer);
180:   return(0);
181: }

183: EXTERN_C_BEGIN
186: int PCCreate_mILU(PC pc)
187: {
188:   PC  base_pc;

192:   pc->ops->apply            = PCApply_mILU;
193:   pc->ops->applyrichardson  = 0;
194:   pc->ops->destroy          = PCDestroy_mILU;
195:   pc->ops->setfromoptions   = 0;
196:   pc->ops->setup            = PCSetup_mILU;
197:   pc->ops->view             = PCView_mILU;

199:   PCCreate(pc->comm,&base_pc);
200:   PCSetType(base_pc,PCILU);
201:   pc->data = (void*)base_pc;

203:   return(0);
204: }
205: EXTERN_C_END