Actual source code: sor.c

petsc-master 2016-07-26
Report Typos and Errors
  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;
 33:   PetscReal      fshift;

 36:   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
 37:   MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);
 38:   MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);
 39:   return(0);
 40: }

 44: 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)
 45: {
 46:   PC_SOR         *jac = (PC_SOR*)pc->data;
 48:   MatSORType     stype = jac->sym;
 49:   PetscReal      fshift;

 52:   PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);
 53:   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
 54:   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
 55:   MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);
 56:   MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);
 57:   *outits = its;
 58:   *reason = PCRICHARDSON_CONVERGED_ITS;
 59:   return(0);
 60: }

 64: PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
 65: {
 66:   PC_SOR         *jac = (PC_SOR*)pc->data;
 68:   PetscBool      flg;

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

 94: PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
 95: {
 96:   PC_SOR         *jac = (PC_SOR*)pc->data;
 97:   MatSORType     sym  = jac->sym;
 98:   const char     *sortype;
100:   PetscBool      iascii;

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


122: /* ------------------------------------------------------------------------------*/
125: static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
126: {
127:   PC_SOR *jac = (PC_SOR*)pc->data;

130:   jac->sym = flag;
131:   return(0);
132: }

136: static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
137: {
138:   PC_SOR *jac = (PC_SOR*)pc->data;

141:   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
142:   jac->omega = omega;
143:   return(0);
144: }

148: static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
149: {
150:   PC_SOR *jac = (PC_SOR*)pc->data;

153:   jac->its  = its;
154:   jac->lits = lits;
155:   return(0);
156: }

160: static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
161: {
162:   PC_SOR *jac = (PC_SOR*)pc->data;

165:   *flag = jac->sym;
166:   return(0);
167: }

171: static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
172: {
173:   PC_SOR *jac = (PC_SOR*)pc->data;

176:   *omega = jac->omega;
177:   return(0);
178: }

182: static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
183: {
184:   PC_SOR *jac = (PC_SOR*)pc->data;

187:   if (its)  *its = jac->its;
188:   if (lits) *lits = jac->lits;
189:   return(0);
190: }

192: /* ------------------------------------------------------------------------------*/
195: /*@
196:    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
197:    each processor.  By default forward relaxation is used.

199:    Logically Collective on PC

201:    Input Parameter:
202: .  pc - the preconditioner context

204:    Output Parameter:
205: .  flag - one of the following
206: .vb
207:     SOR_FORWARD_SWEEP
208:     SOR_BACKWARD_SWEEP
209:     SOR_SYMMETRIC_SWEEP
210:     SOR_LOCAL_FORWARD_SWEEP
211:     SOR_LOCAL_BACKWARD_SWEEP
212:     SOR_LOCAL_SYMMETRIC_SWEEP
213: .ve

215:    Options Database Keys:
216: +  -pc_sor_symmetric - Activates symmetric version
217: .  -pc_sor_backward - Activates backward version
218: .  -pc_sor_local_forward - Activates local forward version
219: .  -pc_sor_local_symmetric - Activates local symmetric version
220: -  -pc_sor_local_backward - Activates local backward version

222:    Notes:
223:    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
224:    which can be chosen with the option
225: .  -pc_type eisenstat - Activates Eisenstat trick

227:    Level: intermediate

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

231: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
232: @*/
233: PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
234: {

239:   PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));
240:   return(0);
241: }

245: /*@
246:    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
247:    (where omega = 1.0 by default).

249:    Logically Collective on PC

251:    Input Parameter:
252: .  pc - the preconditioner context

254:    Output Parameter:
255: .  omega - relaxation coefficient (0 < omega < 2).

257:    Options Database Key:
258: .  -pc_sor_omega <omega> - Sets omega

260:    Level: intermediate

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

264: .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
265: @*/
266: PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
267: {

272:   PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));
273:   return(0);
274: }

278: /*@
279:    PCSORGetIterations - Gets the number of inner iterations to
280:    be used by the SOR preconditioner. The default is 1.

282:    Logically Collective on PC

284:    Input Parameter:
285: .  pc - the preconditioner context

287:    Output Parameter:
288: +  lits - number of local iterations, smoothings over just variables on processor
289: -  its - number of parallel iterations to use; each parallel iteration has lits local iterations

291:    Options Database Key:
292: +  -pc_sor_its <its> - Sets number of iterations
293: -  -pc_sor_lits <lits> - Sets number of local iterations

295:    Level: intermediate

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

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

301: .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
302: @*/
303: PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
304: {

309:   PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));
310:   return(0);
311: }

