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