Actual source code: sor.c

petsc-3.3-p7 2013-05-11
  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) {
 50:     stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
 51:   }
 52:   MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);
 53:   *outits = its;
 54:   *reason = PCRICHARDSON_CONVERGED_ITS;
 55:   return(0);
 56: }

 60: PetscErrorCode PCSetFromOptions_SOR(PC pc)
 61: {
 62:   PC_SOR         *jac = (PC_SOR*)pc->data;
 64:   PetscBool      flg;

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

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

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


122: /* ------------------------------------------------------------------------------*/
123: EXTERN_C_BEGIN
126: PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
127: {
128:   PC_SOR *jac;

131:   jac = (PC_SOR*)pc->data;
132:   jac->sym = flag;
133:   return(0);
134: }
135: EXTERN_C_END

137: EXTERN_C_BEGIN
140: PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
141: {
142:   PC_SOR *jac;

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

152: EXTERN_C_BEGIN
155: PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
156: {
157:   PC_SOR *jac;

160:   jac      = (PC_SOR*)pc->data;
161:   jac->its = its;
162:   jac->lits = lits;
163:   return(0);
164: }
165: EXTERN_C_END

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

175:    Logically Collective on PC

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

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

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

201:    Level: intermediate

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

205: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
206: @*/
207: PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
208: {

214:   PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));
215:   return(0);
216: }

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

224:    Logically 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: {

246:   PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));
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:    Logically 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: {

282:   PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));
283:   return(0);
284: }

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

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

300:    Level: beginner

302:   Concepts: SOR, preconditioners, Gauss-Seidel

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

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

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

314: EXTERN_C_BEGIN
317: PetscErrorCode  PCCreate_SOR(PC pc)
318: {
320:   PC_SOR         *jac;

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

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

338:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
339:                     PCSORSetSymmetric_SOR);
340:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
341:                     PCSORSetOmega_SOR);
342:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
343:                     PCSORSetIterations_SOR);

345:   return(0);
346: }
347: EXTERN_C_END