xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 285fb4e2b69b3de46a0633bd0adc6a7f684caa1e)
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:
287     When run on one processor the number of smoothings is lits*its
288 
289 .keywords: PC, SOR, SSOR, set, iterations
290 
291 .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
292 @*/
293 PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
294 {
295   PetscErrorCode ierr;
296 
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
299   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 /*@
304    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
305    backward, or forward relaxation.  The local variants perform SOR on
306    each processor.  By default forward relaxation is used.
307 
308    Logically Collective on PC
309 
310    Input Parameters:
311 +  pc - the preconditioner context
312 -  flag - one of the following
313 .vb
314     SOR_FORWARD_SWEEP
315     SOR_BACKWARD_SWEEP
316     SOR_SYMMETRIC_SWEEP
317     SOR_LOCAL_FORWARD_SWEEP
318     SOR_LOCAL_BACKWARD_SWEEP
319     SOR_LOCAL_SYMMETRIC_SWEEP
320 .ve
321 
322    Options Database Keys:
323 +  -pc_sor_symmetric - Activates symmetric version
324 .  -pc_sor_backward - Activates backward version
325 .  -pc_sor_local_forward - Activates local forward version
326 .  -pc_sor_local_symmetric - Activates local symmetric version
327 -  -pc_sor_local_backward - Activates local backward version
328 
329    Notes:
330    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
331    which can be chosen with the option
332 .  -pc_type eisenstat - Activates Eisenstat trick
333 
334    Level: intermediate
335 
336 .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
337 
338 .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
339 @*/
340 PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
341 {
342   PetscErrorCode ierr;
343 
344   PetscFunctionBegin;
345   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
346   PetscValidLogicalCollectiveEnum(pc,flag,2);
347   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 /*@
352    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
353    (where omega = 1.0 by default).
354 
355    Logically Collective on PC
356 
357    Input Parameters:
358 +  pc - the preconditioner context
359 -  omega - relaxation coefficient (0 < omega < 2).
360 
361    Options Database Key:
362 .  -pc_sor_omega <omega> - Sets omega
363 
364    Level: intermediate
365 
366 .keywords: PC, SOR, SSOR, set, relaxation, omega
367 
368 .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
369 @*/
370 PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
371 {
372   PetscErrorCode ierr;
373 
374   PetscFunctionBegin;
375   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
376   PetscValidLogicalCollectiveReal(pc,omega,2);
377   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 /*@
382    PCSORSetIterations - Sets the number of inner iterations to
383    be used by the SOR preconditioner. The default is 1.
384 
385    Logically Collective on PC
386 
387    Input Parameters:
388 +  pc - the preconditioner context
389 .  lits - number of local iterations, smoothings over just variables on processor
390 -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
391 
392    Options Database Key:
393 +  -pc_sor_its <its> - Sets number of iterations
394 -  -pc_sor_lits <lits> - Sets number of local iterations
395 
396    Level: intermediate
397 
398    Notes:
399     When run on one processor the number of smoothings is lits*its
400 
401 .keywords: PC, SOR, SSOR, set, iterations
402 
403 .seealso: PCSORSetOmega(), PCSORSetSymmetric()
404 @*/
405 PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
406 {
407   PetscErrorCode ierr;
408 
409   PetscFunctionBegin;
410   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
411   PetscValidLogicalCollectiveInt(pc,its,2);
412   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
413   PetscFunctionReturn(0);
414 }
415 
416 /*MC
417      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
418 
419    Options Database Keys:
420 +  -pc_sor_symmetric - Activates symmetric version
421 .  -pc_sor_backward - Activates backward version
422 .  -pc_sor_forward - Activates forward version
423 .  -pc_sor_local_forward - Activates local forward version
424 .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
425 .  -pc_sor_local_backward - Activates local backward version
426 .  -pc_sor_omega <omega> - Sets omega
427 .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
428 .  -pc_sor_its <its> - Sets number of iterations   (default 1)
429 -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
430 
431    Level: beginner
432 
433   Concepts: SOR, preconditioners, Gauss-Seidel
434 
435    Notes:
436     Only implemented for the AIJ  and SeqBAIJ matrix formats.
437           Not a true parallel SOR, in parallel this implementation corresponds to block
438           Jacobi with SOR on each block.
439 
440           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
441           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
442           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
443           zero pivot.
444 
445           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
446 
447           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
448           the computation is stopped with an error
449 
450           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
451           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
452 
453 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
454            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
455 M*/
456 
457 PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
458 {
459   PetscErrorCode ierr;
460   PC_SOR         *jac;
461 
462   PetscFunctionBegin;
463   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
464 
465   pc->ops->apply           = PCApply_SOR;
466   pc->ops->applytranspose  = PCApplyTranspose_SOR;
467   pc->ops->applyrichardson = PCApplyRichardson_SOR;
468   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
469   pc->ops->setup           = 0;
470   pc->ops->view            = PCView_SOR;
471   pc->ops->destroy         = PCDestroy_SOR;
472   pc->data                 = (void*)jac;
473   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
474   jac->omega               = 1.0;
475   jac->fshift              = 0.0;
476   jac->its                 = 1;
477   jac->lits                = 1;
478 
479   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
480   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
481   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
482   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
483   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
484   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
485   PetscFunctionReturn(0);
486 }
487 
488 
489 
490 
491 
492