Actual source code: sor.c

  1: /*$Id: sor.c,v 1.103 2001/08/21 21:03:17 bsmith Exp $*/

  3: /*
  4:    Defines a  (S)SOR  preconditioner for any Mat implementation
  5: */
 6:  #include src/ksp/pc/pcimpl.h

  8: typedef struct {
  9:   int        its;        /* inner iterations, number of sweeps */
 10:   int        lits;       /* local inner iterations, number of sweeps applied by the local matrix mat->A */
 11:   MatSORType sym;        /* forward, reverse, symmetric etc. */
 12:   PetscReal  omega;
 13: } PC_SOR;

 17: static int PCDestroy_SOR(PC pc)
 18: {
 19:   PC_SOR *jac = (PC_SOR*)pc->data;
 20:   int    ierr;

 23:   PetscFree(jac);
 24:   return(0);
 25: }

 29: static int PCApply_SOR(PC pc,Vec x,Vec y)
 30: {
 31:   PC_SOR *jac = (PC_SOR*)pc->data;
 32:   int    ierr,flag = jac->sym | SOR_ZERO_INITIAL_GUESS;

 35:   MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,0.0,jac->its,jac->lits,y);
 36:   return(0);
 37: }

 41: static int PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its)
 42: {
 43:   PC_SOR *jac = (PC_SOR*)pc->data;
 44:   int    ierr;

 47:   PetscLogInfo(pc,"PCApplyRichardson_SOR: Warning, convergence critera ignored, using %d iterations\n",its);
 48:   its  = its*jac->its;
 49:   MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,0.0,its,jac->lits,y);
 50:   return(0);
 51: }

 55: int PCSetFromOptions_SOR(PC pc)
 56: {
 57:   PC_SOR     *jac = (PC_SOR*)pc->data;
 58:   int        ierr;
 59:   PetscTruth flg;

 62:   PetscOptionsHead("(S)SOR options");
 63:     PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);
 64:     PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);
 65:     PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);
 66:     PetscOptionsLogicalGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);
 67:     if (flg) {PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);}
 68:     PetscOptionsLogicalGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);
 69:     if (flg) {PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);}
 70:     PetscOptionsLogicalGroup("-pc_sor_local_symmetric","use SSOR seperately on each processor","PCSORSetSymmetric",&flg);
 71:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);}
 72:     PetscOptionsLogicalGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);
 73:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);}
 74:     PetscOptionsLogicalGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);
 75:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);}
 76:   PetscOptionsTail();
 77:   return(0);
 78: }

 82: int PCView_SOR(PC pc,PetscViewer viewer)
 83: {
 84:   PC_SOR     *jac = (PC_SOR*)pc->data;
 85:   MatSORType sym = jac->sym;
 86:   const char *sortype;
 87:   int        ierr;
 88:   PetscTruth isascii;

 91:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
 92:   if (isascii) {
 93:     if (sym & SOR_ZERO_INITIAL_GUESS) {PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");}
 94:     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
 95:     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
 96:     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
 97:     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
 98:                                              sortype = "symmetric";
 99:     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
100:     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
101:     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
102:                                              sortype = "local_symmetric";
103:     else if (sym & SOR_LOCAL_FORWARD_SWEEP)  sortype = "local_forward";
104:     else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward";
105:     else                                     sortype = "unknown";
106:     PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %d, omega = %g\n",sortype,jac->its,jac->omega);
107:   } else {
108:     SETERRQ1(1,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
109:   }
110:   return(0);
111: }


114: /* ------------------------------------------------------------------------------*/
115: EXTERN_C_BEGIN
118: int PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
119: {
120:   PC_SOR *jac;

123:   jac = (PC_SOR*)pc->data;
124:   jac->sym = flag;
125:   return(0);
126: }
127: EXTERN_C_END

129: EXTERN_C_BEGIN
132: int PCSORSetOmega_SOR(PC pc,PetscReal omega)
133: {
134:   PC_SOR *jac;

137:   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
138:   jac        = (PC_SOR*)pc->data;
139:   jac->omega = omega;
140:   return(0);
141: }
142: EXTERN_C_END

144: EXTERN_C_BEGIN
147: int PCSORSetIterations_SOR(PC pc,int its,int lits)
148: {
149:   PC_SOR *jac;

152:   jac      = (PC_SOR*)pc->data;
153:   jac->its = its;
154:   jac->lits = lits;
155:   return(0);
156: }
157: EXTERN_C_END

159: /* ------------------------------------------------------------------------------*/
162: /*@
163:    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 
164:    backward, or forward relaxation.  The local variants perform SOR on
165:    each processor.  By default forward relaxation is used.

167:    Collective on PC

169:    Input Parameters:
170: +  pc - the preconditioner context
171: -  flag - one of the following
172: .vb
173:     SOR_FORWARD_SWEEP
174:     SOR_BACKWARD_SWEEP
175:     SOR_SYMMETRIC_SWEEP
176:     SOR_LOCAL_FORWARD_SWEEP
177:     SOR_LOCAL_BACKWARD_SWEEP
178:     SOR_LOCAL_SYMMETRIC_SWEEP
179: .ve

181:    Options Database Keys:
182: +  -pc_sor_symmetric - Activates symmetric version
183: .  -pc_sor_backward - Activates backward version
184: .  -pc_sor_local_forward - Activates local forward version
185: .  -pc_sor_local_symmetric - Activates local symmetric version
186: -  -pc_sor_local_backward - Activates local backward version

188:    Notes: 
189:    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
190:    which can be chosen with the option 
191: .  -pc_type eisenstat - Activates Eisenstat trick

193:    Level: intermediate

195: .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric

197: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
198: @*/
199: int PCSORSetSymmetric(PC pc,MatSORType flag)
200: {
201:   int ierr,(*f)(PC,MatSORType);

205:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);
206:   if (f) {
207:     (*f)(pc,flag);
208:   }
209:   return(0);
210: }

214: /*@
215:    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
216:    (where omega = 1.0 by default).

218:    Collective on PC

220:    Input Parameters:
221: +  pc - the preconditioner context
222: -  omega - relaxation coefficient (0 < omega < 2). 

224:    Options Database Key:
225: .  -pc_sor_omega <omega> - Sets omega

227:    Level: intermediate

229: .keywords: PC, SOR, SSOR, set, relaxation, omega

231: .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
232: @*/
233: int PCSORSetOmega(PC pc,PetscReal omega)
234: {
235:   int ierr,(*f)(PC,PetscReal);

238:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);
239:   if (f) {
240:     (*f)(pc,omega);
241:   }
242:   return(0);
243: }

