Actual source code: sor.c

petsc-3.4.5 2014-06-29
  2: /*
  3:    Defines a  (S)SOR  preconditioner for any Mat implementation
  4: */
  5: #include <petsc-private/pcimpl.h>               /*I "petscpc.h" I*/

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

 17: static PetscErrorCode PCDestroy_SOR(PC pc)
 18: {

 22:   PetscFree(pc->data);
 23:   return(0);
 24: }

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

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

 41: static PetscErrorCode PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
 42: {
 43:   PC_SOR         *jac = (PC_SOR*)pc->data;
 45:   MatSORType     stype = jac->sym;

 48:   PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);
 49:   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
 50:   MatSOR(pc->pmat,b,jac->omega,stype,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:   PetscBool      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:   PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);
 71:   if (flg) {PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);}
 72:   PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);
 73:   if (flg) {PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);}
 74:   PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);
 75:   if (flg) {PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);}
 76:   PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);
 77:   if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);}
 78:   PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);
 79:   if (flg) {PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);}
 80:   PetscOptionsBoolGroupEnd("-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:   PetscBool      iascii;

 97:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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)             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) sortype = "local_symmetric";
107:     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
108:     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
109:     else                                                                     sortype = "unknown";
110:     PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);
111:   }
112:   return(0);
113: }


116: /* ------------------------------------------------------------------------------*/
119: static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
120: {
121:   PC_SOR *jac;

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

131: static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
132: {
133:   PC_SOR *jac;

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

144: static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
145: {
146:   PC_SOR *jac;

149:   jac       = (PC_SOR*)pc->data;
150:   jac->its  = its;
151:   jac->lits = lits;
152:   return(0);
153: }

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

163:    Logically Collective on PC

165:    Input Parameters:
166: +  pc - the preconditioner context
167: -  flag - one of the following
168: .vb
169:     SOR_FORWARD_SWEEP
170:     SOR_BACKWARD_SWEEP
171:     SOR_SYMMETRIC_SWEEP
172:     SOR_LOCAL_FORWARD_SWEEP
173:     SOR_LOCAL_BACKWARD_SWEEP
174:     SOR_LOCAL_SYMMETRIC_SWEEP
175: .ve

177:    Options Database Keys:
178: +  -pc_sor_symmetric - Activates symmetric version
179: .  -pc_sor_backward - Activates backward version
180: .  -pc_sor_local_forward - Activates local forward version
181: .  -pc_sor_local_symmetric - Activates local symmetric version
182: -  -pc_sor_local_backward - Activates local backward version

184:    Notes:
185:    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
186:    which can be chosen with the option
187: .  -pc_type eisenstat - Activates Eisenstat trick

189:    Level: intermediate

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

193: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
194: @*/
195: PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
196: {

202:   PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));
203:   return(0);
204: }

208: /*@
209:    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
210:    (where omega = 1.0 by default).

212:    Logically Collective on PC

214:    Input Parameters:
215: +  pc - the preconditioner context
216: -  omega - relaxation coefficient (0 < omega < 2).

218:    Options Database Key:
219: .  -pc_sor_omega <omega> - Sets omega

221:    Level: intermediate

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

225: .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
226: @*/
227: PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
228: {

234:   PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));
235:   return(0);
236: }

240: /*@
241:    PCSORSetIterations - Sets the number of inner iterations to
242:    be used by the SOR preconditioner. The default is 1.

244:    Logically Collective on PC

246:    Input Parameters:
247: +  pc - the preconditioner context
248: .  lits - number of local iterations, smoothings over just variables on processor
249: -  its - number of parallel iterations to use; each parallel iteration has lits local iterations

251:    Options Database Key:
252: +  -pc_sor_its <its> - Sets number of iterations
253: -  -pc_sor_lits <lits> - Sets number of local iterations

255:    Level: intermediate

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

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

261: .seealso: PCSORSetOmega(), PCSORSetSymmetric()
262: @*/
263: PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
264: {

270:   PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));
271:   return(0);
272: }

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

277:    Options Database Keys:
278: +  -pc_sor_symmetric - Activates symmetric version
279: .  -pc_sor_backward - Activates backward version
280: .  -pc_sor_forward - Activates forward version
281: .  -pc_sor_local_forward - Activates local forward version
282: .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
283: .  -pc_sor_local_backward - Activates local backward version
284: .  -pc_sor_omega <omega> - Sets omega
285: .  -pc_sor_its <its> - Sets number of iterations   (default 1)
286: -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)

288:    Level: beginner

290:   Concepts: SOR, preconditioners, Gauss-Seidel

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

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

298: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
299:            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
300: M*/

304: PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
305: {
307:   PC_SOR         *jac;

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

312:   pc->ops->apply           = PCApply_SOR;
313:   pc->ops->applyrichardson = PCApplyRichardson_SOR;
314:   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
315:   pc->ops->setup           = 0;
316:   pc->ops->view            = PCView_SOR;
317:   pc->ops->destroy         = PCDestroy_SOR;
318:   pc->data                 = (void*)jac;
319:   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
320:   jac->omega               = 1.0;
321:   jac->fshift              = 0.0;
322:   jac->its                 = 1;
323:   jac->lits                = 1;

325:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);
326:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);
327:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);
328:   return(0);
329: }