xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Defines a  (S)SOR  preconditioner for any Mat implementation
44b9ad928SBarry Smith */
5*b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>               /*I "petscpc.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith typedef struct {
8c1ac3661SBarry Smith   PetscInt    its;        /* inner iterations, number of sweeps */
9c1ac3661SBarry Smith   PetscInt    lits;       /* local inner iterations, number of sweeps applied by the local matrix mat->A */
104b9ad928SBarry Smith   MatSORType  sym;        /* forward, reverse, symmetric etc. */
114b9ad928SBarry Smith   PetscReal   omega;
1229c1d7e0SHong Zhang   PetscReal   fshift;
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 {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
22c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
234b9ad928SBarry Smith   PetscFunctionReturn(0);
244b9ad928SBarry Smith }
254b9ad928SBarry Smith 
264b9ad928SBarry Smith #undef __FUNCT__
274b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR"
286849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
294b9ad928SBarry Smith {
304b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
31dfbe8321SBarry Smith   PetscErrorCode ierr;
32c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
334b9ad928SBarry Smith 
344b9ad928SBarry Smith   PetscFunctionBegin;
3541f059aeSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
364b9ad928SBarry Smith   PetscFunctionReturn(0);
374b9ad928SBarry Smith }
384b9ad928SBarry Smith 
394b9ad928SBarry Smith #undef __FUNCT__
404b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR"
41ace3abfcSBarry Smith 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)
424b9ad928SBarry Smith {
434b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
44dfbe8321SBarry Smith   PetscErrorCode ierr;
457319c654SBarry Smith   MatSORType     stype = jac->sym;
464b9ad928SBarry Smith 
474b9ad928SBarry Smith   PetscFunctionBegin;
48ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
497319c654SBarry Smith   if (guesszero) {
50de8939e3SBarry Smith     stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
517319c654SBarry Smith   }
5241f059aeSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
534d0a8057SBarry Smith   *outits = its;
544d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
554b9ad928SBarry Smith   PetscFunctionReturn(0);
564b9ad928SBarry Smith }
574b9ad928SBarry Smith 
584b9ad928SBarry Smith #undef __FUNCT__
594b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR"
60dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc)
614b9ad928SBarry Smith {
624b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
63dfbe8321SBarry Smith   PetscErrorCode ierr;
64ace3abfcSBarry Smith   PetscBool      flg;
654b9ad928SBarry Smith 
664b9ad928SBarry Smith   PetscFunctionBegin;
674b9ad928SBarry Smith   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
684b9ad928SBarry Smith     ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr);
6996fc60bcSBarry Smith     ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr);
704b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
714b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr);
72acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
734b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
74acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
754b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
76acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
77a9510f2eSBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
78acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
794b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
80acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
814b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
82acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
834b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
844b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
854b9ad928SBarry Smith   PetscFunctionReturn(0);
864b9ad928SBarry Smith }
874b9ad928SBarry Smith 
884b9ad928SBarry Smith #undef __FUNCT__
894b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR"
90dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
914b9ad928SBarry Smith {
924b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
934b9ad928SBarry Smith   MatSORType     sym = jac->sym;
942fc52814SBarry Smith   const char     *sortype;
95dfbe8321SBarry Smith   PetscErrorCode ierr;
96ace3abfcSBarry Smith   PetscBool      iascii;
974b9ad928SBarry Smith 
984b9ad928SBarry Smith   PetscFunctionBegin;
992692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10032077d6dSBarry Smith   if (iascii) {
1014b9ad928SBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
1024b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
1034b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
1044b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
1054b9ad928SBarry Smith     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
1064b9ad928SBarry Smith                                              sortype = "symmetric";
1074b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
1084b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
1094b9ad928SBarry Smith     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
1104b9ad928SBarry Smith                                              sortype = "local_symmetric";
1114b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)  sortype = "local_forward";
1124b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward";
1134b9ad928SBarry Smith     else                                     sortype = "unknown";
114cf576665SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);CHKERRQ(ierr);
1154b9ad928SBarry Smith   } else {
11665e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
1174b9ad928SBarry Smith   }
1184b9ad928SBarry Smith   PetscFunctionReturn(0);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1234b9ad928SBarry Smith EXTERN_C_BEGIN
1244b9ad928SBarry Smith #undef __FUNCT__
1254b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR"
1267087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1274b9ad928SBarry Smith {
1284b9ad928SBarry Smith   PC_SOR *jac;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
1314b9ad928SBarry Smith   jac = (PC_SOR*)pc->data;
1324b9ad928SBarry Smith   jac->sym = flag;
1334b9ad928SBarry Smith   PetscFunctionReturn(0);
1344b9ad928SBarry Smith }
1354b9ad928SBarry Smith EXTERN_C_END
1364b9ad928SBarry Smith 
1374b9ad928SBarry Smith EXTERN_C_BEGIN
1384b9ad928SBarry Smith #undef __FUNCT__
1394b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR"
1407087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1414b9ad928SBarry Smith {
1424b9ad928SBarry Smith   PC_SOR *jac;
1434b9ad928SBarry Smith 
1444b9ad928SBarry Smith   PetscFunctionBegin;
145e7e72b3dSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1464b9ad928SBarry Smith   jac        = (PC_SOR*)pc->data;
1474b9ad928SBarry Smith   jac->omega = omega;
1484b9ad928SBarry Smith   PetscFunctionReturn(0);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith EXTERN_C_END
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith EXTERN_C_BEGIN
1534b9ad928SBarry Smith #undef __FUNCT__
1544b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR"
1557087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1564b9ad928SBarry Smith {
1574b9ad928SBarry Smith   PC_SOR *jac;
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith   PetscFunctionBegin;
1604b9ad928SBarry Smith   jac      = (PC_SOR*)pc->data;
1614b9ad928SBarry Smith   jac->its = its;
1624b9ad928SBarry Smith   jac->lits = lits;
1634b9ad928SBarry Smith   PetscFunctionReturn(0);
1644b9ad928SBarry Smith }
1654b9ad928SBarry Smith EXTERN_C_END
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1684b9ad928SBarry Smith #undef __FUNCT__
1694b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric"
1704b9ad928SBarry Smith /*@
1714b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
1724b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
1734b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
1744b9ad928SBarry Smith 
1753f9fe445SBarry Smith    Logically Collective on PC
1764b9ad928SBarry Smith 
1774b9ad928SBarry Smith    Input Parameters:
1784b9ad928SBarry Smith +  pc - the preconditioner context
1794b9ad928SBarry Smith -  flag - one of the following
1804b9ad928SBarry Smith .vb
1814b9ad928SBarry Smith     SOR_FORWARD_SWEEP
1824b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
1834b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
1844b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
1854b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
1864b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
1874b9ad928SBarry Smith .ve
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith    Options Database Keys:
1904b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
1914b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
1924b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
1934b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
1944b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
1954b9ad928SBarry Smith 
1964b9ad928SBarry Smith    Notes:
1974b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
1984b9ad928SBarry Smith    which can be chosen with the option
1994b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
2004b9ad928SBarry Smith 
2014b9ad928SBarry Smith    Level: intermediate
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
2064b9ad928SBarry Smith @*/
2077087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
2084b9ad928SBarry Smith {
2094ac538c5SBarry Smith   PetscErrorCode ierr;
2104b9ad928SBarry Smith 
2114b9ad928SBarry Smith   PetscFunctionBegin;
2120700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
213c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
2144ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
2154b9ad928SBarry Smith   PetscFunctionReturn(0);
2164b9ad928SBarry Smith }
2174b9ad928SBarry Smith 
2184b9ad928SBarry Smith #undef __FUNCT__
2194b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega"
2204b9ad928SBarry Smith /*@
2214b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
2224b9ad928SBarry Smith    (where omega = 1.0 by default).
2234b9ad928SBarry Smith 
2243f9fe445SBarry Smith    Logically Collective on PC
2254b9ad928SBarry Smith 
2264b9ad928SBarry Smith    Input Parameters:
2274b9ad928SBarry Smith +  pc - the preconditioner context
2284b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
2294b9ad928SBarry Smith 
2304b9ad928SBarry Smith    Options Database Key:
2314b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith    Level: intermediate
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
2384b9ad928SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
2404b9ad928SBarry Smith {
2414ac538c5SBarry Smith   PetscErrorCode ierr;
2424b9ad928SBarry Smith 
2434b9ad928SBarry Smith   PetscFunctionBegin;
244c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
245c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
2464ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
2474b9ad928SBarry Smith   PetscFunctionReturn(0);
2484b9ad928SBarry Smith }
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith #undef __FUNCT__
2514b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations"
2524b9ad928SBarry Smith /*@
2534b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
2544b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
2554b9ad928SBarry Smith 
2563f9fe445SBarry Smith    Logically Collective on PC
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith    Input Parameters:
2594b9ad928SBarry Smith +  pc - the preconditioner context
2604b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
2614b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
2624b9ad928SBarry Smith 
2634b9ad928SBarry Smith    Options Database Key:
2644b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
2654b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
2664b9ad928SBarry Smith 
2674b9ad928SBarry Smith    Level: intermediate
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
2744b9ad928SBarry Smith @*/
2757087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
2764b9ad928SBarry Smith {
2774ac538c5SBarry Smith   PetscErrorCode ierr;
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith   PetscFunctionBegin;
2800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
281c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
2824ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
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"
3177087cfbeSBarry Smith PetscErrorCode  PCCreate_SOR(PC pc)
3184b9ad928SBarry Smith {
319dfbe8321SBarry Smith   PetscErrorCode ierr;
3204b9ad928SBarry Smith   PC_SOR         *jac;
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith   PetscFunctionBegin;
32338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr);
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
3264b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
3274b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
3284b9ad928SBarry Smith   pc->ops->setup           = 0;
3294b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
3304b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
3314b9ad928SBarry Smith   pc->data                 = (void*)jac;
332d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
3334b9ad928SBarry Smith   jac->omega               = 1.0;
33496fc60bcSBarry Smith   jac->fshift              = 0.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