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