xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 9d2471e09ceed5d1c3a3f2e8591e84d1de740827)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Defines a  (S)SOR  preconditioner for any Mat implementation
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #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;
333060034bSBarry Smith   PetscReal      fshift;
344b9ad928SBarry Smith 
354b9ad928SBarry Smith   PetscFunctionBegin;
36b0e188deSSatish Balay   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
37422a814eSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
38539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
394b9ad928SBarry Smith   PetscFunctionReturn(0);
404b9ad928SBarry Smith }
414b9ad928SBarry Smith 
424b9ad928SBarry Smith #undef __FUNCT__
43*9d2471e0SBarry Smith #define __FUNCT__ "PCApplyTranspose_SOR"
44*9d2471e0SBarry Smith static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
45*9d2471e0SBarry Smith {
46*9d2471e0SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
47*9d2471e0SBarry Smith   PetscErrorCode ierr;
48*9d2471e0SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
49*9d2471e0SBarry Smith   PetscReal      fshift;
50*9d2471e0SBarry Smith   PetscBool      set,sym;
51*9d2471e0SBarry Smith 
52*9d2471e0SBarry Smith   PetscFunctionBegin;
53*9d2471e0SBarry Smith   ierr = MatIsSymmetricKnown(pc->pmat,&set,&sym);CHKERRQ(ierr);
54*9d2471e0SBarry Smith   if (!set || !sym || (jac->sym != SOR_SYMMETRIC_SWEEP && jac->sym != SOR_LOCAL_SYMMETRIC_SWEEP)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Can only apply transpose of SOR if matrix is symmetric and sweep is symmetric");
55*9d2471e0SBarry Smith   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
56*9d2471e0SBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
57*9d2471e0SBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
58*9d2471e0SBarry Smith   PetscFunctionReturn(0);
59*9d2471e0SBarry Smith }
60*9d2471e0SBarry Smith 
61*9d2471e0SBarry Smith #undef __FUNCT__
624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR"
63ace3abfcSBarry 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)
644b9ad928SBarry Smith {
654b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
677319c654SBarry Smith   MatSORType     stype = jac->sym;
683060034bSBarry Smith   PetscReal      fshift;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
722fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
73b0e188deSSatish Balay   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
74422a814eSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
75539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
764d0a8057SBarry Smith   *outits = its;
774d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
784b9ad928SBarry Smith   PetscFunctionReturn(0);
794b9ad928SBarry Smith }
804b9ad928SBarry Smith 
814b9ad928SBarry Smith #undef __FUNCT__
824b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR"
834416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
844b9ad928SBarry Smith {
854b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
86dfbe8321SBarry Smith   PetscErrorCode ierr;
87ace3abfcSBarry Smith   PetscBool      flg;
884b9ad928SBarry Smith 
894b9ad928SBarry Smith   PetscFunctionBegin;
90e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
9194ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
9294ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
9394ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
9494ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
95acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
964b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
97acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
984b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
99acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
100a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
101acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
1024b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
103acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
1044b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
105acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
1064b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
1074b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1084b9ad928SBarry Smith   PetscFunctionReturn(0);
1094b9ad928SBarry Smith }
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith #undef __FUNCT__
1124b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR"
113dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
1144b9ad928SBarry Smith {
1154b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
1164b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
1172fc52814SBarry Smith   const char     *sortype;
118dfbe8321SBarry Smith   PetscErrorCode ierr;
119ace3abfcSBarry Smith   PetscBool      iascii;
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith   PetscFunctionBegin;
122251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12332077d6dSBarry Smith   if (iascii) {
1244b9ad928SBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
1254b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1264b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1274b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
128db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1294b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1304b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
131db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1324b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1334b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1344b9ad928SBarry Smith     else                                                                     sortype = "unknown";
13557622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  SOR: type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
1364b9ad928SBarry Smith   }
1374b9ad928SBarry Smith   PetscFunctionReturn(0);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith 
1414b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1424b9ad928SBarry Smith #undef __FUNCT__
1434b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR"
144f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1454b9ad928SBarry Smith {
146c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   PetscFunctionBegin;
1494b9ad928SBarry Smith   jac->sym = flag;
1504b9ad928SBarry Smith   PetscFunctionReturn(0);
1514b9ad928SBarry Smith }
1524b9ad928SBarry Smith 
1534b9ad928SBarry Smith #undef __FUNCT__
1544b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR"
155f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1564b9ad928SBarry Smith {
157c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith   PetscFunctionBegin;
160ce94432eSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1614b9ad928SBarry Smith   jac->omega = omega;
1624b9ad928SBarry Smith   PetscFunctionReturn(0);
1634b9ad928SBarry Smith }
1644b9ad928SBarry Smith 
1654b9ad928SBarry Smith #undef __FUNCT__
1664b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR"
167f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1684b9ad928SBarry Smith {
169c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1704b9ad928SBarry Smith 
1714b9ad928SBarry Smith   PetscFunctionBegin;
1724b9ad928SBarry Smith   jac->its  = its;
1734b9ad928SBarry Smith   jac->lits = lits;
1744b9ad928SBarry Smith   PetscFunctionReturn(0);
1754b9ad928SBarry Smith }
1764b9ad928SBarry Smith 
177c60c7ad4SBarry Smith #undef __FUNCT__
178c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric_SOR"
179c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
180c60c7ad4SBarry Smith {
181c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
182c60c7ad4SBarry Smith 
183c60c7ad4SBarry Smith   PetscFunctionBegin;
184c60c7ad4SBarry Smith   *flag = jac->sym;
185c60c7ad4SBarry Smith   PetscFunctionReturn(0);
186c60c7ad4SBarry Smith }
187c60c7ad4SBarry Smith 
188c60c7ad4SBarry Smith #undef __FUNCT__
189c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega_SOR"
190c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
191c60c7ad4SBarry Smith {
192c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
193c60c7ad4SBarry Smith 
194c60c7ad4SBarry Smith   PetscFunctionBegin;
195c60c7ad4SBarry Smith   *omega = jac->omega;
196c60c7ad4SBarry Smith   PetscFunctionReturn(0);
197c60c7ad4SBarry Smith }
198c60c7ad4SBarry Smith 
199c60c7ad4SBarry Smith #undef __FUNCT__
200c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations_SOR"
201c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
202c60c7ad4SBarry Smith {
203c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
204c60c7ad4SBarry Smith 
205c60c7ad4SBarry Smith   PetscFunctionBegin;
206c60c7ad4SBarry Smith   if (its)  *its = jac->its;
207c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
208c60c7ad4SBarry Smith   PetscFunctionReturn(0);
209c60c7ad4SBarry Smith }
210c60c7ad4SBarry Smith 
2114b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
2124b9ad928SBarry Smith #undef __FUNCT__
213c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric"
214c60c7ad4SBarry Smith /*@
2151ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
216c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
217c60c7ad4SBarry Smith 
218c60c7ad4SBarry Smith    Logically Collective on PC
219c60c7ad4SBarry Smith 
220c60c7ad4SBarry Smith    Input Parameter:
221c60c7ad4SBarry Smith .  pc - the preconditioner context
222c60c7ad4SBarry Smith 
223c60c7ad4SBarry Smith    Output Parameter:
224c60c7ad4SBarry Smith .  flag - one of the following
225c60c7ad4SBarry Smith .vb
226c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
227c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
228c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
229c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
230c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
231c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
232c60c7ad4SBarry Smith .ve
233c60c7ad4SBarry Smith 
234c60c7ad4SBarry Smith    Options Database Keys:
235c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
236c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
237c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
238c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
239c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
240c60c7ad4SBarry Smith 
241c60c7ad4SBarry Smith    Notes:
242c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
243c60c7ad4SBarry Smith    which can be chosen with the option
244c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
245c60c7ad4SBarry Smith 
246c60c7ad4SBarry Smith    Level: intermediate
247c60c7ad4SBarry Smith 
248c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
249c60c7ad4SBarry Smith 
250c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
251c60c7ad4SBarry Smith @*/
252c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
253c60c7ad4SBarry Smith {
254c60c7ad4SBarry Smith   PetscErrorCode ierr;
255c60c7ad4SBarry Smith 
256c60c7ad4SBarry Smith   PetscFunctionBegin;
257c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
258c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
259c60c7ad4SBarry Smith   PetscFunctionReturn(0);
260c60c7ad4SBarry Smith }
261c60c7ad4SBarry Smith 
262c60c7ad4SBarry Smith #undef __FUNCT__
263c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega"
264c60c7ad4SBarry Smith /*@
265c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
266c60c7ad4SBarry Smith    (where omega = 1.0 by default).
267c60c7ad4SBarry Smith 
268c60c7ad4SBarry Smith    Logically Collective on PC
269c60c7ad4SBarry Smith 
270c60c7ad4SBarry Smith    Input Parameter:
271c60c7ad4SBarry Smith .  pc - the preconditioner context
272c60c7ad4SBarry Smith 
273c60c7ad4SBarry Smith    Output Parameter:
274c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
275c60c7ad4SBarry Smith 
276c60c7ad4SBarry Smith    Options Database Key:
277c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
278c60c7ad4SBarry Smith 
279c60c7ad4SBarry Smith    Level: intermediate
280c60c7ad4SBarry Smith 
281c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
282c60c7ad4SBarry Smith 
283c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
284c60c7ad4SBarry Smith @*/
285c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
286c60c7ad4SBarry Smith {
287c60c7ad4SBarry Smith   PetscErrorCode ierr;
288c60c7ad4SBarry Smith 
289c60c7ad4SBarry Smith   PetscFunctionBegin;
290c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
291c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
292c60c7ad4SBarry Smith   PetscFunctionReturn(0);
293c60c7ad4SBarry Smith }
294c60c7ad4SBarry Smith 
295c60c7ad4SBarry Smith #undef __FUNCT__
296c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations"
297c60c7ad4SBarry Smith /*@
298c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
299c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
300c60c7ad4SBarry Smith 
301c60c7ad4SBarry Smith    Logically Collective on PC
302c60c7ad4SBarry Smith 
303c60c7ad4SBarry Smith    Input Parameter:
304c60c7ad4SBarry Smith .  pc - the preconditioner context
305c60c7ad4SBarry Smith 
306c60c7ad4SBarry Smith    Output Parameter:
307c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
308c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
309c60c7ad4SBarry Smith 
310c60c7ad4SBarry Smith    Options Database Key:
311c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
312c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
313c60c7ad4SBarry Smith 
314c60c7ad4SBarry Smith    Level: intermediate
315c60c7ad4SBarry Smith 
316c60c7ad4SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
317c60c7ad4SBarry Smith 
318c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
319c60c7ad4SBarry Smith 
320c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
321c60c7ad4SBarry Smith @*/
322c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
323c60c7ad4SBarry Smith {
324c60c7ad4SBarry Smith   PetscErrorCode ierr;
325c60c7ad4SBarry Smith 
326c60c7ad4SBarry Smith   PetscFunctionBegin;
327c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
328c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
329c60c7ad4SBarry Smith   PetscFunctionReturn(0);
330c60c7ad4SBarry Smith }
331c60c7ad4SBarry Smith 
332c60c7ad4SBarry Smith #undef __FUNCT__
3334b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric"
3344b9ad928SBarry Smith /*@
3354b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
3364b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
3374b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
3384b9ad928SBarry Smith 
3393f9fe445SBarry Smith    Logically Collective on PC
3404b9ad928SBarry Smith 
3414b9ad928SBarry Smith    Input Parameters:
3424b9ad928SBarry Smith +  pc - the preconditioner context
3434b9ad928SBarry Smith -  flag - one of the following
3444b9ad928SBarry Smith .vb
3454b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3464b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3474b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3484b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3494b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3504b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3514b9ad928SBarry Smith .ve
3524b9ad928SBarry Smith 
3534b9ad928SBarry Smith    Options Database Keys:
3544b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3554b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3564b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3574b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3584b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith    Notes:
3614b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3624b9ad928SBarry Smith    which can be chosen with the option
3634b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Level: intermediate
3664b9ad928SBarry Smith 
3674b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3704b9ad928SBarry Smith @*/
3717087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3724b9ad928SBarry Smith {
3734ac538c5SBarry Smith   PetscErrorCode ierr;
3744b9ad928SBarry Smith 
3754b9ad928SBarry Smith   PetscFunctionBegin;
3760700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
377c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3784ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3794b9ad928SBarry Smith   PetscFunctionReturn(0);
3804b9ad928SBarry Smith }
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith #undef __FUNCT__
3834b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega"
3844b9ad928SBarry Smith /*@
3854b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3864b9ad928SBarry Smith    (where omega = 1.0 by default).
3874b9ad928SBarry Smith 
3883f9fe445SBarry Smith    Logically Collective on PC
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith    Input Parameters:
3914b9ad928SBarry Smith +  pc - the preconditioner context
3924b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith    Options Database Key:
3954b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Level: intermediate
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
4024b9ad928SBarry Smith @*/
4037087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
4044b9ad928SBarry Smith {
4054ac538c5SBarry Smith   PetscErrorCode ierr;
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith   PetscFunctionBegin;
408c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
409c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
4104ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
4114b9ad928SBarry Smith   PetscFunctionReturn(0);
4124b9ad928SBarry Smith }
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith #undef __FUNCT__
4154b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations"
4164b9ad928SBarry Smith /*@
4174b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
4184b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
4194b9ad928SBarry Smith 
4203f9fe445SBarry Smith    Logically Collective on PC
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith    Input Parameters:
4234b9ad928SBarry Smith +  pc - the preconditioner context
4244b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
4254b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Options Database Key:
4284b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
4294b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith    Level: intermediate
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
4384b9ad928SBarry Smith @*/
4397087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
4404b9ad928SBarry Smith {
4414ac538c5SBarry Smith   PetscErrorCode ierr;
4424b9ad928SBarry Smith 
4434b9ad928SBarry Smith   PetscFunctionBegin;
4440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
445c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
4464ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
4474b9ad928SBarry Smith   PetscFunctionReturn(0);
4484b9ad928SBarry Smith }
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith /*MC
4514b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4524b9ad928SBarry Smith 
4534b9ad928SBarry Smith    Options Database Keys:
4544b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4554b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
456a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4574b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
458a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4594b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4604b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
461422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
462a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
463a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith    Level: beginner
4664b9ad928SBarry Smith 
4674b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
4684b9ad928SBarry Smith 
46937a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
4704b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4714b9ad928SBarry Smith           Jacobi with SOR on each block.
4724b9ad928SBarry Smith 
473422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
474422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
475422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
476422a814eSBarry Smith           zero pivot.
477422a814eSBarry Smith 
47837a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
47937a17b4dSBarry Smith 
480422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
481422a814eSBarry Smith           the computation is stopped with an error
482422a814eSBarry Smith 
4834b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4844b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
4854b9ad928SBarry Smith M*/
4864b9ad928SBarry Smith 
4874b9ad928SBarry Smith #undef __FUNCT__
4884b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR"
4898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4904b9ad928SBarry Smith {
491dfbe8321SBarry Smith   PetscErrorCode ierr;
4924b9ad928SBarry Smith   PC_SOR         *jac;
4934b9ad928SBarry Smith 
4944b9ad928SBarry Smith   PetscFunctionBegin;
495b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4964b9ad928SBarry Smith 
4974b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
498*9d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4994b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
5004b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
5014b9ad928SBarry Smith   pc->ops->setup           = 0;
5024b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
5034b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
5044b9ad928SBarry Smith   pc->data                 = (void*)jac;
505d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
5064b9ad928SBarry Smith   jac->omega               = 1.0;
50796fc60bcSBarry Smith   jac->fshift              = 0.0;
5084b9ad928SBarry Smith   jac->its                 = 1;
5094b9ad928SBarry Smith   jac->lits                = 1;
5104b9ad928SBarry Smith 
511bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
512bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
513bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
514c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
515c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
516c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
5174b9ad928SBarry Smith   PetscFunctionReturn(0);
5184b9ad928SBarry Smith }
5194b9ad928SBarry Smith 
5204b9ad928SBarry Smith 
5214b9ad928SBarry Smith 
5224b9ad928SBarry Smith 
5234b9ad928SBarry Smith 
524