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