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,PetscInt *outits,PCRichardsonConvergedReason *reason)
 44: {
 45:   PC_SOR         *jac = (PC_SOR*)pc->data;

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

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

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

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

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


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

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

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

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

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

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

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

173:    Collective on PC

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

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

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

199:    Level: intermediate

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

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

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

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

224:    Collective on PC

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

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

233:    Level: intermediate

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

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

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

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

257:    Collective on PC

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

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

268:    Level: intermediate

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

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

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

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

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

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

303:    Level: beginner

305:   Concepts: SOR, preconditioners, Gauss-Seidel

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

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

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

320: PetscErrorCode  PCCreate_SOR(PC pc)
321: {
323:   PC_SOR         *jac;

326:   PetscNewLog(pc,PC_SOR,&jac);

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: }