xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 5e71baeff2f3138f93cd4f5927dfd596eb8325cc)
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 static PetscErrorCode PCDestroy_SOR(PC pc)
16 {
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   ierr = PetscFree(pc->data);CHKERRQ(ierr);
21   PetscFunctionReturn(0);
22 }
23 
24 static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
25 {
26   PC_SOR         *jac = (PC_SOR*)pc->data;
27   PetscErrorCode ierr;
28   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
29   PetscReal      fshift;
30 
31   PetscFunctionBegin;
32   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
33   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
34   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
39 {
40   PC_SOR         *jac = (PC_SOR*)pc->data;
41   PetscErrorCode ierr;
42   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
43   PetscReal      fshift;
44   PetscBool      set,sym;
45 
46   PetscFunctionBegin;
47   ierr = MatIsSymmetricKnown(pc->pmat,&set,&sym);CHKERRQ(ierr);
48   if (!set || !sym || (jac->sym != SOR_SYMMETRIC_SWEEP && jac->sym != SOR_LOCAL_SYMMETRIC_SWEEP)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Can only apply transpose of SOR if matrix is symmetric and sweep is symmetric");
49   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
50   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
51   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 
55 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)
56 {
57   PC_SOR         *jac = (PC_SOR*)pc->data;
58   PetscErrorCode ierr;
59   MatSORType     stype = jac->sym;
60   PetscReal      fshift;
61 
62   PetscFunctionBegin;
63   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
64   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
65   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
66   ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
67   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
68   *outits = its;
69   *reason = PCRICHARDSON_CONVERGED_ITS;
70   PetscFunctionReturn(0);
71 }
72 
73 PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
74 {
75   PC_SOR         *jac = (PC_SOR*)pc->data;
76   PetscErrorCode ierr;
77   PetscBool      flg;
78 
79   PetscFunctionBegin;
80   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
81   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
82   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
83   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
84   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
85   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
86   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
87   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
88   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
89   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
90   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
91   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
92   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
93   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
94   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
95   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
96   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
97   ierr = PetscOptionsTail();CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
102 {
103   PC_SOR         *jac = (PC_SOR*)pc->data;
104   MatSORType     sym  = jac->sym;
105   const char     *sortype;
106   PetscErrorCode ierr;
107   PetscBool      iascii;
108 
109   PetscFunctionBegin;
110   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
111   if (iascii) {
112     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  zero initial guess\n");CHKERRQ(ierr);}
113     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
114     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
115     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
116     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
117     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
118     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
119     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
120     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
121     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
122     else                                                                     sortype = "unknown";
123     ierr = PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
124   }
125   PetscFunctionReturn(0);
126 }
127 
128 
129 /* ------------------------------------------------------------------------------*/
130 static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
131 {
132   PC_SOR *jac = (PC_SOR*)pc->data;
133 
134   PetscFunctionBegin;
135   jac->sym = flag;
136   PetscFunctionReturn(0);
137 }
138 
139 static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
140 {
141   PC_SOR *jac = (PC_SOR*)pc->data;
142 
143   PetscFunctionBegin;
144   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
145   jac->omega = omega;
146   PetscFunctionReturn(0);
147 }
148 
149 static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
150 {
151   PC_SOR *jac = (PC_SOR*)pc->data;
152 
153   PetscFunctionBegin;
154   jac->its  = its;
155   jac->lits = lits;
156   PetscFunctionReturn(0);
157 }
158 
159 static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
160 {
161   PC_SOR *jac = (PC_SOR*)pc->data;
162 
163   PetscFunctionBegin;
164   *flag = jac->sym;
165   PetscFunctionReturn(0);
166 }
167 
168 static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
169 {
170   PC_SOR *jac = (PC_SOR*)pc->data;
171 
172   PetscFunctionBegin;
173   *omega = jac->omega;
174   PetscFunctionReturn(0);
175 }
176 
177 static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
178 {
179   PC_SOR *jac = (PC_SOR*)pc->data;
180 
181   PetscFunctionBegin;
182   if (its)  *its = jac->its;
183   if (lits) *lits = jac->lits;
184   PetscFunctionReturn(0);
185 }
186 
187 /* ------------------------------------------------------------------------------*/
188 /*@
189    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
190    each processor.  By default forward relaxation is used.
191 
192    Logically Collective on PC
193 
194    Input Parameter:
195 .  pc - the preconditioner context
196 
197    Output Parameter:
198 .  flag - one of the following
199 .vb
200     SOR_FORWARD_SWEEP
201     SOR_BACKWARD_SWEEP
202     SOR_SYMMETRIC_SWEEP
203     SOR_LOCAL_FORWARD_SWEEP
204     SOR_LOCAL_BACKWARD_SWEEP
205     SOR_LOCAL_SYMMETRIC_SWEEP
206 .ve
207 
208    Options Database Keys:
209 +  -pc_sor_symmetric - Activates symmetric version
210 .  -pc_sor_backward - Activates backward version
211 .  -pc_sor_local_forward - Activates local forward version
212 .  -pc_sor_local_symmetric - Activates local symmetric version
213 -  -pc_sor_local_backward - Activates local backward version
214 
215    Notes:
216    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
217    which can be chosen with the option
218 .  -pc_type eisenstat - Activates Eisenstat trick
219 
220    Level: intermediate
221 
222 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
223 
224 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
225 @*/
226 PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
227 {
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
232   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 
236 /*@
237    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
238    (where omega = 1.0 by default).
239 
240    Logically Collective on PC
241 
242    Input Parameter:
243 .  pc - the preconditioner context
244 
245    Output Parameter:
246 .  omega - relaxation coefficient (0 < omega < 2).
247 
248    Options Database Key:
249 .  -pc_sor_omega <omega> - Sets omega
250 
251    Level: intermediate
252 
253 .keywords: PC, SOR, SSOR, set, relaxation, omega
254 
255 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
256 @*/
257 PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
258 {
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
263   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
264   PetscFunctionReturn(0);
265 }
266 
267 /*@
268    PCSORGetIterations - Gets the number of inner iterations to
269    be used by the SOR preconditioner. The default is 1.
270 
271    Logically Collective on PC
272 
273    Input Parameter:
274 .  pc - the preconditioner context
275 
276    Output Parameter:
277 +  lits - number of local iterations, smoothings over just variables on processor
278 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
279 
280    Options Database Key:
281 +  -pc_sor_its <its> - Sets number of iterations
282 -  -pc_sor_lits <lits> - Sets number of local iterations
283 
284    Level: intermediate
285 
286    Notes: When run on one processor the number of smoothings is lits*its
287 
288 .keywords: PC, SOR, SSOR, set, iterations
289 
290 .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
291 @*/
292 PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
293 {
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
298   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
299   PetscFunctionReturn(0);
300 }
301 
302 /*@
303    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
304    backward, or forward relaxation.  The local variants perform SOR on
305    each processor.  By default forward relaxation is used.
306 
307    Logically Collective on PC
308 
309    Input Parameters:
310 +  pc - the preconditioner context
311 -  flag - one of the following
312 .vb
313     SOR_FORWARD_SWEEP
314     SOR_BACKWARD_SWEEP
315     SOR_SYMMETRIC_SWEEP
316     SOR_LOCAL_FORWARD_SWEEP
317     SOR_LOCAL_BACKWARD_SWEEP
318     SOR_LOCAL_SYMMETRIC_SWEEP
319 .ve
320 
321    Options Database Keys:
322 +  -pc_sor_symmetric - Activates symmetric version
323 .  -pc_sor_backward - Activates backward version
324 .  -pc_sor_local_forward - Activates local forward version
325 .  -pc_sor_local_symmetric - Activates local symmetric version
326 -  -pc_sor_local_backward - Activates local backward version
327 
328    Notes:
329    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
330    which can be chosen with the option
331 .  -pc_type eisenstat - Activates Eisenstat trick
332 
333    Level: intermediate
334 
335 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
336 
337 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
338 @*/
339 PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
340 {
341   PetscErrorCode ierr;
342 
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
345   PetscValidLogicalCollectiveEnum(pc,flag,2);
346   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
347   PetscFunctionReturn(0);
348 }
349 
350 /*@
351    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
352    (where omega = 1.0 by default).
353 
354    Logically Collective on PC
355 
356    Input Parameters:
357 +  pc - the preconditioner context
358 -  omega - relaxation coefficient (0 < omega < 2).
359 
360    Options Database Key:
361 .  -pc_sor_omega <omega> - Sets omega
362 
363    Level: intermediate
364 
365 .keywords: PC, SOR, SSOR, set, relaxation, omega
366 
367 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
368 @*/
369 PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
370 {
371   PetscErrorCode ierr;
372 
373   PetscFunctionBegin;
374   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
375   PetscValidLogicalCollectiveReal(pc,omega,2);
376   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 /*@
381    PCSORSetIterations - Sets the number of inner iterations to
382    be used by the SOR preconditioner. The default is 1.
383 
384    Logically Collective on PC
385 
386    Input Parameters:
387 +  pc - the preconditioner context
388 .  lits - number of local iterations, smoothings over just variables on processor
389 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
390 
391    Options Database Key:
392 +  -pc_sor_its <its> - Sets number of iterations
393 -  -pc_sor_lits <lits> - Sets number of local iterations
394 
395    Level: intermediate
396 
397    Notes: When run on one processor the number of smoothings is lits*its
398 
399 .keywords: PC, SOR, SSOR, set, iterations
400 
401 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
402 @*/
403 PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
404 {
405   PetscErrorCode ierr;
406 
407   PetscFunctionBegin;
408   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
409   PetscValidLogicalCollectiveInt(pc,its,2);
410   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 
414 /*MC
415      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
416 
417    Options Database Keys:
418 +  -pc_sor_symmetric - Activates symmetric version
419 .  -pc_sor_backward - Activates backward version
420 .  -pc_sor_forward - Activates forward version
421 .  -pc_sor_local_forward - Activates local forward version
422 .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
423 .  -pc_sor_local_backward - Activates local backward version
424 .  -pc_sor_omega <omega> - Sets omega
425 .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
426 .  -pc_sor_its <its> - Sets number of iterations   (default 1)
427 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
428 
429    Level: beginner
430 
431   Concepts: SOR, preconditioners, Gauss-Seidel
432 
433    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
434           Not a true parallel SOR, in parallel this implementation corresponds to block
435           Jacobi with SOR on each block.
436 
437           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
438           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
439           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
440           zero pivot.
441 
442           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
443 
444           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
445           the computation is stopped with an error
446 
447           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
448           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
449 
450 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
451            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
452 M*/
453 
454 PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
455 {
456   PetscErrorCode ierr;
457   PC_SOR         *jac;
458 
459   PetscFunctionBegin;
460   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
461 
462   pc->ops->apply           = PCApply_SOR;
463   pc->ops->applytranspose  = PCApplyTranspose_SOR;
464   pc->ops->applyrichardson = PCApplyRichardson_SOR;
465   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
466   pc->ops->setup           = 0;
467   pc->ops->view            = PCView_SOR;
468   pc->ops->destroy         = PCDestroy_SOR;
469   pc->data                 = (void*)jac;
470   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
471   jac->omega               = 1.0;
472   jac->fshift              = 0.0;
473   jac->its                 = 1;
474   jac->lits                = 1;
475 
476   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
477   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
478   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
479   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
480   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
481   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 
486 
487 
488 
489