Actual source code: sor.c

  1: #define PETSCKSP_DLL

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

  8: typedef struct {
  9:   PetscInt    its;        /* inner iterations, number of sweeps */
 10:   PetscInt    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:   PetscReal   fshift;
 14: } PC_SOR;

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

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

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

 37:   MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);
 38:   return(0);
 39: }

 43: static PetscErrorCode PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
 44: {
 45:   PC_SOR         *jac = (PC_SOR*)pc->data;

 49:   PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);
 50:   its  = its*jac->its;
 51:   MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,jac->fshift,its,jac->lits,y);
 52:   return(0);
 53: }

 57: PetscErrorCode PCSetFromOptions_SOR(PC pc)
 58: {
 59:   PC_SOR         *jac = (PC_SOR*)pc->data;
 61:   PetscTruth     flg;

 64:   PetscOptionsHead("(S)SOR options");
 65:     PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);
 66:     PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);
 67:     PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);
 68:     PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);
 69:     PetscOptionsTruthGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);
 70:     if (flg) {PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);}
 71:     PetscOptionsTruthGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);
 72:     if (flg) {PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);}
 73:     PetscOptionsTruthGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);
 74:     if (flg) {PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);}
 75:     PetscOptionsTruthGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);
 76:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);}
 77:     PetscOptionsTruthGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);
 78:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);}
 79:     PetscOptionsTruthGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);
 80:     if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);}
 81:   PetscOptionsTail();
 82:   return(0);
 83: }

 87: PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
 88: {
 89:   PC_SOR         *jac = (PC_SOR*)pc->data;
 90:   MatSORType     sym = jac->sym;
 91:   const char     *sortype;
 93:   PetscTruth     iascii;

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


119: /* ------------------------------------------------------------------------------*/
123: PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
124: {
125:   PC_SOR *jac;

128:   jac = (PC_SOR*)pc->data;
129:   jac->sym = flag;
130:   return(0);
131: }

137: PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
138: {
139:   PC_SOR *jac;

142:   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
143:   jac        = (PC_SOR*)pc->data;
144:   jac->omega = omega;
145:   return(0);
146: }

152: PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
153: {
154:   PC_SOR *jac;

157:   jac      = (PC_SOR*)pc->data;
158:   jac->its = its;
159:   jac->lits = lits;
160:   return(0);
161: }

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

172:    Collective on PC

174:    Input Parameters:
175: +  pc - the preconditioner context
176: -  flag - one of the following
177: .vb
178:     SOR_FORWARD_SWEEP
179:     SOR_BACKWARD_SWEEP
180:     SOR_SYMMETRIC_SWEEP
181:     SOR_LOCAL_FORWARD_SWEEP
182:     SOR_LOCAL_BACKWARD_SWEEP
183:     SOR_LOCAL_SYMMETRIC_SWEEP
184: .ve

186:    Options Database Keys:
187: +  -pc_sor_symmetric - Activates symmetric version
188: .  -pc_sor_backward - Activates backward version
189: .  -pc_sor_local_forward - Activates local forward version
190: .  -pc_sor_local_symmetric - Activates local symmetric version
191: -  -pc_sor_local_backward - Activates local backward version

193:    Notes: 
194:    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
195:    which can be chosen with the option 
196: .  -pc_type eisenstat - Activates Eisenstat trick

198:    Level: intermediate

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

202: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
203: @*/
204: PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
205: {
206:   PetscErrorCode ierr,(*f)(PC,MatSORType);

210:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);
211:   if (f) {
212:     (*f)(pc,flag);
213:   }
214:   return(0);
215: }

219: /*@
220:    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
221:    (where omega = 1.0 by default).

223:    Collective on PC

225:    Input Parameters:
226: +  pc - the preconditioner context
227: -  omega - relaxation coefficient (0 < omega < 2). 

229:    Options Database Key:
230: .  -pc_sor_omega <omega> - Sets omega

232:    Level: intermediate

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

236: .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
237: @*/
238: PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
239: {
240:   PetscErrorCode ierr,(*f)(PC,PetscReal);

243:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);
244:   if (f) {
245:     (*f)(pc,omega);
246:   }
247:   return(0);
248: }

252: /*@
253:    PCSORSetIterations - Sets the number of inner iterations to 
254:    be used by the SOR preconditioner. The default is 1.

256:    Collective on PC

258:    Input Parameters:
259: +  pc - the preconditioner context
260: .  lits - number of local iterations, smoothings over just variables on processor
261: -  its - number of parallel iterations to use; each parallel iteration has lits local iterations

263:    Options Database Key:
264: +  -pc_sor_its <its> - Sets number of iterations
265: -  -pc_sor_lits <lits> - Sets number of local iterations

267:    Level: intermediate

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

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

273: .seealso: PCSORSetOmega(), PCSORSetSymmetric()
274: @*/
275: PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
276: {
277:   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt);

281:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);
282:   if (f) {
283:     (*f)(pc,its,lits);
284:   }
285:   return(0);
286: }

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

291:    Options Database Keys:
292: +  -pc_sor_symmetric - Activates symmetric version
293: .  -pc_sor_backward - Activates backward version
294: .  -pc_sor_forward - Activates forward version
295: .  -pc_sor_local_forward - Activates local forward version
296: .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
297: .  -pc_sor_local_backward - Activates local backward version
298: .  -pc_sor_omega <omega> - Sets omega
299: .  -pc_sor_its <its> - Sets number of iterations   (default 1)
300: -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)

302:    Level: beginner

304:   Concepts: SOR, preconditioners, Gauss-Seidel

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

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

312: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
313:            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
314: M*/

319: PetscErrorCode  PCCreate_SOR(PC pc)
320: {
322:   PC_SOR         *jac;

325:   PetscNew(PC_SOR,&jac);
326:   PetscLogObjectMemory(pc,sizeof(PC_SOR));

328:   pc->ops->apply           = PCApply_SOR;
329:   pc->ops->applyrichardson = PCApplyRichardson_SOR;
330:   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
331:   pc->ops->setup           = 0;
332:   pc->ops->view            = PCView_SOR;
333:   pc->ops->destroy         = PCDestroy_SOR;
334:   pc->data                 = (void*)jac;
335:   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
336:   jac->omega               = 1.0;
337:   jac->fshift              = 0.0;
338:   jac->its                 = 1;
339:   jac->lits                = 1;

341:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
342:                     PCSORSetSymmetric_SOR);
343:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
344:                     PCSORSetOmega_SOR);
345:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
346:                     PCSORSetIterations_SOR);

348:   return(0);
349: }