Actual source code: eisen.c

  1: /*$Id: eisen.c,v 1.113 2001/08/06 21:16:28 bsmith Exp $*/

  3: /*
  4:    Defines a  Eisenstat trick SSOR  preconditioner. This uses about 
  5:  %50 of the usual amount of floating point ops used for SSOR + Krylov 
  6:  method. But it requires actually solving the preconditioned problem 
  7:  with both left and right preconditioning. 
  8: */
 9:  #include src/sles/pc/pcimpl.h

 11: typedef struct {
 12:   Mat        shell,A;
 13:   Vec        b,diag;     /* temporary storage for true right hand side */
 14:   PetscReal  omega;
 15:   PetscTruth usediag;    /* indicates preconditioner should include diagonal scaling*/
 16: } PC_Eisenstat;


 19: static int PCMult_Eisenstat(Mat mat,Vec b,Vec x)
 20: {
 21:   int          ierr;
 22:   PC           pc;
 23:   PC_Eisenstat *eis;

 26:   MatShellGetContext(mat,(void **)&pc);
 27:   eis = (PC_Eisenstat*)pc->data;
 28:   MatRelax(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);
 29:   return(0);
 30: }

 32: static int PCApply_Eisenstat(PC pc,Vec x,Vec y)
 33: {
 34:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
 35:   int          ierr;

 38:   if (eis->usediag)  {VecPointwiseMult(x,eis->diag,y);}
 39:   else               {VecCopy(x,y);}
 40:   return(0);
 41: }

 43: static int PCPre_Eisenstat(PC pc,KSP ksp,Vec x,Vec b)
 44: {
 45:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
 46:   PetscTruth   nonzero;
 47:   int          ierr;

 50:   if (pc->mat != pc->pmat) SETERRQ(PETSC_ERR_SUP,"Cannot have different mat and pmat");
 51: 
 52:   /* swap shell matrix and true matrix */
 53:   eis->A    = pc->mat;
 54:   pc->mat   = eis->shell;

 56:   if (!eis->b) {
 57:     VecDuplicate(b,&eis->b);
 58:     PetscLogObjectParent(pc,eis->b);
 59:   }
 60: 
 61:   /* save true b, other option is to swap pointers */
 62:   VecCopy(b,eis->b);

 64:   /* if nonzero initial guess, modify x */
 65:   KSPGetInitialGuessNonzero(ksp,&nonzero);
 66:   if (nonzero) {
 67:     MatRelax(eis->A,x,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);
 68:   }

 70:   /* modify b by (L + D)^{-1} */
 71:     MatRelax(eis->A,b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS |
 72:                                         SOR_FORWARD_SWEEP),0.0,1,1,b);
 73:   return(0);
 74: }

 76: static int PCPost_Eisenstat(PC pc,KSP ksp,Vec x,Vec b)
 77: {
 78:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
 79:   int          ierr;

 82:     MatRelax(eis->A,x,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS |
 83:                                  SOR_BACKWARD_SWEEP),0.0,1,1,x);
 84:   pc->mat = eis->A;
 85:   /* get back true b */
 86:   VecCopy(eis->b,b);
 87:   return(0);
 88: }

 90: static int PCDestroy_Eisenstat(PC pc)
 91: {
 92:   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
 93:   int          ierr;

 96:   if (eis->b)     {VecDestroy(eis->b);}
 97:   if (eis->shell) {MatDestroy(eis->shell);}
 98:   if (eis->diag)  {VecDestroy(eis->diag);}
 99:   PetscFree(eis);
100:   return(0);
101: }

103: static int PCSetFromOptions_Eisenstat(PC pc)
104: {
105:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
106:   int        ierr;
107:   PetscTruth flg;

110:   PetscOptionsHead("Eisenstat SSOR options");
111:     PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,0);
112:     PetscOptionsName("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",&flg);
113:     if (flg) {
114:       PCEisenstatNoDiagonalScaling(pc);
115:     }
116:   PetscOptionsTail();
117:   return(0);
118: }

120: static int PCView_Eisenstat(PC pc,PetscViewer viewer)
121: {
122:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
123:   int          ierr;
124:   PetscTruth   isascii;

127:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
128:   if (isascii) {
129:     PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %gn",eis->omega);
130:     if (eis->usediag) {
131:       PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)n");
132:     } else {
133:       PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scalingn");
134:     }
135:   } else {
136:     SETERRQ1(1,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
137:   }
138:   return(0);
139: }

