xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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 
156849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc)
164b9ad928SBarry Smith {
17dfbe8321SBarry Smith   PetscErrorCode ierr;
184b9ad928SBarry Smith 
194b9ad928SBarry Smith   PetscFunctionBegin;
20c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
214b9ad928SBarry Smith   PetscFunctionReturn(0);
224b9ad928SBarry Smith }
234b9ad928SBarry Smith 
246849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
254b9ad928SBarry Smith {
264b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
27dfbe8321SBarry Smith   PetscErrorCode ierr;
28c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
293060034bSBarry Smith   PetscReal      fshift;
304b9ad928SBarry Smith 
314b9ad928SBarry Smith   PetscFunctionBegin;
32b0e188deSSatish Balay   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
33422a814eSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
34539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
354b9ad928SBarry Smith   PetscFunctionReturn(0);
364b9ad928SBarry Smith }
374b9ad928SBarry Smith 
389d2471e0SBarry Smith static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
399d2471e0SBarry Smith {
409d2471e0SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
419d2471e0SBarry Smith   PetscErrorCode ierr;
429d2471e0SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
439d2471e0SBarry Smith   PetscReal      fshift;
449d2471e0SBarry Smith   PetscBool      set,sym;
459d2471e0SBarry Smith 
469d2471e0SBarry Smith   PetscFunctionBegin;
479d2471e0SBarry Smith   ierr = MatIsSymmetricKnown(pc->pmat,&set,&sym);CHKERRQ(ierr);
489d2471e0SBarry 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");
499d2471e0SBarry Smith   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
509d2471e0SBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
519d2471e0SBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
529d2471e0SBarry Smith   PetscFunctionReturn(0);
539d2471e0SBarry Smith }
549d2471e0SBarry Smith 
55ace3abfcSBarry 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)
564b9ad928SBarry Smith {
574b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
58dfbe8321SBarry Smith   PetscErrorCode ierr;
597319c654SBarry Smith   MatSORType     stype = jac->sym;
603060034bSBarry Smith   PetscReal      fshift;
614b9ad928SBarry Smith 
624b9ad928SBarry Smith   PetscFunctionBegin;
63ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
642fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
65b0e188deSSatish Balay   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
66422a814eSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
67539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
684d0a8057SBarry Smith   *outits = its;
694d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
704b9ad928SBarry Smith   PetscFunctionReturn(0);
714b9ad928SBarry Smith }
724b9ad928SBarry Smith 
734416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
744b9ad928SBarry Smith {
754b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
76dfbe8321SBarry Smith   PetscErrorCode ierr;
77ace3abfcSBarry Smith   PetscBool      flg;
784b9ad928SBarry Smith 
794b9ad928SBarry Smith   PetscFunctionBegin;
80e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
8194ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
8294ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
8394ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
8494ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
85acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
864b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
87acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
884b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
89acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
90a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
91acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
924b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
93acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
944b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
95acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
964b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
974b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
984b9ad928SBarry Smith   PetscFunctionReturn(0);
994b9ad928SBarry Smith }
1004b9ad928SBarry Smith 
101dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
1024b9ad928SBarry Smith {
1034b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
1044b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
1052fc52814SBarry Smith   const char     *sortype;
106dfbe8321SBarry Smith   PetscErrorCode ierr;
107ace3abfcSBarry Smith   PetscBool      iascii;
1084b9ad928SBarry Smith 
1094b9ad928SBarry Smith   PetscFunctionBegin;
110251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11132077d6dSBarry Smith   if (iascii) {
112efd4aadfSBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  zero initial guess\n");CHKERRQ(ierr);}
1134b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1144b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1154b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
116db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1174b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1184b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
119db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1204b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1214b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1224b9ad928SBarry Smith     else                                                                     sortype = "unknown";
123efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
1244b9ad928SBarry Smith   }
1254b9ad928SBarry Smith   PetscFunctionReturn(0);
1264b9ad928SBarry Smith }
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
130f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1314b9ad928SBarry Smith {
132c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   PetscFunctionBegin;
1354b9ad928SBarry Smith   jac->sym = flag;
1364b9ad928SBarry Smith   PetscFunctionReturn(0);
1374b9ad928SBarry Smith }
1384b9ad928SBarry Smith 
139f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1404b9ad928SBarry Smith {
141c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith   PetscFunctionBegin;
144ce94432eSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1454b9ad928SBarry Smith   jac->omega = omega;
1464b9ad928SBarry Smith   PetscFunctionReturn(0);
1474b9ad928SBarry Smith }
1484b9ad928SBarry Smith 
149f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1504b9ad928SBarry Smith {
151c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1524b9ad928SBarry Smith 
1534b9ad928SBarry Smith   PetscFunctionBegin;
1544b9ad928SBarry Smith   jac->its  = its;
1554b9ad928SBarry Smith   jac->lits = lits;
1564b9ad928SBarry Smith   PetscFunctionReturn(0);
1574b9ad928SBarry Smith }
1584b9ad928SBarry Smith 
159c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
160c60c7ad4SBarry Smith {
161c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
162c60c7ad4SBarry Smith 
163c60c7ad4SBarry Smith   PetscFunctionBegin;
164c60c7ad4SBarry Smith   *flag = jac->sym;
165c60c7ad4SBarry Smith   PetscFunctionReturn(0);
166c60c7ad4SBarry Smith }
167c60c7ad4SBarry Smith 
168c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
169c60c7ad4SBarry Smith {
170c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
171c60c7ad4SBarry Smith 
172c60c7ad4SBarry Smith   PetscFunctionBegin;
173c60c7ad4SBarry Smith   *omega = jac->omega;
174c60c7ad4SBarry Smith   PetscFunctionReturn(0);
175c60c7ad4SBarry Smith }
176c60c7ad4SBarry Smith 
177c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
178c60c7ad4SBarry Smith {
179c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
180c60c7ad4SBarry Smith 
181c60c7ad4SBarry Smith   PetscFunctionBegin;
182c60c7ad4SBarry Smith   if (its)  *its = jac->its;
183c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
184c60c7ad4SBarry Smith   PetscFunctionReturn(0);
185c60c7ad4SBarry Smith }
186c60c7ad4SBarry Smith 
1874b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
188c60c7ad4SBarry Smith /*@
1891ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
190c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
191c60c7ad4SBarry Smith 
192c60c7ad4SBarry Smith    Logically Collective on PC
193c60c7ad4SBarry Smith 
194c60c7ad4SBarry Smith    Input Parameter:
195c60c7ad4SBarry Smith .  pc - the preconditioner context
196c60c7ad4SBarry Smith 
197c60c7ad4SBarry Smith    Output Parameter:
198c60c7ad4SBarry Smith .  flag - one of the following
199c60c7ad4SBarry Smith .vb
200c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
201c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
202c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
203c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
204c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
205c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
206c60c7ad4SBarry Smith .ve
207c60c7ad4SBarry Smith 
208c60c7ad4SBarry Smith    Options Database Keys:
209c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
210c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
211c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
212c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
213c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
214c60c7ad4SBarry Smith 
215c60c7ad4SBarry Smith    Notes:
216c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
217c60c7ad4SBarry Smith    which can be chosen with the option
218c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
219c60c7ad4SBarry Smith 
220c60c7ad4SBarry Smith    Level: intermediate
221c60c7ad4SBarry Smith 
222c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
223c60c7ad4SBarry Smith 
224c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
225c60c7ad4SBarry Smith @*/
226c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
227c60c7ad4SBarry Smith {
228c60c7ad4SBarry Smith   PetscErrorCode ierr;
229c60c7ad4SBarry Smith 
230c60c7ad4SBarry Smith   PetscFunctionBegin;
231c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
232c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
233c60c7ad4SBarry Smith   PetscFunctionReturn(0);
234c60c7ad4SBarry Smith }
235c60c7ad4SBarry Smith 
236c60c7ad4SBarry Smith /*@
237c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
238c60c7ad4SBarry Smith    (where omega = 1.0 by default).
239c60c7ad4SBarry Smith 
240c60c7ad4SBarry Smith    Logically Collective on PC
241c60c7ad4SBarry Smith 
242c60c7ad4SBarry Smith    Input Parameter:
243c60c7ad4SBarry Smith .  pc - the preconditioner context
244c60c7ad4SBarry Smith 
245c60c7ad4SBarry Smith    Output Parameter:
246c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
247c60c7ad4SBarry Smith 
248c60c7ad4SBarry Smith    Options Database Key:
249c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
250c60c7ad4SBarry Smith 
251c60c7ad4SBarry Smith    Level: intermediate
252c60c7ad4SBarry Smith 
253c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
254c60c7ad4SBarry Smith 
255c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
256c60c7ad4SBarry Smith @*/
257c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
258c60c7ad4SBarry Smith {
259c60c7ad4SBarry Smith   PetscErrorCode ierr;
260c60c7ad4SBarry Smith 
261c60c7ad4SBarry Smith   PetscFunctionBegin;
262c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
263c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
264c60c7ad4SBarry Smith   PetscFunctionReturn(0);
265c60c7ad4SBarry Smith }
266c60c7ad4SBarry Smith 
267c60c7ad4SBarry Smith /*@
268c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
269c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
270c60c7ad4SBarry Smith 
271c60c7ad4SBarry Smith    Logically Collective on PC
272c60c7ad4SBarry Smith 
273c60c7ad4SBarry Smith    Input Parameter:
274c60c7ad4SBarry Smith .  pc - the preconditioner context
275c60c7ad4SBarry Smith 
276c60c7ad4SBarry Smith    Output Parameter:
277c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
278c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
279c60c7ad4SBarry Smith 
280c60c7ad4SBarry Smith    Options Database Key:
281c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
282c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
283c60c7ad4SBarry Smith 
284c60c7ad4SBarry Smith    Level: intermediate
285c60c7ad4SBarry Smith 
286*95452b02SPatrick Sanan    Notes:
287*95452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
288c60c7ad4SBarry Smith 
289c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
290c60c7ad4SBarry Smith 
291c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
292c60c7ad4SBarry Smith @*/
293c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
294c60c7ad4SBarry Smith {
295c60c7ad4SBarry Smith   PetscErrorCode ierr;
296c60c7ad4SBarry Smith 
297c60c7ad4SBarry Smith   PetscFunctionBegin;
298c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
299c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
300c60c7ad4SBarry Smith   PetscFunctionReturn(0);
301c60c7ad4SBarry Smith }
302c60c7ad4SBarry Smith 
3034b9ad928SBarry Smith /*@
3044b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
3054b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
3064b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
3074b9ad928SBarry Smith 
3083f9fe445SBarry Smith    Logically Collective on PC
3094b9ad928SBarry Smith 
3104b9ad928SBarry Smith    Input Parameters:
3114b9ad928SBarry Smith +  pc - the preconditioner context
3124b9ad928SBarry Smith -  flag - one of the following
3134b9ad928SBarry Smith .vb
3144b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3154b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3164b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3174b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3184b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3194b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3204b9ad928SBarry Smith .ve
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Options Database Keys:
3234b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3244b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3254b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3264b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3274b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith    Notes:
3304b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3314b9ad928SBarry Smith    which can be chosen with the option
3324b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3334b9ad928SBarry Smith 
3344b9ad928SBarry Smith    Level: intermediate
3354b9ad928SBarry Smith 
3364b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3394b9ad928SBarry Smith @*/
3407087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3414b9ad928SBarry Smith {
3424ac538c5SBarry Smith   PetscErrorCode ierr;
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith   PetscFunctionBegin;
3450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
346c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3474ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3484b9ad928SBarry Smith   PetscFunctionReturn(0);
3494b9ad928SBarry Smith }
3504b9ad928SBarry Smith 
3514b9ad928SBarry Smith /*@
3524b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3534b9ad928SBarry Smith    (where omega = 1.0 by default).
3544b9ad928SBarry Smith 
3553f9fe445SBarry Smith    Logically Collective on PC
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith    Input Parameters:
3584b9ad928SBarry Smith +  pc - the preconditioner context
3594b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3604b9ad928SBarry Smith 
3614b9ad928SBarry Smith    Options Database Key:
3624b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3634b9ad928SBarry Smith 
3644b9ad928SBarry Smith    Level: intermediate
3654b9ad928SBarry Smith 
3664b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
3694b9ad928SBarry Smith @*/
3707087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3714b9ad928SBarry Smith {
3724ac538c5SBarry Smith   PetscErrorCode ierr;
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith   PetscFunctionBegin;
375c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
376c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3774ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
3784b9ad928SBarry Smith   PetscFunctionReturn(0);
3794b9ad928SBarry Smith }
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith /*@
3824b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3834b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3844b9ad928SBarry Smith 
3853f9fe445SBarry Smith    Logically Collective on PC
3864b9ad928SBarry Smith 
3874b9ad928SBarry Smith    Input Parameters:
3884b9ad928SBarry Smith +  pc - the preconditioner context
3894b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
3904b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Options Database Key:
3934b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
3944b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith    Level: intermediate
3974b9ad928SBarry Smith 
398*95452b02SPatrick Sanan    Notes:
399*95452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
4044b9ad928SBarry Smith @*/
4057087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
4064b9ad928SBarry Smith {
4074ac538c5SBarry Smith   PetscErrorCode ierr;
4084b9ad928SBarry Smith 
4094b9ad928SBarry Smith   PetscFunctionBegin;
4100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
411c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
4124ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
4134b9ad928SBarry Smith   PetscFunctionReturn(0);
4144b9ad928SBarry Smith }
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith /*MC
4174b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith    Options Database Keys:
4204b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4214b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
422a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4234b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
424a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4254b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4264b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
427422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
428a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
429a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith    Level: beginner
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
4344b9ad928SBarry Smith 
435*95452b02SPatrick Sanan    Notes:
436*95452b02SPatrick Sanan     Only implemented for the AIJ  and SeqBAIJ matrix formats.
4374b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4384b9ad928SBarry Smith           Jacobi with SOR on each block.
4394b9ad928SBarry Smith 
440422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
441422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
442422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
443422a814eSBarry Smith           zero pivot.
444422a814eSBarry Smith 
44537a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
44637a17b4dSBarry Smith 
447422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
448422a814eSBarry Smith           the computation is stopped with an error
449422a814eSBarry Smith 
45067991ba4SBarry Smith           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
45167991ba4SBarry Smith           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
45267991ba4SBarry Smith 
4534b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4544b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
4554b9ad928SBarry Smith M*/
4564b9ad928SBarry Smith 
4578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4584b9ad928SBarry Smith {
459dfbe8321SBarry Smith   PetscErrorCode ierr;
4604b9ad928SBarry Smith   PC_SOR         *jac;
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith   PetscFunctionBegin;
463b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4669d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4674b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4684b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4694b9ad928SBarry Smith   pc->ops->setup           = 0;
4704b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4714b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4724b9ad928SBarry Smith   pc->data                 = (void*)jac;
473d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4744b9ad928SBarry Smith   jac->omega               = 1.0;
47596fc60bcSBarry Smith   jac->fshift              = 0.0;
4764b9ad928SBarry Smith   jac->its                 = 1;
4774b9ad928SBarry Smith   jac->lits                = 1;
4784b9ad928SBarry Smith 
479bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
480bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
481bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
482c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
483c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
484c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
4854b9ad928SBarry Smith   PetscFunctionReturn(0);
4864b9ad928SBarry Smith }
4874b9ad928SBarry Smith 
4884b9ad928SBarry Smith 
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith 
4914b9ad928SBarry Smith 
492