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

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

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

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

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

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

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

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

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

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

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


117: /* ------------------------------------------------------------------------------*/
121: PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
122: {
123:   PC_SOR *jac;

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

135: PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega_SOR(PC pc,PetscReal omega)
136: {
137:   PC_SOR *jac;

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

150: PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
151: {
152:   PC_SOR *jac;

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

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

170:    Collective on PC

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

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

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

196:    Level: intermediate

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

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

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

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

221:    Collective on PC

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

227:    Options Database Key:
228: .  -pc_sor_omega <omega> - Sets omega

230:    Level: intermediate

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

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

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

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

254:    Collective on PC

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

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

265:    Level: intermediate

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

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

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

279:   PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);
280:   if (f) {
281:     (*f)(pc,its,lits);
282:   }
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*/

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

323:   PetscNew(PC_SOR,&jac);
324:   PetscLogObjectMemory(pc,sizeof(PC_SOR));

326:   pc->ops->apply           = PCApply_SOR;
327:   pc->ops->applyrichardson = PCApplyRichardson_SOR;
328:   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
329:   pc->ops->setup           = 0;
330:   pc->ops->view            = PCView_SOR;
331:   pc->ops->destroy         = PCDestroy_SOR;
332:   pc->data                 = (void*)jac;
333:   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
334:   jac->omega               = 1.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: }