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