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