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