xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision a6dfd86ebdf75fa024070af26d51b62fd16650a3)
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)             sortype = "symmetric";
106     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
107     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
108     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
109     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
110     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
111     else                                                                     sortype = "unknown";
112     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);CHKERRQ(ierr);
113   }
114   PetscFunctionReturn(0);
115 }
116 
117 
118 /* ------------------------------------------------------------------------------*/
119 EXTERN_C_BEGIN
120 #undef __FUNCT__
121 #define __FUNCT__ "PCSORSetSymmetric_SOR"
122 PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
123 {
124   PC_SOR *jac;
125 
126   PetscFunctionBegin;
127   jac = (PC_SOR*)pc->data;
128   jac->sym = flag;
129   PetscFunctionReturn(0);
130 }
131 EXTERN_C_END
132 
133 EXTERN_C_BEGIN
134 #undef __FUNCT__
135 #define __FUNCT__ "PCSORSetOmega_SOR"
136 PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
137 {
138   PC_SOR *jac;
139 
140   PetscFunctionBegin;
141   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
142   jac        = (PC_SOR*)pc->data;
143   jac->omega = omega;
144   PetscFunctionReturn(0);
145 }
146 EXTERN_C_END
147 
148 EXTERN_C_BEGIN
149 #undef __FUNCT__
150 #define __FUNCT__ "PCSORSetIterations_SOR"
151 PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
152 {
153   PC_SOR *jac;
154 
155   PetscFunctionBegin;
156   jac      = (PC_SOR*)pc->data;
157   jac->its = its;
158   jac->lits = lits;
159   PetscFunctionReturn(0);
160 }
161 EXTERN_C_END
162 
163 /* ------------------------------------------------------------------------------*/
164 #undef __FUNCT__
165 #define __FUNCT__ "PCSORSetSymmetric"
166 /*@
167    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
168    backward, or forward relaxation.  The local variants perform SOR on
169    each processor.  By default forward relaxation is used.
170 
171    Logically Collective on PC
172 
173    Input Parameters:
174 +  pc - the preconditioner context
175 -  flag - one of the following
176 .vb
177     SOR_FORWARD_SWEEP
178     SOR_BACKWARD_SWEEP
179     SOR_SYMMETRIC_SWEEP
180     SOR_LOCAL_FORWARD_SWEEP
181     SOR_LOCAL_BACKWARD_SWEEP
182     SOR_LOCAL_SYMMETRIC_SWEEP
183 .ve
184 
185    Options Database Keys:
186 +  -pc_sor_symmetric - Activates symmetric version
187 .  -pc_sor_backward - Activates backward version
188 .  -pc_sor_local_forward - Activates local forward version
189 .  -pc_sor_local_symmetric - Activates local symmetric version
190 -  -pc_sor_local_backward - Activates local backward version
191 
192    Notes:
193    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
194    which can be chosen with the option
195 .  -pc_type eisenstat - Activates Eisenstat trick
196 
197    Level: intermediate
198 
199 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
200 
201 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
202 @*/
203 PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
204 {
205   PetscErrorCode ierr;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
209   PetscValidLogicalCollectiveEnum(pc,flag,2);
210   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "PCSORSetOmega"
216 /*@
217    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
218    (where omega = 1.0 by default).
219 
220    Logically Collective on PC
221 
222    Input Parameters:
223 +  pc - the preconditioner context
224 -  omega - relaxation coefficient (0 < omega < 2).
225 
226    Options Database Key:
227 .  -pc_sor_omega <omega> - Sets omega
228 
229    Level: intermediate
230 
231 .keywords: PC, SOR, SSOR, set, relaxation, omega
232 
233 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
234 @*/
235 PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
236 {
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
241   PetscValidLogicalCollectiveReal(pc,omega,2);
242   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "PCSORSetIterations"
248 /*@
249    PCSORSetIterations - Sets the number of inner iterations to
250    be used by the SOR preconditioner. The default is 1.
251 
252    Logically Collective on PC
253 
254    Input Parameters:
255 +  pc - the preconditioner context
256 .  lits - number of local iterations, smoothings over just variables on processor
257 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
258 
259    Options Database Key:
260 +  -pc_sor_its <its> - Sets number of iterations
261 -  -pc_sor_lits <lits> - Sets number of local iterations
262 
263    Level: intermediate
264 
265    Notes: When run on one processor the number of smoothings is lits*its
266 
267 .keywords: PC, SOR, SSOR, set, iterations
268 
269 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
270 @*/
271 PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
272 {
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
277   PetscValidLogicalCollectiveInt(pc,its,2);
278   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
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_forward - Activates forward version
289 .  -pc_sor_local_forward - Activates local forward version
290 .  -pc_sor_local_symmetric - Activates local symmetric version  (default 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   (default 1)
294 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
295 
296    Level: beginner
297 
298   Concepts: SOR, preconditioners, Gauss-Seidel
299 
300    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
301           Not a true parallel SOR, in parallel this implementation corresponds to block
302           Jacobi with SOR on each block.
303 
304           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
305 
306 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
307            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
308 M*/
309 
310 EXTERN_C_BEGIN
311 #undef __FUNCT__
312 #define __FUNCT__ "PCCreate_SOR"
313 PetscErrorCode  PCCreate_SOR(PC pc)
314 {
315   PetscErrorCode ierr;
316   PC_SOR         *jac;
317 
318   PetscFunctionBegin;
319   ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr);
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_LOCAL_SYMMETRIC_SWEEP;
329   jac->omega               = 1.0;
330   jac->fshift              = 0.0;
331   jac->its                 = 1;
332   jac->lits                = 1;
333 
334   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
335                     PCSORSetSymmetric_SOR);CHKERRQ(ierr);
336   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
337                     PCSORSetOmega_SOR);CHKERRQ(ierr);
338   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
339                     PCSORSetIterations_SOR);CHKERRQ(ierr);
340 
341   PetscFunctionReturn(0);
342 }
343 EXTERN_C_END
344 
345 
346 
347 
348 
349