xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 3e7ff0edd3573be01c8c0fa32db97c3db8fa5c8d)
1 
2 /*
3    Defines a  (S)SOR  preconditioner for any Mat implementation
4 */
5 #include <petsc-private/pcimpl.h>               /*I "petscpc.h" I*/
6 
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;
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "PCDestroy_SOR"
17 static PetscErrorCode PCDestroy_SOR(PC pc)
18 {
19   PetscErrorCode ierr;
20 
21   PetscFunctionBegin;
22   ierr = PetscFree(pc->data);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PCApply_SOR"
28 static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
29 {
30   PC_SOR         *jac = (PC_SOR*)pc->data;
31   PetscErrorCode ierr;
32   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
33 
34   PetscFunctionBegin;
35   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
36   PetscFunctionReturn(0);
37 }
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "PCApplyRichardson_SOR"
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;
44   PetscErrorCode ierr;
45   MatSORType     stype = jac->sym;
46 
47   PetscFunctionBegin;
48   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
49   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
50   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
51   *outits = its;
52   *reason = PCRICHARDSON_CONVERGED_ITS;
53   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "PCSetFromOptions_SOR"
58 PetscErrorCode PCSetFromOptions_SOR(PC pc)
59 {
60   PC_SOR         *jac = (PC_SOR*)pc->data;
61   PetscErrorCode ierr;
62   PetscBool      flg;
63 
64   PetscFunctionBegin;
65   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
66   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
67   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
68   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
69   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
70   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
71   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
72   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
73   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
74   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
75   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
76   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
77   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
78   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
79   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
80   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
81   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
82   ierr = PetscOptionsTail();CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "PCView_SOR"
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;
93   PetscErrorCode ierr;
94   PetscBool      iascii;
95 
96   PetscFunctionBegin;
97   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
98   if (iascii) {
99     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
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     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 
116 /* ------------------------------------------------------------------------------*/
117 #undef __FUNCT__
118 #define __FUNCT__ "PCSORSetSymmetric_SOR"
119 static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
120 {
121   PC_SOR *jac = (PC_SOR*)pc->data;
122 
123   PetscFunctionBegin;
124   jac->sym = flag;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "PCSORSetOmega_SOR"
130 static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
131 {
132   PC_SOR *jac = (PC_SOR*)pc->data;
133 
134   PetscFunctionBegin;
135   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
136   jac->omega = omega;
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "PCSORSetIterations_SOR"
142 static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
143 {
144   PC_SOR *jac = (PC_SOR*)pc->data;
145 
146   PetscFunctionBegin;
147   jac->its  = its;
148   jac->lits = lits;
149   PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "PCSORGetSymmetric_SOR"
154 static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
155 {
156   PC_SOR *jac = (PC_SOR*)pc->data;
157 
158   PetscFunctionBegin;
159   *flag = jac->sym;
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "PCSORGetOmega_SOR"
165 static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
166 {
167   PC_SOR *jac = (PC_SOR*)pc->data;
168 
169   PetscFunctionBegin;
170   *omega = jac->omega;
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "PCSORGetIterations_SOR"
176 static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
177 {
178   PC_SOR *jac = (PC_SOR*)pc->data;
179 
180   PetscFunctionBegin;
181   if (its)  *its = jac->its;
182   if (lits) *lits = jac->lits;
183   PetscFunctionReturn(0);
184 }
185 
186 /* ------------------------------------------------------------------------------*/
187 #undef __FUNCT__
188 #define __FUNCT__ "PCSORGetSymmetric"
189 /*@
190    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
191    each processor.  By default forward relaxation is used.
192 
193    Logically Collective on PC
194 
195    Input Parameter:
196 .  pc - the preconditioner context
197 
198    Output Parameter:
199 .  flag - one of the following
200 .vb
201     SOR_FORWARD_SWEEP
202     SOR_BACKWARD_SWEEP
203     SOR_SYMMETRIC_SWEEP
204     SOR_LOCAL_FORWARD_SWEEP
205     SOR_LOCAL_BACKWARD_SWEEP
206     SOR_LOCAL_SYMMETRIC_SWEEP
207 .ve
208 
209    Options Database Keys:
210 +  -pc_sor_symmetric - Activates symmetric version
211 .  -pc_sor_backward - Activates backward version
212 .  -pc_sor_local_forward - Activates local forward version
213 .  -pc_sor_local_symmetric - Activates local symmetric version
214 -  -pc_sor_local_backward - Activates local backward version
215 
216    Notes:
217    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
218    which can be chosen with the option
219 .  -pc_type eisenstat - Activates Eisenstat trick
220 
221    Level: intermediate
222 
223 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
224 
225 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
226 @*/
227 PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
228 {
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
233   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "PCSORGetOmega"
239 /*@
240    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
241    (where omega = 1.0 by default).
242 
243    Logically Collective on PC
244 
245    Input Parameter:
246 .  pc - the preconditioner context
247 
248    Output Parameter:
249 .  omega - relaxation coefficient (0 < omega < 2).
250 
251    Options Database Key:
252 .  -pc_sor_omega <omega> - Sets omega
253 
254    Level: intermediate
255 
256 .keywords: PC, SOR, SSOR, set, relaxation, omega
257 
258 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
259 @*/
260 PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
261 {
262   PetscErrorCode ierr;
263 
264   PetscFunctionBegin;
265   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
266   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "PCSORGetIterations"
272 /*@
273    PCSORGetIterations - Gets the number of inner iterations to
274    be used by the SOR preconditioner. The default is 1.
275 
276    Logically Collective on PC
277 
278    Input Parameter:
279 .  pc - the preconditioner context
280 
281    Output Parameter:
282 +  lits - number of local iterations, smoothings over just variables on processor
283 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
284 
285    Options Database Key:
286 +  -pc_sor_its <its> - Sets number of iterations
287 -  -pc_sor_lits <lits> - Sets number of local iterations
288 
289    Level: intermediate
290 
291    Notes: When run on one processor the number of smoothings is lits*its
292 
293 .keywords: PC, SOR, SSOR, set, iterations
294 
295 .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
296 @*/
297 PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
298 {
299   PetscErrorCode ierr;
300 
301   PetscFunctionBegin;
302   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
303   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "PCSORSetSymmetric"
309 /*@
310    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
311    backward, or forward relaxation.  The local variants perform SOR on
312    each processor.  By default forward relaxation is used.
313 
314    Logically Collective on PC
315 
316    Input Parameters:
317 +  pc - the preconditioner context
318 -  flag - one of the following
319 .vb
320     SOR_FORWARD_SWEEP
321     SOR_BACKWARD_SWEEP
322     SOR_SYMMETRIC_SWEEP
323     SOR_LOCAL_FORWARD_SWEEP
324     SOR_LOCAL_BACKWARD_SWEEP
325     SOR_LOCAL_SYMMETRIC_SWEEP
326 .ve
327 
328    Options Database Keys:
329 +  -pc_sor_symmetric - Activates symmetric version
330 .  -pc_sor_backward - Activates backward version
331 .  -pc_sor_local_forward - Activates local forward version
332 .  -pc_sor_local_symmetric - Activates local symmetric version
333 -  -pc_sor_local_backward - Activates local backward version
334 
335    Notes:
336    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
337    which can be chosen with the option
338 .  -pc_type eisenstat - Activates Eisenstat trick
339 
340    Level: intermediate
341 
342 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
343 
344 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
345 @*/
346 PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
347 {
348   PetscErrorCode ierr;
349 
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
352   PetscValidLogicalCollectiveEnum(pc,flag,2);
353   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "PCSORSetOmega"
359 /*@
360    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
361    (where omega = 1.0 by default).
362 
363    Logically Collective on PC
364 
365    Input Parameters:
366 +  pc - the preconditioner context
367 -  omega - relaxation coefficient (0 < omega < 2).
368 
369    Options Database Key:
370 .  -pc_sor_omega <omega> - Sets omega
371 
372    Level: intermediate
373 
374 .keywords: PC, SOR, SSOR, set, relaxation, omega
375 
376 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
377 @*/
378 PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
379 {
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384   PetscValidLogicalCollectiveReal(pc,omega,2);
385   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "PCSORSetIterations"
391 /*@
392    PCSORSetIterations - Sets the number of inner iterations to
393    be used by the SOR preconditioner. The default is 1.
394 
395    Logically Collective on PC
396 
397    Input Parameters:
398 +  pc - the preconditioner context
399 .  lits - number of local iterations, smoothings over just variables on processor
400 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
401 
402    Options Database Key:
403 +  -pc_sor_its <its> - Sets number of iterations
404 -  -pc_sor_lits <lits> - Sets number of local iterations
405 
406    Level: intermediate
407 
408    Notes: When run on one processor the number of smoothings is lits*its
409 
410 .keywords: PC, SOR, SSOR, set, iterations
411 
412 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
413 @*/
414 PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
415 {
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
420   PetscValidLogicalCollectiveInt(pc,its,2);
421   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
422   PetscFunctionReturn(0);
423 }
424 
425 /*MC
426      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
427 
428    Options Database Keys:
429 +  -pc_sor_symmetric - Activates symmetric version
430 .  -pc_sor_backward - Activates backward version
431 .  -pc_sor_forward - Activates forward version
432 .  -pc_sor_local_forward - Activates local forward version
433 .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
434 .  -pc_sor_local_backward - Activates local backward version
435 .  -pc_sor_omega <omega> - Sets omega
436 .  -pc_sor_its <its> - Sets number of iterations   (default 1)
437 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
438 
439    Level: beginner
440 
441   Concepts: SOR, preconditioners, Gauss-Seidel
442 
443    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
444           Not a true parallel SOR, in parallel this implementation corresponds to block
445           Jacobi with SOR on each block.
446 
447           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
448 
449 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
450            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
451 M*/
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "PCCreate_SOR"
455 PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
456 {
457   PetscErrorCode ierr;
458   PC_SOR         *jac;
459 
460   PetscFunctionBegin;
461   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
462 
463   pc->ops->apply           = PCApply_SOR;
464   pc->ops->applyrichardson = PCApplyRichardson_SOR;
465   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
466   pc->ops->setup           = 0;
467   pc->ops->view            = PCView_SOR;
468   pc->ops->destroy         = PCDestroy_SOR;
469   pc->data                 = (void*)jac;
470   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
471   jac->omega               = 1.0;
472   jac->fshift              = 0.0;
473   jac->its                 = 1;
474   jac->lits                = 1;
475 
476   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
477   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
478   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
479   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
480   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
481   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 
486 
487 
488 
489