xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
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;
1329c1d7e0SHong Zhang   PetscReal   fshift;
144b9ad928SBarry Smith } PC_SOR;
154b9ad928SBarry Smith 
164b9ad928SBarry Smith #undef __FUNCT__
174b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR"
186849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc)
194b9ad928SBarry Smith {
204b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
21dfbe8321SBarry Smith   PetscErrorCode ierr;
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   PetscFunctionBegin;
244b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
254b9ad928SBarry Smith   PetscFunctionReturn(0);
264b9ad928SBarry Smith }
274b9ad928SBarry Smith 
284b9ad928SBarry Smith #undef __FUNCT__
294b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR"
306849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
314b9ad928SBarry Smith {
324b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
33dfbe8321SBarry Smith   PetscErrorCode ierr;
34c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
354b9ad928SBarry Smith 
364b9ad928SBarry Smith   PetscFunctionBegin;
3741f059aeSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
384b9ad928SBarry Smith   PetscFunctionReturn(0);
394b9ad928SBarry Smith }
404b9ad928SBarry Smith 
414b9ad928SBarry Smith #undef __FUNCT__
424b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR"
43ace3abfcSBarry 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)
444b9ad928SBarry Smith {
454b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
46dfbe8321SBarry Smith   PetscErrorCode ierr;
477319c654SBarry Smith   MatSORType     stype = jac->sym;
484b9ad928SBarry Smith 
494b9ad928SBarry Smith   PetscFunctionBegin;
50ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
517319c654SBarry Smith   if (guesszero) {
52de8939e3SBarry Smith     stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
537319c654SBarry Smith   }
5441f059aeSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
554d0a8057SBarry Smith   *outits = its;
564d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
574b9ad928SBarry Smith   PetscFunctionReturn(0);
584b9ad928SBarry Smith }
594b9ad928SBarry Smith 
604b9ad928SBarry Smith #undef __FUNCT__
614b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR"
62dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc)
634b9ad928SBarry Smith {
644b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
65dfbe8321SBarry Smith   PetscErrorCode ierr;
66ace3abfcSBarry Smith   PetscBool      flg;
674b9ad928SBarry Smith 
684b9ad928SBarry Smith   PetscFunctionBegin;
694b9ad928SBarry Smith   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
704b9ad928SBarry Smith     ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr);
7196fc60bcSBarry Smith     ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr);
724b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr);
734b9ad928SBarry Smith     ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr);
74*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
754b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
76*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
774b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
78*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
79a9510f2eSBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
80*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
814b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
82*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
834b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
84*acfcf0e5SJed Brown     ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
854b9ad928SBarry Smith     if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
864b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
874b9ad928SBarry Smith   PetscFunctionReturn(0);
884b9ad928SBarry Smith }
894b9ad928SBarry Smith 
904b9ad928SBarry Smith #undef __FUNCT__
914b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR"
92dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
934b9ad928SBarry Smith {
944b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
954b9ad928SBarry Smith   MatSORType     sym = jac->sym;
962fc52814SBarry Smith   const char     *sortype;
97dfbe8321SBarry Smith   PetscErrorCode ierr;
98ace3abfcSBarry Smith   PetscBool      iascii;
994b9ad928SBarry Smith 
1004b9ad928SBarry Smith   PetscFunctionBegin;
1012692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10232077d6dSBarry Smith   if (iascii) {
1034b9ad928SBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
1044b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)              sortype = "apply_upper";
1054b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)         sortype = "apply_lower";
1064b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)            sortype = "Eisenstat";
1074b9ad928SBarry Smith     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)
1084b9ad928SBarry Smith                                              sortype = "symmetric";
1094b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)       sortype = "backward";
1104b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)        sortype = "forward";
1114b9ad928SBarry Smith     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP)
1124b9ad928SBarry Smith                                              sortype = "local_symmetric";
1134b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)  sortype = "local_forward";
1144b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward";
1154b9ad928SBarry Smith     else                                     sortype = "unknown";
116cf576665SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);CHKERRQ(ierr);
1174b9ad928SBarry Smith   } else {
11865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name);
1194b9ad928SBarry Smith   }
1204b9ad928SBarry Smith   PetscFunctionReturn(0);
1214b9ad928SBarry Smith }
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith 
1244b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1254b9ad928SBarry Smith EXTERN_C_BEGIN
1264b9ad928SBarry Smith #undef __FUNCT__
1274b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR"
128dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1294b9ad928SBarry Smith {
1304b9ad928SBarry Smith   PC_SOR *jac;
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith   PetscFunctionBegin;
1334b9ad928SBarry Smith   jac = (PC_SOR*)pc->data;
1344b9ad928SBarry Smith   jac->sym = flag;
1354b9ad928SBarry Smith   PetscFunctionReturn(0);
1364b9ad928SBarry Smith }
1374b9ad928SBarry Smith EXTERN_C_END
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith EXTERN_C_BEGIN
1404b9ad928SBarry Smith #undef __FUNCT__
1414b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR"
142dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega_SOR(PC pc,PetscReal omega)
1434b9ad928SBarry Smith {
1444b9ad928SBarry Smith   PC_SOR *jac;
1454b9ad928SBarry Smith 
1464b9ad928SBarry Smith   PetscFunctionBegin;
147e7e72b3dSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1484b9ad928SBarry Smith   jac        = (PC_SOR*)pc->data;
1494b9ad928SBarry Smith   jac->omega = omega;
1504b9ad928SBarry Smith   PetscFunctionReturn(0);
1514b9ad928SBarry Smith }
1524b9ad928SBarry Smith EXTERN_C_END
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith EXTERN_C_BEGIN
1554b9ad928SBarry Smith #undef __FUNCT__
1564b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR"
157dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1584b9ad928SBarry Smith {
1594b9ad928SBarry Smith   PC_SOR *jac;
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith   PetscFunctionBegin;
1624b9ad928SBarry Smith   jac      = (PC_SOR*)pc->data;
1634b9ad928SBarry Smith   jac->its = its;
1644b9ad928SBarry Smith   jac->lits = lits;
1654b9ad928SBarry Smith   PetscFunctionReturn(0);
1664b9ad928SBarry Smith }
1674b9ad928SBarry Smith EXTERN_C_END
1684b9ad928SBarry Smith 
1694b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1704b9ad928SBarry Smith #undef __FUNCT__
1714b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric"
1724b9ad928SBarry Smith /*@
1734b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
1744b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
1754b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
1764b9ad928SBarry Smith 
1773f9fe445SBarry Smith    Logically Collective on PC
1784b9ad928SBarry Smith 
1794b9ad928SBarry Smith    Input Parameters:
1804b9ad928SBarry Smith +  pc - the preconditioner context
1814b9ad928SBarry Smith -  flag - one of the following
1824b9ad928SBarry Smith .vb
1834b9ad928SBarry Smith     SOR_FORWARD_SWEEP
1844b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
1854b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
1864b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
1874b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
1884b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
1894b9ad928SBarry Smith .ve
1904b9ad928SBarry Smith 
1914b9ad928SBarry Smith    Options Database Keys:
1924b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
1934b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
1944b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
1954b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
1964b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
1974b9ad928SBarry Smith 
1984b9ad928SBarry Smith    Notes:
1994b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
2004b9ad928SBarry Smith    which can be chosen with the option
2014b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith    Level: intermediate
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
2064b9ad928SBarry Smith 
2074b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
2084b9ad928SBarry Smith @*/
209dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC pc,MatSORType flag)
2104b9ad928SBarry Smith {
2114ac538c5SBarry Smith   PetscErrorCode ierr;
2124b9ad928SBarry Smith 
2134b9ad928SBarry Smith   PetscFunctionBegin;
2140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
215c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
2164ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
2174b9ad928SBarry Smith   PetscFunctionReturn(0);
2184b9ad928SBarry Smith }
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith #undef __FUNCT__
2214b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega"
2224b9ad928SBarry Smith /*@
2234b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
2244b9ad928SBarry Smith    (where omega = 1.0 by default).
2254b9ad928SBarry Smith 
2263f9fe445SBarry Smith    Logically Collective on PC
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith    Input Parameters:
2294b9ad928SBarry Smith +  pc - the preconditioner context
2304b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith    Options Database Key:
2334b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith    Level: intermediate
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
2384b9ad928SBarry Smith 
2394b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
2404b9ad928SBarry Smith @*/
241dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC pc,PetscReal omega)
2424b9ad928SBarry Smith {
2434ac538c5SBarry Smith   PetscErrorCode ierr;
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith   PetscFunctionBegin;
246c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
247c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
2484ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
2494b9ad928SBarry Smith   PetscFunctionReturn(0);
2504b9ad928SBarry Smith }
2514b9ad928SBarry Smith 
2524b9ad928SBarry Smith #undef __FUNCT__
2534b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations"
2544b9ad928SBarry Smith /*@
2554b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
2564b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
2574b9ad928SBarry Smith 
2583f9fe445SBarry Smith    Logically Collective on PC
2594b9ad928SBarry Smith 
2604b9ad928SBarry Smith    Input Parameters:
2614b9ad928SBarry Smith +  pc - the preconditioner context
2624b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
2634b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Options Database Key:
2664b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
2674b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith    Level: intermediate
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
2764b9ad928SBarry Smith @*/
277dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
2784b9ad928SBarry Smith {
2794ac538c5SBarry Smith   PetscErrorCode ierr;
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith   PetscFunctionBegin;
2820700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
283c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
2844ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
2854b9ad928SBarry Smith   PetscFunctionReturn(0);
2864b9ad928SBarry Smith }
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith /*MC
2894b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
2904b9ad928SBarry Smith 
2914b9ad928SBarry Smith    Options Database Keys:
2924b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
2934b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
294a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
2954b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
296a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
2974b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
2984b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
299a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
300a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith    Level: beginner
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
3054b9ad928SBarry Smith 
30637a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
3074b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
3084b9ad928SBarry Smith           Jacobi with SOR on each block.
3094b9ad928SBarry Smith 
31037a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
31137a17b4dSBarry Smith 
3124b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
3134b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
3144b9ad928SBarry Smith M*/
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith EXTERN_C_BEGIN
3174b9ad928SBarry Smith #undef __FUNCT__
3184b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR"
319dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SOR(PC pc)
3204b9ad928SBarry Smith {
321dfbe8321SBarry Smith   PetscErrorCode ierr;
3224b9ad928SBarry Smith   PC_SOR         *jac;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
32538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr);
3264b9ad928SBarry Smith 
3274b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
3284b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
3294b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
3304b9ad928SBarry Smith   pc->ops->setup           = 0;
3314b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
3324b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
3334b9ad928SBarry Smith   pc->data                 = (void*)jac;
334d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
3354b9ad928SBarry Smith   jac->omega               = 1.0;
33696fc60bcSBarry Smith   jac->fshift              = 0.0;
3374b9ad928SBarry Smith   jac->its                 = 1;
3384b9ad928SBarry Smith   jac->lits                = 1;
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR",
3414b9ad928SBarry Smith                     PCSORSetSymmetric_SOR);CHKERRQ(ierr);
3424b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR",
3434b9ad928SBarry Smith                     PCSORSetOmega_SOR);CHKERRQ(ierr);
3444b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR",
3454b9ad928SBarry Smith                     PCSORSetIterations_SOR);CHKERRQ(ierr);
3464b9ad928SBarry Smith 
3474b9ad928SBarry Smith   PetscFunctionReturn(0);
3484b9ad928SBarry Smith }
3494b9ad928SBarry Smith EXTERN_C_END
3504b9ad928SBarry Smith 
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith 
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith 
355