xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 1ff2113e2c5156446e3465ca400c359a0cfbe633)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Defines a  (S)SOR  preconditioner for any Mat implementation
44b9ad928SBarry Smith */
5b45d2f2cSJed 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);
492fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
5041f059aeSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
514d0a8057SBarry Smith   *outits = its;
524d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
534b9ad928SBarry Smith   PetscFunctionReturn(0);
544b9ad928SBarry Smith }
554b9ad928SBarry Smith 
564b9ad928SBarry Smith #undef __FUNCT__
574b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR"
58dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc)
594b9ad928SBarry Smith {
604b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
61dfbe8321SBarry Smith   PetscErrorCode ierr;
62ace3abfcSBarry Smith   PetscBool      flg;
634b9ad928SBarry Smith 
644b9ad928SBarry Smith   PetscFunctionBegin;
654b9ad928SBarry Smith   ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr);
6694ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
6794ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
6894ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
6994ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
70acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
714b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
72acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
734b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
74acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
75a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
76acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
774b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
78acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
794b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
80acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
814b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
824b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
834b9ad928SBarry Smith   PetscFunctionReturn(0);
844b9ad928SBarry Smith }
854b9ad928SBarry Smith 
864b9ad928SBarry Smith #undef __FUNCT__
874b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR"
88dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
894b9ad928SBarry Smith {
904b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
914b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
922fc52814SBarry Smith   const char     *sortype;
93dfbe8321SBarry Smith   PetscErrorCode ierr;
94ace3abfcSBarry Smith   PetscBool      iascii;
954b9ad928SBarry Smith 
964b9ad928SBarry Smith   PetscFunctionBegin;
97251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9832077d6dSBarry Smith   if (iascii) {
994b9ad928SBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  SOR:  zero initial guess\n");CHKERRQ(ierr);}
1004b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1014b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1024b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
103db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1044b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1054b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
106db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1074b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1084b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1094b9ad928SBarry Smith     else                                                                     sortype = "unknown";
11057622a8eSBarry 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);
1114b9ad928SBarry Smith   }
1124b9ad928SBarry Smith   PetscFunctionReturn(0);
1134b9ad928SBarry Smith }
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith 
1164b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1174b9ad928SBarry Smith #undef __FUNCT__
1184b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR"
119f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1204b9ad928SBarry Smith {
121c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith   PetscFunctionBegin;
1244b9ad928SBarry Smith   jac->sym = flag;
1254b9ad928SBarry Smith   PetscFunctionReturn(0);
1264b9ad928SBarry Smith }
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith #undef __FUNCT__
1294b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR"
130f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1314b9ad928SBarry Smith {
132c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   PetscFunctionBegin;
135ce94432eSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1364b9ad928SBarry Smith   jac->omega = omega;
1374b9ad928SBarry Smith   PetscFunctionReturn(0);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith #undef __FUNCT__
1414b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR"
142f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1434b9ad928SBarry Smith {
144c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1454b9ad928SBarry Smith 
1464b9ad928SBarry Smith   PetscFunctionBegin;
1474b9ad928SBarry Smith   jac->its  = its;
1484b9ad928SBarry Smith   jac->lits = lits;
1494b9ad928SBarry Smith   PetscFunctionReturn(0);
1504b9ad928SBarry Smith }
1514b9ad928SBarry Smith 
152c60c7ad4SBarry Smith #undef __FUNCT__
153c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric_SOR"
154c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
155c60c7ad4SBarry Smith {
156c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
157c60c7ad4SBarry Smith 
158c60c7ad4SBarry Smith   PetscFunctionBegin;
159c60c7ad4SBarry Smith   *flag = jac->sym;
160c60c7ad4SBarry Smith   PetscFunctionReturn(0);
161c60c7ad4SBarry Smith }
162c60c7ad4SBarry Smith 
163c60c7ad4SBarry Smith #undef __FUNCT__
164c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega_SOR"
165c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
166c60c7ad4SBarry Smith {
167c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
168c60c7ad4SBarry Smith 
169c60c7ad4SBarry Smith   PetscFunctionBegin;
170c60c7ad4SBarry Smith   *omega = jac->omega;
171c60c7ad4SBarry Smith   PetscFunctionReturn(0);
172c60c7ad4SBarry Smith }
173c60c7ad4SBarry Smith 
174c60c7ad4SBarry Smith #undef __FUNCT__
175c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations_SOR"
176c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
177c60c7ad4SBarry Smith {
178c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
179c60c7ad4SBarry Smith 
180c60c7ad4SBarry Smith   PetscFunctionBegin;
181c60c7ad4SBarry Smith   if (its)  *its = jac->its;
182c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
183c60c7ad4SBarry Smith   PetscFunctionReturn(0);
184c60c7ad4SBarry Smith }
185c60c7ad4SBarry Smith 
1864b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
1874b9ad928SBarry Smith #undef __FUNCT__
188c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric"
189c60c7ad4SBarry Smith /*@
190*1ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
191c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
192c60c7ad4SBarry Smith 
193c60c7ad4SBarry Smith    Logically Collective on PC
194c60c7ad4SBarry Smith 
195c60c7ad4SBarry Smith    Input Parameter:
196c60c7ad4SBarry Smith .  pc - the preconditioner context
197c60c7ad4SBarry Smith 
198c60c7ad4SBarry Smith    Output Parameter:
199c60c7ad4SBarry Smith .  flag - one of the following
200c60c7ad4SBarry Smith .vb
201c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
202c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
203c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
204c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
205c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
206c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
207c60c7ad4SBarry Smith .ve
208c60c7ad4SBarry Smith 
209c60c7ad4SBarry Smith    Options Database Keys:
210c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
211c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
212c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
213c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
214c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
215c60c7ad4SBarry Smith 
216c60c7ad4SBarry Smith    Notes:
217c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
218c60c7ad4SBarry Smith    which can be chosen with the option
219c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
220c60c7ad4SBarry Smith 
221c60c7ad4SBarry Smith    Level: intermediate
222c60c7ad4SBarry Smith 
223c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
224c60c7ad4SBarry Smith 
225c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
226c60c7ad4SBarry Smith @*/
227c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
228c60c7ad4SBarry Smith {
229c60c7ad4SBarry Smith   PetscErrorCode ierr;
230c60c7ad4SBarry Smith 
231c60c7ad4SBarry Smith   PetscFunctionBegin;
232c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
233c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
234c60c7ad4SBarry Smith   PetscFunctionReturn(0);
235c60c7ad4SBarry Smith }
236c60c7ad4SBarry Smith 
237c60c7ad4SBarry Smith #undef __FUNCT__
238c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega"
239c60c7ad4SBarry Smith /*@
240c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
241c60c7ad4SBarry Smith    (where omega = 1.0 by default).
242c60c7ad4SBarry Smith 
243c60c7ad4SBarry Smith    Logically Collective on PC
244c60c7ad4SBarry Smith 
245c60c7ad4SBarry Smith    Input Parameter:
246c60c7ad4SBarry Smith .  pc - the preconditioner context
247c60c7ad4SBarry Smith 
248c60c7ad4SBarry Smith    Output Parameter:
249c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
250c60c7ad4SBarry Smith 
251c60c7ad4SBarry Smith    Options Database Key:
252c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
253c60c7ad4SBarry Smith 
254c60c7ad4SBarry Smith    Level: intermediate
255c60c7ad4SBarry Smith 
256c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
257c60c7ad4SBarry Smith 
258c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
259c60c7ad4SBarry Smith @*/
260c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
261c60c7ad4SBarry Smith {
262c60c7ad4SBarry Smith   PetscErrorCode ierr;
263c60c7ad4SBarry Smith 
264c60c7ad4SBarry Smith   PetscFunctionBegin;
265c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
266c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
267c60c7ad4SBarry Smith   PetscFunctionReturn(0);
268c60c7ad4SBarry Smith }
269c60c7ad4SBarry Smith 
270c60c7ad4SBarry Smith #undef __FUNCT__
271c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations"
272c60c7ad4SBarry Smith /*@
273c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
274c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
275c60c7ad4SBarry Smith 
276c60c7ad4SBarry Smith    Logically Collective on PC
277c60c7ad4SBarry Smith 
278c60c7ad4SBarry Smith    Input Parameter:
279c60c7ad4SBarry Smith .  pc - the preconditioner context
280c60c7ad4SBarry Smith 
281c60c7ad4SBarry Smith    Output Parameter:
282c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
283c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
284c60c7ad4SBarry Smith 
285c60c7ad4SBarry Smith    Options Database Key:
286c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
287c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
288c60c7ad4SBarry Smith 
289c60c7ad4SBarry Smith    Level: intermediate
290c60c7ad4SBarry Smith 
291c60c7ad4SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
292c60c7ad4SBarry Smith 
293c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
294c60c7ad4SBarry Smith 
295c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
296c60c7ad4SBarry Smith @*/
297c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
298c60c7ad4SBarry Smith {
299c60c7ad4SBarry Smith   PetscErrorCode ierr;
300c60c7ad4SBarry Smith 
301c60c7ad4SBarry Smith   PetscFunctionBegin;
302c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
303c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
304c60c7ad4SBarry Smith   PetscFunctionReturn(0);
305c60c7ad4SBarry Smith }
306c60c7ad4SBarry Smith 
307c60c7ad4SBarry Smith #undef __FUNCT__
3084b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric"
3094b9ad928SBarry Smith /*@
3104b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
3114b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
3124b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
3134b9ad928SBarry Smith 
3143f9fe445SBarry Smith    Logically Collective on PC
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith    Input Parameters:
3174b9ad928SBarry Smith +  pc - the preconditioner context
3184b9ad928SBarry Smith -  flag - one of the following
3194b9ad928SBarry Smith .vb
3204b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3214b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3224b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3234b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3244b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3254b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3264b9ad928SBarry Smith .ve
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith    Options Database Keys:
3294b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3304b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3314b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3324b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3334b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3344b9ad928SBarry Smith 
3354b9ad928SBarry Smith    Notes:
3364b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3374b9ad928SBarry Smith    which can be chosen with the option
3384b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith    Level: intermediate
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3454b9ad928SBarry Smith @*/
3467087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3474b9ad928SBarry Smith {
3484ac538c5SBarry Smith   PetscErrorCode ierr;
3494b9ad928SBarry Smith 
3504b9ad928SBarry Smith   PetscFunctionBegin;
3510700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
352c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3534ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3544b9ad928SBarry Smith   PetscFunctionReturn(0);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith #undef __FUNCT__
3584b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega"
3594b9ad928SBarry Smith /*@
3604b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3614b9ad928SBarry Smith    (where omega = 1.0 by default).
3624b9ad928SBarry Smith 
3633f9fe445SBarry Smith    Logically Collective on PC
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Input Parameters:
3664b9ad928SBarry Smith +  pc - the preconditioner context
3674b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith    Options Database Key:
3704b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3714b9ad928SBarry Smith 
3724b9ad928SBarry Smith    Level: intermediate
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
3774b9ad928SBarry Smith @*/
3787087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3794b9ad928SBarry Smith {
3804ac538c5SBarry Smith   PetscErrorCode ierr;
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith   PetscFunctionBegin;
383c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
384c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3854ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
3864b9ad928SBarry Smith   PetscFunctionReturn(0);
3874b9ad928SBarry Smith }
3884b9ad928SBarry Smith 
3894b9ad928SBarry Smith #undef __FUNCT__
3904b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations"
3914b9ad928SBarry Smith /*@
3924b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3934b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3944b9ad928SBarry Smith 
3953f9fe445SBarry Smith    Logically Collective on PC
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Input Parameters:
3984b9ad928SBarry Smith +  pc - the preconditioner context
3994b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
4004b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith    Options Database Key:
4034b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
4044b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
4054b9ad928SBarry Smith 
4064b9ad928SBarry Smith    Level: intermediate
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
4114b9ad928SBarry Smith 
4124b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
4134b9ad928SBarry Smith @*/
4147087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
4154b9ad928SBarry Smith {
4164ac538c5SBarry Smith   PetscErrorCode ierr;
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith   PetscFunctionBegin;
4190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
420c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
4214ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
4224b9ad928SBarry Smith   PetscFunctionReturn(0);
4234b9ad928SBarry Smith }
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith /*MC
4264b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith    Options Database Keys:
4294b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4304b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
431a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4324b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
433a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4344b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4354b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
436a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
437a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4384b9ad928SBarry Smith 
4394b9ad928SBarry Smith    Level: beginner
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
4424b9ad928SBarry Smith 
44337a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
4444b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4454b9ad928SBarry Smith           Jacobi with SOR on each block.
4464b9ad928SBarry Smith 
44737a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
44837a17b4dSBarry Smith 
4494b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4504b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
4514b9ad928SBarry Smith M*/
4524b9ad928SBarry Smith 
4534b9ad928SBarry Smith #undef __FUNCT__
4544b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR"
4558cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4564b9ad928SBarry Smith {
457dfbe8321SBarry Smith   PetscErrorCode ierr;
4584b9ad928SBarry Smith   PC_SOR         *jac;
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith   PetscFunctionBegin;
461b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4624b9ad928SBarry Smith 
4634b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4644b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4654b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4664b9ad928SBarry Smith   pc->ops->setup           = 0;
4674b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4684b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4694b9ad928SBarry Smith   pc->data                 = (void*)jac;
470d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4714b9ad928SBarry Smith   jac->omega               = 1.0;
47296fc60bcSBarry Smith   jac->fshift              = 0.0;
4734b9ad928SBarry Smith   jac->its                 = 1;
4744b9ad928SBarry Smith   jac->lits                = 1;
4754b9ad928SBarry Smith 
476bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
477bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
478bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
479c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
480c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
481c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
4824b9ad928SBarry Smith   PetscFunctionReturn(0);
4834b9ad928SBarry Smith }
4844b9ad928SBarry Smith 
4854b9ad928SBarry Smith 
4864b9ad928SBarry Smith 
4874b9ad928SBarry Smith 
4884b9ad928SBarry Smith 
489