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