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