xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 4d615ea068c82feafef4a984334a49e3f0dfdaf7) !
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,PetscInt *outits,PCRichardsonConvergedReason *reason)
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   ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
51   *outits = its;
52   *reason = PCRICHARDSON_CONVERGED_ITS;
53   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "PCSetFromOptions_SOR"
58 PetscErrorCode PCSetFromOptions_SOR(PC pc)
59 {
60   PC_SOR         *jac = (PC_SOR*)pc->data;
61   PetscErrorCode ierr;
62   PetscTruth     flg;
63 
64   PetscFunctionBegin;
65   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
66     ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr);
67     ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr);
68     ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
69     ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr);
70     ierr = PetscOptionsTruthGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
71     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
72     ierr = PetscOptionsTruthGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
73     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
74     ierr = PetscOptionsTruthGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
75     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
76     ierr = PetscOptionsTruthGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
77     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
78     ierr = PetscOptionsTruthGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
79     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
80     ierr = PetscOptionsTruthGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
81     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
82   ierr = PetscOptionsTail();CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "PCView_SOR"
88 PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
89 {
90   PC_SOR         *jac = (PC_SOR*)pc->data;
91   MatSORType     sym = jac->sym;
92   const char     *sortype;
93   PetscErrorCode ierr;
94   PetscTruth     iascii;
95 
96   PetscFunctionBegin;
97   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
98   if (iascii) {
99     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
100     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
101     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
102     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
103     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
104                                              sortype = "symmetric";
105     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
106     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
107     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
108                                              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   } else {
114     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 
120 /* ------------------------------------------------------------------------------*/
121 EXTERN_C_BEGIN
122 #undef __FUNCT__
123 #define __FUNCT__ "PCSORSetSymmetric_SOR"
124 PetscErrorCode PETSCKSP_DLLEXPORT 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 PETSCKSP_DLLEXPORT PCSORSetOmega_SOR(PC pc,PetscReal omega)
139 {
140   PC_SOR *jac;
141 
142   PetscFunctionBegin;
143   if (omega >= 2.0 || omega <= 0.0) SETERRQ(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 PETSCKSP_DLLEXPORT 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    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 PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC pc,MatSORType flag)
206 {
207   PetscErrorCode ierr,(*f)(PC,MatSORType);
208 
209   PetscFunctionBegin;
210   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
211   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr);
212   if (f) {
213     ierr = (*f)(pc,flag);CHKERRQ(ierr);
214   }
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    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 PETSCKSP_DLLEXPORT PCSORSetOmega(PC pc,PetscReal omega)
240 {
241   PetscErrorCode ierr,(*f)(PC,PetscReal);
242 
243   PetscFunctionBegin;
244   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr);
245   if (f) {
246     ierr = (*f)(pc,omega);CHKERRQ(ierr);
247   }
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "PCSORSetIterations"
253 /*@
254    PCSORSetIterations - Sets the number of inner iterations to
255    be used by the SOR preconditioner. The default is 1.
256 
257    Collective on PC
258 
259    Input Parameters:
260 +  pc - the preconditioner context
261 .  lits - number of local iterations, smoothings over just variables on processor
262 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
263 
264    Options Database Key:
265 +  -pc_sor_its <its> - Sets number of iterations
266 -  -pc_sor_lits <lits> - Sets number of local iterations
267 
268    Level: intermediate
269 
270    Notes: When run on one processor the number of smoothings is lits*its
271 
272 .keywords: PC, SOR, SSOR, set, iterations
273 
274 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
275 @*/
276 PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
277 {
278   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt);
279 
280   PetscFunctionBegin;
281   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
282   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
283   if (f) {
284     ierr = (*f)(pc,its,lits);CHKERRQ(ierr);
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 /*MC
290      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
291 
292    Options Database Keys:
293 +  -pc_sor_symmetric - Activates symmetric version
294 .  -pc_sor_backward - Activates backward version
295 .  -pc_sor_forward - Activates forward version
296 .  -pc_sor_local_forward - Activates local forward version
297 .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
298 .  -pc_sor_local_backward - Activates local backward version
299 .  -pc_sor_omega <omega> - Sets omega
300 .  -pc_sor_its <its> - Sets number of iterations   (default 1)
301 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
302 
303    Level: beginner
304 
305   Concepts: SOR, preconditioners, Gauss-Seidel
306 
307    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
308           Not a true parallel SOR, in parallel this implementation corresponds to block
309           Jacobi with SOR on each block.
310 
311           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
312 
313 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
314            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
315 M*/
316 
317 EXTERN_C_BEGIN
318 #undef __FUNCT__
319 #define __FUNCT__ "PCCreate_SOR"
320 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SOR(PC pc)
321 {
322   PetscErrorCode ierr;
323   PC_SOR         *jac;
324 
325   PetscFunctionBegin;
326   ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr);
327 
328   pc->ops->apply           = PCApply_SOR;
329   pc->ops->applyrichardson = PCApplyRichardson_SOR;
330   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
331   pc->ops->setup           = 0;
332   pc->ops->view            = PCView_SOR;
333   pc->ops->destroy         = PCDestroy_SOR;
334   pc->data                 = (void*)jac;
335   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
336   jac->omega               = 1.0;
337   jac->fshift              = 0.0;
338   jac->its                 = 1;
339   jac->lits                = 1;
340 
341   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
342                     PCSORSetSymmetric_SOR);CHKERRQ(ierr);
343   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
344                     PCSORSetOmega_SOR);CHKERRQ(ierr);
345   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
346                     PCSORSetIterations_SOR);CHKERRQ(ierr);
347 
348   PetscFunctionReturn(0);
349 }
350 EXTERN_C_END
351 
352 
353 
354 
355 
356