247: /*@
248:    PCSORSetIterations - Sets the number of inner iterations to 
249:    be used by the SOR preconditioner. The default is 1.

251:    Collective on PC

253:    Input Parameters:
254: +  pc - the preconditioner context
255: .  lits - number of local iterations, smoothings over just variables on processor
256: -  its - number of parallel iterations to use; each parallel iteration has lits local iterations

258:    Options Database Key:
259: +  -pc_sor_its <its> - Sets number of iterations
260: -  -pc_sor_lits <lits> - Sets number of local iterations

262:    Level: intermediate

264:    Notes: When run on one processor the number of smoothings is lits*its

266: .keywords: PC, SOR, SSOR, set, iterations

268: .seealso: PCSORSetOmega(), PCSORSetSymmetric()
269: @*/
270: int PCSORSetIterations(PC pc,int its,int lits)
271: {
272:   int ierr,(*f)(PC,int,int);

276:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);
277:   if (f) {
278:     (*f)(pc,its,lits);
279:   }
280:   return(0);
281: }

283: /*MC
284:      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning

286:    Options Database Keys:
287: +  -pc_sor_symmetric - Activates symmetric version
288: .  -pc_sor_backward - Activates backward version
289: .  -pc_sor_local_forward - Activates local forward version
290: .  -pc_sor_local_symmetric - Activates local symmetric version
291: .  -pc_sor_local_backward - Activates local backward version
292: .  -pc_sor_omega <omega> - Sets omega
293: .  -pc_sor_its <its> - Sets number of iterations
294: -  -pc_sor_lits <lits> - Sets number of local iterations

296:    Level: beginner

298:   Concepts: SOR, preconditioners, Gauss-Seidel

300:    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
301:           Not a true parallel SOR, in parallel this implementation corresponds to block
302:           Jacobi with SOR on each block.

304:           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.

306: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
307:            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
308: M*/

310: EXTERN_C_BEGIN
313: int PCCreate_SOR(PC pc)
314: {
315:   int    ierr;
316:   PC_SOR *jac;

319:   PetscNew(PC_SOR,&jac);
320:   PetscLogObjectMemory(pc,sizeof(PC_SOR));

322:   pc->ops->apply           = PCApply_SOR;
323:   pc->ops->applyrichardson = PCApplyRichardson_SOR;
324:   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
325:   pc->ops->setup           = 0;
326:   pc->ops->view            = PCView_SOR;
327:   pc->ops->destroy         = PCDestroy_SOR;
328:   pc->data           = (void*)jac;
329:   jac->sym           = SOR_FORWARD_SWEEP;
330:   jac->omega         = 1.0;
331:   jac->its           = 1;
332:   jac->lits          = 1;

334:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
335:                     PCSORSetSymmetric_SOR);
336:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
337:                     PCSORSetOmega_SOR);
338:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
339:                     PCSORSetIterations_SOR);

341:   return(0);
342: }
343: EXTERN_C_END