xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision ae15b995b5732fffd2de5a75cf61ef7190c6fef1)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith    Defines a  (S)SOR  preconditioner for any Mat implementation
54b9ad928SBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"               /*I "petscpc.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith typedef struct {
9c1ac3661SBarry Smith   PetscInt   its;        /* inner iterations, number of sweeps */
10c1ac3661SBarry Smith   PetscInt   lits;       /* local inner iterations, number of sweeps applied by the local matrix mat->A */
114b9ad928SBarry Smith   MatSORType sym;        /* forward, reverse, symmetric etc. */
124b9ad928SBarry Smith   PetscReal  omega;
134b9ad928SBarry Smith } PC_SOR;
144b9ad928SBarry Smith 
154b9ad928SBarry Smith #undef __FUNCT__
164b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR"
176849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc)
184b9ad928SBarry Smith {
194b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
20dfbe8321SBarry Smith   PetscErrorCode ierr;
214b9ad928SBarry Smith 
224b9ad928SBarry Smith   PetscFunctionBegin;
234b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
244b9ad928SBarry Smith   PetscFunctionReturn(0);
254b9ad928SBarry Smith }
264b9ad928SBarry Smith 
274b9ad928SBarry Smith #undef __FUNCT__
284b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR"
296849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
304b9ad928SBarry Smith {
314b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
32dfbe8321SBarry Smith   PetscErrorCode ierr;
33c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
344b9ad928SBarry Smith 
354b9ad928SBarry Smith   PetscFunctionBegin;
364b9ad928SBarry Smith   ierr = MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,0.0,jac->its,jac->lits,y);CHKERRQ(ierr);
374b9ad928SBarry Smith   PetscFunctionReturn(0);
384b9ad928SBarry Smith }
394b9ad928SBarry Smith 
404b9ad928SBarry Smith #undef __FUNCT__
414b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR"
4270441072SBarry Smith static PetscErrorCode PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
434b9ad928SBarry Smith {
444b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
45dfbe8321SBarry Smith   PetscErrorCode ierr;
464b9ad928SBarry Smith 
474b9ad928SBarry Smith   PetscFunctionBegin;
48*ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
494b9ad928SBarry Smith   its  = its*jac->its;
504b9ad928SBarry Smith   ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,0.0,its,jac->lits,y);CHKERRQ(ierr);
514b9ad928SBarry Smith   PetscFunctionReturn(0);
524b9ad928SBarry Smith }
534b9ad928SBarry Smith 
544b9ad928SBarry Smith #undef __FUNCT__
554b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR"
56dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc)
574b9ad928SBarry Smith {
584b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
59dfbe8321SBarry Smith   PetscErrorCode ierr;
604b9ad928SBarry Smith   PetscTruth     flg;
614b9ad928SBarry Smith 
624b9ad928SBarry Smith   PetscFunctionBegin;
634b9ad928SBarry Smith   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
644b9ad928SBarry Smith     ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr);
654b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
664b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr);
674dc4c822SBarry Smith     ierr = PetscOptionsTruthGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
684b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
694dc4c822SBarry Smith     ierr = PetscOptionsTruthGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
704b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
71a9510f2eSBarry Smith     ierr = PetscOptionsTruthGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
72a9510f2eSBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
73a8c7a070SBarry Smith     ierr = PetscOptionsTruthGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
744b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
754dc4c822SBarry Smith     ierr = PetscOptionsTruthGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
764b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
774dc4c822SBarry Smith     ierr = PetscOptionsTruthGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
784b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
794b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
804b9ad928SBarry Smith   PetscFunctionReturn(0);
814b9ad928SBarry Smith }
824b9ad928SBarry Smith 
834b9ad928SBarry Smith #undef __FUNCT__
844b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR"
85dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
864b9ad928SBarry Smith {
874b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
884b9ad928SBarry Smith   MatSORType     sym = jac->sym;
892fc52814SBarry Smith   const char     *sortype;
90dfbe8321SBarry Smith   PetscErrorCode ierr;
9132077d6dSBarry Smith   PetscTruth     iascii;
924b9ad928SBarry Smith 
934b9ad928SBarry Smith   PetscFunctionBegin;
9432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
9532077d6dSBarry Smith   if (iascii) {
964b9ad928SBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
974b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
984b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
994b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
1004b9ad928SBarry Smith     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
1014b9ad928SBarry Smith                                              sortype = "symmetric";
1024b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
1034b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
1044b9ad928SBarry Smith     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
1054b9ad928SBarry Smith                                              sortype = "local_symmetric";
1064b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)  sortype = "local_forward";
1074b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward";
1084b9ad928SBarry Smith     else                                     sortype = "unknown";
109a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, omega = %G\n",sortype,jac->its,jac->omega);CHKERRQ(ierr);
1104b9ad928SBarry Smith   } else {
11179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
1124b9ad928SBarry Smith   }
1134b9ad928SBarry Smith   PetscFunctionReturn(0);
1144b9ad928SBarry Smith }
1154b9ad928SBarry Smith 
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1184b9ad928SBarry Smith EXTERN_C_BEGIN
1194b9ad928SBarry Smith #undef __FUNCT__
1204b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR"
121dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1224b9ad928SBarry Smith {
1234b9ad928SBarry Smith   PC_SOR *jac;
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith   PetscFunctionBegin;
1264b9ad928SBarry Smith   jac = (PC_SOR*)pc->data;
1274b9ad928SBarry Smith   jac->sym = flag;
1284b9ad928SBarry Smith   PetscFunctionReturn(0);
1294b9ad928SBarry Smith }
1304b9ad928SBarry Smith EXTERN_C_END
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith EXTERN_C_BEGIN
1334b9ad928SBarry Smith #undef __FUNCT__
1344b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR"
135dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega_SOR(PC pc,PetscReal omega)
1364b9ad928SBarry Smith {
1374b9ad928SBarry Smith   PC_SOR *jac;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
1404b9ad928SBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1414b9ad928SBarry Smith   jac        = (PC_SOR*)pc->data;
1424b9ad928SBarry Smith   jac->omega = omega;
1434b9ad928SBarry Smith   PetscFunctionReturn(0);
1444b9ad928SBarry Smith }
1454b9ad928SBarry Smith EXTERN_C_END
1464b9ad928SBarry Smith 
1474b9ad928SBarry Smith EXTERN_C_BEGIN
1484b9ad928SBarry Smith #undef __FUNCT__
1494b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR"
150dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1514b9ad928SBarry Smith {
1524b9ad928SBarry Smith   PC_SOR *jac;
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith   PetscFunctionBegin;
1554b9ad928SBarry Smith   jac      = (PC_SOR*)pc->data;
1564b9ad928SBarry Smith   jac->its = its;
1574b9ad928SBarry Smith   jac->lits = lits;
1584b9ad928SBarry Smith   PetscFunctionReturn(0);
1594b9ad928SBarry Smith }
1604b9ad928SBarry Smith EXTERN_C_END
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1634b9ad928SBarry Smith #undef __FUNCT__
1644b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric"
1654b9ad928SBarry Smith /*@
1664b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
1674b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
1684b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
1694b9ad928SBarry Smith 
1704b9ad928SBarry Smith    Collective on PC
1714b9ad928SBarry Smith 
1724b9ad928SBarry Smith    Input Parameters:
1734b9ad928SBarry Smith +  pc - the preconditioner context
1744b9ad928SBarry Smith -  flag - one of the following
1754b9ad928SBarry Smith .vb
1764b9ad928SBarry Smith     SOR_FORWARD_SWEEP
1774b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
1784b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
1794b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
1804b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
1814b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
1824b9ad928SBarry Smith .ve
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith    Options Database Keys:
1854b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
1864b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
1874b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
1884b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
1894b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
1904b9ad928SBarry Smith 
1914b9ad928SBarry Smith    Notes:
1924b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
1934b9ad928SBarry Smith    which can be chosen with the option
1944b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
1954b9ad928SBarry Smith 
1964b9ad928SBarry Smith    Level: intermediate
1974b9ad928SBarry Smith 
1984b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
2014b9ad928SBarry Smith @*/
202dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC pc,MatSORType flag)
2034b9ad928SBarry Smith {
204dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,MatSORType);
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith   PetscFunctionBegin;
2074482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
2084b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr);
2094b9ad928SBarry Smith   if (f) {
2104b9ad928SBarry Smith     ierr = (*f)(pc,flag);CHKERRQ(ierr);
2114b9ad928SBarry Smith   }
2124b9ad928SBarry Smith   PetscFunctionReturn(0);
2134b9ad928SBarry Smith }
2144b9ad928SBarry Smith 
2154b9ad928SBarry Smith #undef __FUNCT__
2164b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega"
2174b9ad928SBarry Smith /*@
2184b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
2194b9ad928SBarry Smith    (where omega = 1.0 by default).
2204b9ad928SBarry Smith 
2214b9ad928SBarry Smith    Collective on PC
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith    Input Parameters:
2244b9ad928SBarry Smith +  pc - the preconditioner context
2254b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith    Options Database Key:
2284b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
2294b9ad928SBarry Smith 
2304b9ad928SBarry Smith    Level: intermediate
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
2334b9ad928SBarry Smith 
2344b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
2354b9ad928SBarry Smith @*/
236dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC pc,PetscReal omega)
2374b9ad928SBarry Smith {
238dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   PetscFunctionBegin;
2414b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr);
2424b9ad928SBarry Smith   if (f) {
2434b9ad928SBarry Smith     ierr = (*f)(pc,omega);CHKERRQ(ierr);
2444b9ad928SBarry Smith   }
2454b9ad928SBarry Smith   PetscFunctionReturn(0);
2464b9ad928SBarry Smith }
2474b9ad928SBarry Smith 
2484b9ad928SBarry Smith #undef __FUNCT__
2494b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations"
2504b9ad928SBarry Smith /*@
2514b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
2524b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
2534b9ad928SBarry Smith 
2544b9ad928SBarry Smith    Collective on PC
2554b9ad928SBarry Smith 
2564b9ad928SBarry Smith    Input Parameters:
2574b9ad928SBarry Smith +  pc - the preconditioner context
2584b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
2594b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
2604b9ad928SBarry Smith 
2614b9ad928SBarry Smith    Options Database Key:
2624b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
2634b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Level: intermediate
2664b9ad928SBarry Smith 
2674b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
2724b9ad928SBarry Smith @*/
273dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
2744b9ad928SBarry Smith {
275c1ac3661SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt);
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith   PetscFunctionBegin;
2784482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
2794b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
2804b9ad928SBarry Smith   if (f) {
2814b9ad928SBarry Smith     ierr = (*f)(pc,its,lits);CHKERRQ(ierr);
2824b9ad928SBarry Smith   }
2834b9ad928SBarry Smith   PetscFunctionReturn(0);
2844b9ad928SBarry Smith }
2854b9ad928SBarry Smith 
2864b9ad928SBarry Smith /*MC
2874b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
2884b9ad928SBarry Smith 
2894b9ad928SBarry Smith    Options Database Keys:
2904b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
2914b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
292a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
2934b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
294a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
2954b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
2964b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
297a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
298a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith    Level: beginner
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
3034b9ad928SBarry Smith 
30437a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
3054b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
3064b9ad928SBarry Smith           Jacobi with SOR on each block.
3074b9ad928SBarry Smith 
30837a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
30937a17b4dSBarry Smith 
3104b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
3114b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
3124b9ad928SBarry Smith M*/
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith EXTERN_C_BEGIN
3154b9ad928SBarry Smith #undef __FUNCT__
3164b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR"
317dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SOR(PC pc)
3184b9ad928SBarry Smith {
319dfbe8321SBarry Smith   PetscErrorCode ierr;
3204b9ad928SBarry Smith   PC_SOR         *jac;
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith   PetscFunctionBegin;
3234b9ad928SBarry Smith   ierr = PetscNew(PC_SOR,&jac);CHKERRQ(ierr);
32452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_SOR));CHKERRQ(ierr);
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
3274b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
3284b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
3294b9ad928SBarry Smith   pc->ops->setup           = 0;
3304b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
3314b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
3324b9ad928SBarry Smith   pc->data                 = (void*)jac;
333d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
3344b9ad928SBarry Smith   jac->omega               = 1.0;
3354b9ad928SBarry Smith   jac->its                 = 1;
3364b9ad928SBarry Smith   jac->lits                = 1;
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
3394b9ad928SBarry Smith                     PCSORSetSymmetric_SOR);CHKERRQ(ierr);
3404b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
3414b9ad928SBarry Smith                     PCSORSetOmega_SOR);CHKERRQ(ierr);
3424b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
3434b9ad928SBarry Smith                     PCSORSetIterations_SOR);CHKERRQ(ierr);
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith   PetscFunctionReturn(0);
3464b9ad928SBarry Smith }
3474b9ad928SBarry Smith EXTERN_C_END
3484b9ad928SBarry Smith 
3494b9ad928SBarry Smith 
3504b9ad928SBarry Smith 
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith 
353