xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 047240e14af00aad1ef65e96f6fface8924f7f7e)
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,0);CHKERRQ(ierr);
67   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr);
68   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
69   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);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,jac->omega);CHKERRQ(ierr);
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 
116 /* ------------------------------------------------------------------------------*/
117 EXTERN_C_BEGIN
118 #undef __FUNCT__
119 #define __FUNCT__ "PCSORSetSymmetric_SOR"
120 PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
121 {
122   PC_SOR *jac;
123 
124   PetscFunctionBegin;
125   jac      = (PC_SOR*)pc->data;
126   jac->sym = flag;
127   PetscFunctionReturn(0);
128 }
129 EXTERN_C_END
130 
131 EXTERN_C_BEGIN
132 #undef __FUNCT__
133 #define __FUNCT__ "PCSORSetOmega_SOR"
134 PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
135 {
136   PC_SOR *jac;
137 
138   PetscFunctionBegin;
139   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
140   jac        = (PC_SOR*)pc->data;
141   jac->omega = omega;
142   PetscFunctionReturn(0);
143 }
144 EXTERN_C_END
145 
146 EXTERN_C_BEGIN
147 #undef __FUNCT__
148 #define __FUNCT__ "PCSORSetIterations_SOR"
149 PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
150 {
151   PC_SOR *jac;
152 
153   PetscFunctionBegin;
154   jac       = (PC_SOR*)pc->data;
155   jac->its  = its;
156   jac->lits = lits;
157   PetscFunctionReturn(0);
158 }
159 EXTERN_C_END
160 
161 /* ------------------------------------------------------------------------------*/
162 #undef __FUNCT__
163 #define __FUNCT__ "PCSORSetSymmetric"
164 /*@
165    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
166    backward, or forward relaxation.  The local variants perform SOR on
167    each processor.  By default forward relaxation is used.
168 
169    Logically Collective on PC
170 
171    Input Parameters:
172 +  pc - the preconditioner context
173 -  flag - one of the following
174 .vb
175     SOR_FORWARD_SWEEP
176     SOR_BACKWARD_SWEEP
177     SOR_SYMMETRIC_SWEEP
178     SOR_LOCAL_FORWARD_SWEEP
179     SOR_LOCAL_BACKWARD_SWEEP
180     SOR_LOCAL_SYMMETRIC_SWEEP
181 .ve
182 
183    Options Database Keys:
184 +  -pc_sor_symmetric - Activates symmetric version
185 .  -pc_sor_backward - Activates backward version
186 .  -pc_sor_local_forward - Activates local forward version
187 .  -pc_sor_local_symmetric - Activates local symmetric version
188 -  -pc_sor_local_backward - Activates local backward version
189 
190    Notes:
191    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
192    which can be chosen with the option
193 .  -pc_type eisenstat - Activates Eisenstat trick
194 
195    Level: intermediate
196 
197 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
198 
199 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
200 @*/
201 PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
202 {
203   PetscErrorCode ierr;
204 
205   PetscFunctionBegin;
206   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
207   PetscValidLogicalCollectiveEnum(pc,flag,2);
208   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
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    Logically 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 PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
234 {
235   PetscErrorCode ierr;
236 
237   PetscFunctionBegin;
238   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
239   PetscValidLogicalCollectiveReal(pc,omega,2);
240   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
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    Logically 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,PetscInt its,PetscInt lits)
270 {
271   PetscErrorCode ierr;
272 
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
275   PetscValidLogicalCollectiveInt(pc,its,2);
276   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 /*MC
281      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
282 
283    Options Database Keys:
284 +  -pc_sor_symmetric - Activates symmetric version
285 .  -pc_sor_backward - Activates backward version
286 .  -pc_sor_forward - Activates forward version
287 .  -pc_sor_local_forward - Activates local forward version
288 .  -pc_sor_local_symmetric - Activates local symmetric version  (default 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   (default 1)
292 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
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 PetscErrorCode  PCCreate_SOR(PC pc)
312 {
313   PetscErrorCode ierr;
314   PC_SOR         *jac;
315 
316   PetscFunctionBegin;
317   ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr);
318 
319   pc->ops->apply           = PCApply_SOR;
320   pc->ops->applyrichardson = PCApplyRichardson_SOR;
321   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
322   pc->ops->setup           = 0;
323   pc->ops->view            = PCView_SOR;
324   pc->ops->destroy         = PCDestroy_SOR;
325   pc->data                 = (void*)jac;
326   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
327   jac->omega               = 1.0;
328   jac->fshift              = 0.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   PetscFunctionReturn(0);
339 }
340 EXTERN_C_END
341 
342 
343 
344 
345 
346