315: /*@
316:    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
317:    backward, or forward relaxation.  The local variants perform SOR on
318:    each processor.  By default forward relaxation is used.

320:    Logically Collective on PC

322:    Input Parameters:
323: +  pc - the preconditioner context
324: -  flag - one of the following
325: .vb
326:     SOR_FORWARD_SWEEP
327:     SOR_BACKWARD_SWEEP
328:     SOR_SYMMETRIC_SWEEP
329:     SOR_LOCAL_FORWARD_SWEEP
330:     SOR_LOCAL_BACKWARD_SWEEP
331:     SOR_LOCAL_SYMMETRIC_SWEEP
332: .ve

334:    Options Database Keys:
335: +  -pc_sor_symmetric - Activates symmetric version
336: .  -pc_sor_backward - Activates backward version
337: .  -pc_sor_local_forward - Activates local forward version
338: .  -pc_sor_local_symmetric - Activates local symmetric version
339: -  -pc_sor_local_backward - Activates local backward version

341:    Notes:
342:    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
343:    which can be chosen with the option
344: .  -pc_type eisenstat - Activates Eisenstat trick

346:    Level: intermediate

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

350: .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
351: @*/
352: PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
353: {

359:   PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));
360:   return(0);
361: }

365: /*@
366:    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
367:    (where omega = 1.0 by default).

369:    Logically Collective on PC

371:    Input Parameters:
372: +  pc - the preconditioner context
373: -  omega - relaxation coefficient (0 < omega < 2).

375:    Options Database Key:
376: .  -pc_sor_omega <omega> - Sets omega

378:    Level: intermediate

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

382: .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
383: @*/
384: PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
385: {

391:   PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));
392:   return(0);
393: }

397: /*@
398:    PCSORSetIterations - Sets the number of inner iterations to
399:    be used by the SOR preconditioner. The default is 1.

401:    Logically Collective on PC

403:    Input Parameters:
404: +  pc - the preconditioner context
405: .  lits - number of local iterations, smoothings over just variables on processor
406: -  its - number of parallel iterations to use; each parallel iteration has lits local iterations

408:    Options Database Key:
409: +  -pc_sor_its <its> - Sets number of iterations
410: -  -pc_sor_lits <lits> - Sets number of local iterations

412:    Level: intermediate

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

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

418: .seealso: PCSORSetOmega(), PCSORSetSymmetric()
419: @*/
420: PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
421: {

427:   PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));
428:   return(0);
429: }

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

434:    Options Database Keys:
435: +  -pc_sor_symmetric - Activates symmetric version
436: .  -pc_sor_backward - Activates backward version
437: .  -pc_sor_forward - Activates forward version
438: .  -pc_sor_local_forward - Activates local forward version
439: .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
440: .  -pc_sor_local_backward - Activates local backward version
441: .  -pc_sor_omega <omega> - Sets omega
442: .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
443: .  -pc_sor_its <its> - Sets number of iterations   (default 1)
444: -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)

446:    Level: beginner

448:   Concepts: SOR, preconditioners, Gauss-Seidel

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

454:           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
455:           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
456:           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the 
457:           zero pivot.

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

461:           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected 
462:           the computation is stopped with an error

464: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
465:            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
466: M*/

470: PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
471: {
473:   PC_SOR         *jac;

476:   PetscNewLog(pc,&jac);

478:   pc->ops->apply           = PCApply_SOR;
479:   pc->ops->applyrichardson = PCApplyRichardson_SOR;
480:   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
481:   pc->ops->setup           = 0;
482:   pc->ops->view            = PCView_SOR;
483:   pc->ops->destroy         = PCDestroy_SOR;
484:   pc->data                 = (void*)jac;
485:   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
486:   jac->omega               = 1.0;
487:   jac->fshift              = 0.0;
488:   jac->its                 = 1;
489:   jac->lits                = 1;

491:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);
492:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);
493:   PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);
494:   PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);
495:   PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);
496:   PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);
497:   return(0);
498: }