xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 9af31e4ad595286b4e2df8194fee047feeccfe42)
1 /*
2    Defines a  (S)SOR  preconditioner for any Mat implementation
3 */
4 #include "src/ksp/pc/pcimpl.h"               /*I "petscpc.h" I*/
5 
6 typedef struct {
7   int        its;        /* inner iterations, number of sweeps */
8   int        lits;       /* local inner iterations, number of sweeps applied by the local matrix mat->A */
9   MatSORType sym;        /* forward, reverse, symmetric etc. */
10   PetscReal  omega;
11 } PC_SOR;
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "PCDestroy_SOR"
15 static int PCDestroy_SOR(PC pc)
16 {
17   PC_SOR *jac = (PC_SOR*)pc->data;
18   int    ierr;
19 
20   PetscFunctionBegin;
21   ierr = PetscFree(jac);CHKERRQ(ierr);
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNCT__
26 #define __FUNCT__ "PCApply_SOR"
27 static int PCApply_SOR(PC pc,Vec x,Vec y)
28 {
29   PC_SOR *jac = (PC_SOR*)pc->data;
30   int    ierr,flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
31 
32   PetscFunctionBegin;
33   ierr = MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,0.0,jac->its,jac->lits,y);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "PCApplyRichardson_SOR"
39 static int PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its)
40 {
41   PC_SOR *jac = (PC_SOR*)pc->data;
42   int    ierr;
43 
44   PetscFunctionBegin;
45   PetscLogInfo(pc,"PCApplyRichardson_SOR: Warning, convergence critera ignored, using %d iterations\n",its);
46   its  = its*jac->its;
47   ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,0.0,its,jac->lits,y);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "PCSetFromOptions_SOR"
53 int PCSetFromOptions_SOR(PC pc)
54 {
55   PC_SOR     *jac = (PC_SOR*)pc->data;
56   int        ierr;
57   PetscTruth flg;
58 
59   PetscFunctionBegin;
60   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
61     ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr);
62     ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
63     ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr);
64     ierr = PetscOptionsLogicalGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
65     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
66     ierr = PetscOptionsLogicalGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
67     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
68     ierr = PetscOptionsLogicalGroup("-pc_sor_local_symmetric","use SSOR seperately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
69     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
70     ierr = PetscOptionsLogicalGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
71     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
72     ierr = PetscOptionsLogicalGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
73     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
74   ierr = PetscOptionsTail();CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "PCView_SOR"
80 int PCView_SOR(PC pc,PetscViewer viewer)
81 {
82   PC_SOR     *jac = (PC_SOR*)pc->data;
83   MatSORType sym = jac->sym;
84   const char *sortype;
85   int        ierr;
86   PetscTruth iascii;
87 
88   PetscFunctionBegin;
89   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
90   if (iascii) {
91     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
92     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
93     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
94     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
95     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
96                                              sortype = "symmetric";
97     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
98     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
99     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
100                                              sortype = "local_symmetric";
101     else if (sym & SOR_LOCAL_FORWARD_SWEEP)  sortype = "local_forward";
102     else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward";
103     else                                     sortype = "unknown";
104     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %d, omega = %g\n",sortype,jac->its,jac->omega);CHKERRQ(ierr);
105   } else {
106     SETERRQ1(1,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
107   }
108   PetscFunctionReturn(0);
109 }
110 
111 
112 /* ------------------------------------------------------------------------------*/
113 EXTERN_C_BEGIN
114 #undef __FUNCT__
115 #define __FUNCT__ "PCSORSetSymmetric_SOR"
116 int PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
117 {
118   PC_SOR *jac;
119 
120   PetscFunctionBegin;
121   jac = (PC_SOR*)pc->data;
122   jac->sym = flag;
123   PetscFunctionReturn(0);
124 }
125 EXTERN_C_END
126 
127 EXTERN_C_BEGIN
128 #undef __FUNCT__
129 #define __FUNCT__ "PCSORSetOmega_SOR"
130 int PCSORSetOmega_SOR(PC pc,PetscReal omega)
131 {
132   PC_SOR *jac;
133 
134   PetscFunctionBegin;
135   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
136   jac        = (PC_SOR*)pc->data;
137   jac->omega = omega;
138   PetscFunctionReturn(0);
139 }
140 EXTERN_C_END
141 
142 EXTERN_C_BEGIN
143 #undef __FUNCT__
144 #define __FUNCT__ "PCSORSetIterations_SOR"
145 int PCSORSetIterations_SOR(PC pc,int its,int lits)
146 {
147   PC_SOR *jac;
148 
149   PetscFunctionBegin;
150   jac      = (PC_SOR*)pc->data;
151   jac->its = its;
152   jac->lits = lits;
153   PetscFunctionReturn(0);
154 }
155 EXTERN_C_END
156 
157 /* ------------------------------------------------------------------------------*/
158 #undef __FUNCT__
159 #define __FUNCT__ "PCSORSetSymmetric"
160 /*@
161    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
162    backward, or forward relaxation.  The local variants perform SOR on
163    each processor.  By default forward relaxation is used.
164 
165    Collective on PC
166 
167    Input Parameters:
168 +  pc - the preconditioner context
169 -  flag - one of the following
170 .vb
171     SOR_FORWARD_SWEEP
172     SOR_BACKWARD_SWEEP
173     SOR_SYMMETRIC_SWEEP
174     SOR_LOCAL_FORWARD_SWEEP
175     SOR_LOCAL_BACKWARD_SWEEP
176     SOR_LOCAL_SYMMETRIC_SWEEP
177 .ve
178 
179    Options Database Keys:
180 +  -pc_sor_symmetric - Activates symmetric version
181 .  -pc_sor_backward - Activates backward version
182 .  -pc_sor_local_forward - Activates local forward version
183 .  -pc_sor_local_symmetric - Activates local symmetric version
184 -  -pc_sor_local_backward - Activates local backward version
185 
186    Notes:
187    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
188    which can be chosen with the option
189 .  -pc_type eisenstat - Activates Eisenstat trick
190 
191    Level: intermediate
192 
193 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
194 
195 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
196 @*/
197 int PCSORSetSymmetric(PC pc,MatSORType flag)
198 {
199   int ierr,(*f)(PC,MatSORType);
200 
201   PetscFunctionBegin;
202   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
203   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr);
204   if (f) {
205     ierr = (*f)(pc,flag);CHKERRQ(ierr);
206   }
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "PCSORSetOmega"
212 /*@
213    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
214    (where omega = 1.0 by default).
215 
216    Collective on PC
217 
218    Input Parameters:
219 +  pc - the preconditioner context
220 -  omega - relaxation coefficient (0 < omega < 2).
221 
222    Options Database Key:
223 .  -pc_sor_omega <omega> - Sets omega
224 
225    Level: intermediate
226 
227 .keywords: PC, SOR, SSOR, set, relaxation, omega
228 
229 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
230 @*/
231 int PCSORSetOmega(PC pc,PetscReal omega)
232 {
233   int ierr,(*f)(PC,PetscReal);
234 
235   PetscFunctionBegin;
236   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr);
237   if (f) {
238     ierr = (*f)(pc,omega);CHKERRQ(ierr);
239   }
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "PCSORSetIterations"
245 /*@
246    PCSORSetIterations - Sets the number of inner iterations to
247    be used by the SOR preconditioner. The default is 1.
248 
249    Collective on PC
250 
251    Input Parameters:
252 +  pc - the preconditioner context
253 .  lits - number of local iterations, smoothings over just variables on processor
254 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
255 
256    Options Database Key:
257 +  -pc_sor_its <its> - Sets number of iterations
258 -  -pc_sor_lits <lits> - Sets number of local iterations
259 
260    Level: intermediate
261 
262    Notes: When run on one processor the number of smoothings is lits*its
263 
264 .keywords: PC, SOR, SSOR, set, iterations
265 
266 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
267 @*/
268 int PCSORSetIterations(PC pc,int its,int lits)
269 {
270   int ierr,(*f)(PC,int,int);
271 
272   PetscFunctionBegin;
273   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
274   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
275   if (f) {
276     ierr = (*f)(pc,its,lits);CHKERRQ(ierr);
277   }
278   PetscFunctionReturn(0);
279 }
280 
281 /*MC
282      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
283 
284    Options Database Keys:
285 +  -pc_sor_symmetric - Activates symmetric version
286 .  -pc_sor_backward - Activates backward version
287 .  -pc_sor_local_forward - Activates local forward version
288 .  -pc_sor_local_symmetric - Activates local symmetric version
289 .  -pc_sor_local_backward - Activates local backward version
290 .  -pc_sor_omega <omega> - Sets omega
291 .  -pc_sor_its <its> - Sets number of iterations
292 -  -pc_sor_lits <lits> - Sets number of local iterations
293 
294    Level: beginner
295 
296   Concepts: SOR, preconditioners, Gauss-Seidel
297 
298    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
299           Not a true parallel SOR, in parallel this implementation corresponds to block
300           Jacobi with SOR on each block.
301 
302           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
303 
304 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
305            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
306 M*/
307 
308 EXTERN_C_BEGIN
309 #undef __FUNCT__
310 #define __FUNCT__ "PCCreate_SOR"
311 int PCCreate_SOR(PC pc)
312 {
313   int    ierr;
314   PC_SOR *jac;
315 
316   PetscFunctionBegin;
317   ierr = PetscNew(PC_SOR,&jac);CHKERRQ(ierr);
318   PetscLogObjectMemory(pc,sizeof(PC_SOR));
319 
320   pc->ops->apply           = PCApply_SOR;
321   pc->ops->applyrichardson = PCApplyRichardson_SOR;
322   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
323   pc->ops->setup           = 0;
324   pc->ops->view            = PCView_SOR;
325   pc->ops->destroy         = PCDestroy_SOR;
326   pc->data           = (void*)jac;
327   jac->sym           = SOR_FORWARD_SWEEP;
328   jac->omega         = 1.0;
329   jac->its           = 1;
330   jac->lits          = 1;
331 
332   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
333                     PCSORSetSymmetric_SOR);CHKERRQ(ierr);
334   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
335                     PCSORSetOmega_SOR);CHKERRQ(ierr);
336   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
337                     PCSORSetIterations_SOR);CHKERRQ(ierr);
338 
339   PetscFunctionReturn(0);
340 }
341 EXTERN_C_END
342 
343 
344 
345 
346 
347