141: static int PCSetUp_Eisenstat(PC pc)
142: {
143:   int          ierr,M,N,m,n;
144:   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;

147:   if (!pc->setupcalled) {
148:     MatGetSize(pc->mat,&M,&N);
149:     MatGetLocalSize(pc->mat,&m,&n);
150:     MatCreateShell(pc->comm,m,N,M,N,(void*)pc,&eis->shell);
151:     PetscLogObjectParent(pc,eis->shell);
152:     MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);
153:   }
154:   if (!eis->usediag) return(0);
155:   if (!pc->setupcalled) {
156:     VecDuplicate(pc->vec,&eis->diag);
157:     PetscLogObjectParent(pc,eis->diag);
158:   }
159:   MatGetDiagonal(pc->pmat,eis->diag);
160:   return(0);
161: }

163: /* --------------------------------------------------------------------*/

165: EXTERN_C_BEGIN
166: int PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
167: {
168:   PC_Eisenstat  *eis;

171:   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
172:   eis = (PC_Eisenstat*)pc->data;
173:   eis->omega = omega;
174:   return(0);
175: }
176: EXTERN_C_END

178: EXTERN_C_BEGIN
179: int PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
180: {
181:   PC_Eisenstat *eis;

184:   eis = (PC_Eisenstat*)pc->data;
185:   eis->usediag = PETSC_FALSE;
186:   return(0);
187: }
188: EXTERN_C_END

190: /*@ 
191:    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
192:    to use with Eisenstat's trick (where omega = 1.0 by default).

194:    Collective on PC

196:    Input Parameters:
197: +  pc - the preconditioner context
198: -  omega - relaxation coefficient (0 < omega < 2)

200:    Options Database Key:
201: .  -pc_eisenstat_omega <omega> - Sets omega

203:    Notes: 
204:    The Eisenstat trick implementation of SSOR requires about 50% of the
205:    usual amount of floating point operations used for SSOR + Krylov method;
206:    however, the preconditioned problem must be solved with both left 
207:    and right preconditioning.

209:    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 
210:    which can be chosen with the database options
211: $    -pc_type  sor  -pc_sor_symmetric

213:    Level: intermediate

215: .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega

217: .seealso: PCSORSetOmega()
218: @*/
219: int PCEisenstatSetOmega(PC pc,PetscReal omega)
220: {
221:   int ierr,(*f)(PC,PetscReal);

225:   PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);
226:   if (f) {
227:     (*f)(pc,omega);
228:   }
229:   return(0);
230: }

232: /*@
233:    PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
234:    not to do additional diagonal preconditioning. For matrices with a constant 
235:    along the diagonal, this may save a small amount of work.

237:    Collective on PC

239:    Input Parameter:
240: .  pc - the preconditioner context

242:    Options Database Key:
243: .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()

245:    Level: intermediate

247:    Note:
248:      If you use the SLESSetDiagonalScaling() or -sles_diagonal_scale option then you will
249:    likley want to use this routine since it will save you some unneeded flops.

251: .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR

253: .seealso: PCEisenstatSetOmega()
254: @*/
255: int PCEisenstatNoDiagonalScaling(PC pc)
256: {
257:   int ierr,(*f)(PC);

261:   PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);
262:   if (f) {
263:     (*f)(pc);
264:   }
265:   return(0);
266: }

268: /* --------------------------------------------------------------------*/

270: EXTERN_C_BEGIN
271: int PCCreate_Eisenstat(PC pc)
272: {
273:   int          ierr;
274:   PC_Eisenstat *eis;

277:   PetscNew(PC_Eisenstat,&eis);
278:   PetscLogObjectMemory(pc,sizeof(PC_Eisenstat));

280:   pc->ops->apply           = PCApply_Eisenstat;
281:   pc->ops->presolve        = PCPre_Eisenstat;
282:   pc->ops->postsolve       = PCPost_Eisenstat;
283:   pc->ops->applyrichardson = 0;
284:   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
285:   pc->ops->destroy         = PCDestroy_Eisenstat;
286:   pc->ops->view            = PCView_Eisenstat;
287:   pc->ops->setup           = PCSetUp_Eisenstat;

289:   pc->data           = (void*)eis;
290:   eis->omega         = 1.0;
291:   eis->b             = 0;
292:   eis->diag          = 0;
293:   eis->usediag       = PETSC_TRUE;

295:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
296:                     PCEisenstatSetOmega_Eisenstat);
297:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
298:                     "PCEisenstatNoDiagonalScaling_Eisenstat",
299:                     PCEisenstatNoDiagonalScaling_Eisenstat);
300:  return(0);
301: }
302: EXTERN_C_END