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