Actual source code: eisen.c

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

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


 20: static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
 21: {
 23:   PC             pc;
 24:   PC_Eisenstat   *eis;

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

 35: static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
 36: {
 37:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
 39:   PetscBool      hasop;

 42:   if (eis->usediag)  {
 43:     MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
 44:     if (hasop) {
 45:       MatMultDiagonalBlock(pc->pmat,x,y);
 46:     } else {
 47:       VecPointwiseMult(y,x,eis->diag);
 48:     }
 49:   } else {VecCopy(x,y);}
 50:   return(0);
 51: }

 55: static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
 56: {
 57:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
 58:   PetscBool      nonzero;

 62:   if (pc->mat != pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot have different mat and pmat");
 63: 
 64:   /* swap shell matrix and true matrix */
 65:   eis->A    = pc->mat;
 66:   pc->mat   = eis->shell;

 68:   if (!eis->b) {
 69:     VecDuplicate(b,&eis->b);
 70:     PetscLogObjectParent(pc,eis->b);
 71:   }
 72: 

 74:   /* if nonzero initial guess, modify x */
 75:   KSPGetInitialGuessNonzero(ksp,&nonzero);
 76:   if (nonzero) {
 77:     VecCopy(x,eis->b);
 78:     MatSOR(eis->A,eis->b,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);
 79:   }

 81:   /* save true b, other option is to swap pointers */
 82:   VecCopy(b,eis->b);

 84:   /* modify b by (L + D/omega)^{-1} */
 85:     MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);
 86:   return(0);
 87: }

 91: static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
 92: {
 93:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;

 97:   /* get back true b */
 98:   VecCopy(eis->b,b);

100:   /* modify x by (U + D/omega)^{-1} */
101:   VecCopy(x,eis->b);
102:   MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);
103:   pc->mat = eis->A;
104:   return(0);
105: }

109: static PetscErrorCode PCReset_Eisenstat(PC pc)
110: {
111:   PC_Eisenstat   *eis = (PC_Eisenstat *)pc->data;

115:   VecDestroy(&eis->b);
116:   MatDestroy(&eis->shell);
117:   VecDestroy(&eis->diag);
118:   return(0);
119: }

123: static PetscErrorCode PCDestroy_Eisenstat(PC pc)
124: {

128:   PCReset_Eisenstat(pc);
129:   PetscFree(pc->data);
130:   return(0);
131: }

135: static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
136: {
137:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
139:   PetscBool      flg = PETSC_FALSE;

142:   PetscOptionsHead("Eisenstat SSOR options");
143:     PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);
144:     PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);
145:     if (flg) {
146:       PCEisenstatNoDiagonalScaling(pc);
147:     }
148:   PetscOptionsTail();
149:   return(0);
150: }

154: static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
155: {
156:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
158:   PetscBool      iascii;

161:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
162:   if (iascii) {
163:     PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);
164:     if (eis->usediag) {
165:       PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");
166:     } else {
167:       PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");
168:     }
169:   } else {
170:     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
171:   }
172:   return(0);
173: }

177: static PetscErrorCode PCSetUp_Eisenstat(PC pc)
178: {
180:   PetscInt       M,N,m,n;
181:   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;

184:   if (!pc->setupcalled) {
185:     MatGetSize(pc->mat,&M,&N);
186:     MatGetLocalSize(pc->mat,&m,&n);
187:     MatCreate(((PetscObject)pc)->comm,&eis->shell);
188:     MatSetSizes(eis->shell,m,n,M,N);
189:     MatSetType(eis->shell,MATSHELL);
190:     MatShellSetContext(eis->shell,(void*)pc);
191:     PetscLogObjectParent(pc,eis->shell);
192:     MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);
193:   }
194:   if (!eis->usediag) return(0);
195:   if (!pc->setupcalled) {
196:     MatGetVecs(pc->pmat,&eis->diag,0);
197:     PetscLogObjectParent(pc,eis->diag);
198:   }
199:   MatGetDiagonal(pc->pmat,eis->diag);
200:   return(0);
201: }

203: /* --------------------------------------------------------------------*/

208: PetscErrorCode  PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
209: {
210:   PC_Eisenstat  *eis;

213:   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
214:   eis = (PC_Eisenstat*)pc->data;
215:   eis->omega = omega;
216:   return(0);
217: }

223: PetscErrorCode  PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
224: {
225:   PC_Eisenstat *eis;

228:   eis = (PC_Eisenstat*)pc->data;
229:   eis->usediag = PETSC_FALSE;
230:   return(0);
231: }

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

240:    Logically Collective on PC

242:    Input Parameters:
243: +  pc - the preconditioner context
244: -  omega - relaxation coefficient (0 < omega < 2)

246:    Options Database Key:
247: .  -pc_eisenstat_omega <omega> - Sets omega

249:    Notes: 
250:    The Eisenstat trick implementation of SSOR requires about 50% of the
251:    usual amount of floating point operations used for SSOR + Krylov method;
252:    however, the preconditioned problem must be solved with both left 
253:    and right preconditioning.

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

259:    Level: intermediate

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

263: .seealso: PCSORSetOmega()
264: @*/
265: PetscErrorCode  PCEisenstatSetOmega(PC pc,PetscReal omega)
266: {

272:   PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));
273:   return(0);
274: }

278: /*@
279:    PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
280:    not to do additional diagonal preconditioning. For matrices with a constant 
281:    along the diagonal, this may save a small amount of work.

283:    Logically Collective on PC

285:    Input Parameter:
286: .  pc - the preconditioner context

288:    Options Database Key:
289: .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()

291:    Level: intermediate

293:    Note:
294:      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
295:    likley want to use this routine since it will save you some unneeded flops.

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

299: .seealso: PCEisenstatSetOmega()
300: @*/
301: PetscErrorCode  PCEisenstatNoDiagonalScaling(PC pc)
302: {

307:   PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));
308:   return(0);
309: }

311: /* --------------------------------------------------------------------*/

313: /*MC
314:      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
315:            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.

317:    Options Database Keys:
318: +  -pc_eisenstat_omega <omega> - Sets omega
319: -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()

321:    Level: beginner

323:   Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick

325:    Notes: Only implemented for the SeqAIJ matrix format.
326:           Not a true parallel SOR, in parallel this implementation corresponds to block
327:           Jacobi with SOR on each block.

329: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
330:            PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
331: M*/

336: PetscErrorCode  PCCreate_Eisenstat(PC pc)
337: {
339:   PC_Eisenstat   *eis;

342:   PetscNewLog(pc,PC_Eisenstat,&eis);

344:   pc->ops->apply           = PCApply_Eisenstat;
345:   pc->ops->presolve        = PCPreSolve_Eisenstat;
346:   pc->ops->postsolve       = PCPostSolve_Eisenstat;
347:   pc->ops->applyrichardson = 0;
348:   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
349:   pc->ops->destroy         = PCDestroy_Eisenstat;
350:   pc->ops->reset           = PCReset_Eisenstat;
351:   pc->ops->view            = PCView_Eisenstat;
352:   pc->ops->setup           = PCSetUp_Eisenstat;

354:   pc->data           = (void*)eis;
355:   eis->omega         = 1.0;
356:   eis->b             = 0;
357:   eis->diag          = 0;
358:   eis->usediag       = PETSC_TRUE;

360:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
361:                     PCEisenstatSetOmega_Eisenstat);
362:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
363:                     "PCEisenstatNoDiagonalScaling_Eisenstat",
364:                     PCEisenstatNoDiagonalScaling_Eisenstat);
365:  return(0);
366: }