xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 1230317d47e50dc0bffca96e8ee76ec386405dd1)
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;
294b9ad928SBarry Smith 
304b9ad928SBarry Smith   PetscFunctionBegin;
31*1230317dSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
32539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
334b9ad928SBarry Smith   PetscFunctionReturn(0);
344b9ad928SBarry Smith }
354b9ad928SBarry Smith 
369d2471e0SBarry Smith static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
379d2471e0SBarry Smith {
389d2471e0SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
399d2471e0SBarry Smith   PetscErrorCode ierr;
409d2471e0SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
419d2471e0SBarry Smith   PetscBool      set,sym;
429d2471e0SBarry Smith 
439d2471e0SBarry Smith   PetscFunctionBegin;
449d2471e0SBarry Smith   ierr = MatIsSymmetricKnown(pc->pmat,&set,&sym);CHKERRQ(ierr);
459d2471e0SBarry 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");
46*1230317dSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
479d2471e0SBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
489d2471e0SBarry Smith   PetscFunctionReturn(0);
499d2471e0SBarry Smith }
509d2471e0SBarry Smith 
51ace3abfcSBarry 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)
524b9ad928SBarry Smith {
534b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
557319c654SBarry Smith   MatSORType     stype = jac->sym;
563060034bSBarry Smith   PetscReal      fshift;
574b9ad928SBarry Smith 
584b9ad928SBarry Smith   PetscFunctionBegin;
59ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
602fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
61b0e188deSSatish Balay   fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0);
62422a814eSBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
63539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
644d0a8057SBarry Smith   *outits = its;
654d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
664b9ad928SBarry Smith   PetscFunctionReturn(0);
674b9ad928SBarry Smith }
684b9ad928SBarry Smith 
694416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
704b9ad928SBarry Smith {
714b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
72dfbe8321SBarry Smith   PetscErrorCode ierr;
73ace3abfcSBarry Smith   PetscBool      flg;
744b9ad928SBarry Smith 
754b9ad928SBarry Smith   PetscFunctionBegin;
76e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
7794ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
7894ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
7994ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
8094ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
81acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
824b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
83acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
844b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
85acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
86a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
87acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
884b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
89acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
904b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
91acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
924b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
934b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
944b9ad928SBarry Smith   PetscFunctionReturn(0);
954b9ad928SBarry Smith }
964b9ad928SBarry Smith 
97dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
984b9ad928SBarry Smith {
994b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
1004b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
1012fc52814SBarry Smith   const char     *sortype;
102dfbe8321SBarry Smith   PetscErrorCode ierr;
103ace3abfcSBarry Smith   PetscBool      iascii;
1044b9ad928SBarry Smith 
1054b9ad928SBarry Smith   PetscFunctionBegin;
106251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10732077d6dSBarry Smith   if (iascii) {
108efd4aadfSBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  zero initial guess\n");CHKERRQ(ierr);}
1094b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1104b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1114b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
112db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1134b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1144b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
115db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1164b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1174b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1184b9ad928SBarry Smith     else                                                                     sortype = "unknown";
119efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
1204b9ad928SBarry Smith   }
1214b9ad928SBarry Smith   PetscFunctionReturn(0);
1224b9ad928SBarry Smith }
1234b9ad928SBarry Smith 
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
126f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1274b9ad928SBarry Smith {
128c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
1314b9ad928SBarry Smith   jac->sym = flag;
1324b9ad928SBarry Smith   PetscFunctionReturn(0);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
135f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1364b9ad928SBarry Smith {
137c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
140ce94432eSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1414b9ad928SBarry Smith   jac->omega = omega;
1424b9ad928SBarry Smith   PetscFunctionReturn(0);
1434b9ad928SBarry Smith }
1444b9ad928SBarry Smith 
145f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1464b9ad928SBarry Smith {
147c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1484b9ad928SBarry Smith 
1494b9ad928SBarry Smith   PetscFunctionBegin;
1504b9ad928SBarry Smith   jac->its  = its;
1514b9ad928SBarry Smith   jac->lits = lits;
1524b9ad928SBarry Smith   PetscFunctionReturn(0);
1534b9ad928SBarry Smith }
1544b9ad928SBarry Smith 
155c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
156c60c7ad4SBarry Smith {
157c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
158c60c7ad4SBarry Smith 
159c60c7ad4SBarry Smith   PetscFunctionBegin;
160c60c7ad4SBarry Smith   *flag = jac->sym;
161c60c7ad4SBarry Smith   PetscFunctionReturn(0);
162c60c7ad4SBarry Smith }
163c60c7ad4SBarry Smith 
164c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
165c60c7ad4SBarry Smith {
166c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
167c60c7ad4SBarry Smith 
168c60c7ad4SBarry Smith   PetscFunctionBegin;
169c60c7ad4SBarry Smith   *omega = jac->omega;
170c60c7ad4SBarry Smith   PetscFunctionReturn(0);
171c60c7ad4SBarry Smith }
172c60c7ad4SBarry Smith 
173c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
174c60c7ad4SBarry Smith {
175c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
176c60c7ad4SBarry Smith 
177c60c7ad4SBarry Smith   PetscFunctionBegin;
178c60c7ad4SBarry Smith   if (its)  *its = jac->its;
179c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
180c60c7ad4SBarry Smith   PetscFunctionReturn(0);
181c60c7ad4SBarry Smith }
182c60c7ad4SBarry Smith 
1834b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
184c60c7ad4SBarry Smith /*@
1851ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
186c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
187c60c7ad4SBarry Smith 
188c60c7ad4SBarry Smith    Logically Collective on PC
189c60c7ad4SBarry Smith 
190c60c7ad4SBarry Smith    Input Parameter:
191c60c7ad4SBarry Smith .  pc - the preconditioner context
192c60c7ad4SBarry Smith 
193c60c7ad4SBarry Smith    Output Parameter:
194c60c7ad4SBarry Smith .  flag - one of the following
195c60c7ad4SBarry Smith .vb
196c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
197c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
198c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
199c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
200c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
201c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
202c60c7ad4SBarry Smith .ve
203c60c7ad4SBarry Smith 
204c60c7ad4SBarry Smith    Options Database Keys:
205c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
206c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
207c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
208c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
209c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
210c60c7ad4SBarry Smith 
211c60c7ad4SBarry Smith    Notes:
212c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
213c60c7ad4SBarry Smith    which can be chosen with the option
214c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
215c60c7ad4SBarry Smith 
216c60c7ad4SBarry Smith    Level: intermediate
217c60c7ad4SBarry Smith 
218c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
219c60c7ad4SBarry Smith 
220c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
221c60c7ad4SBarry Smith @*/
222c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
223c60c7ad4SBarry Smith {
224c60c7ad4SBarry Smith   PetscErrorCode ierr;
225c60c7ad4SBarry Smith 
226c60c7ad4SBarry Smith   PetscFunctionBegin;
227c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
228c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
229c60c7ad4SBarry Smith   PetscFunctionReturn(0);
230c60c7ad4SBarry Smith }
231c60c7ad4SBarry Smith 
232c60c7ad4SBarry Smith /*@
233c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
234c60c7ad4SBarry Smith    (where omega = 1.0 by default).
235c60c7ad4SBarry Smith 
236c60c7ad4SBarry Smith    Logically Collective on PC
237c60c7ad4SBarry Smith 
238c60c7ad4SBarry Smith    Input Parameter:
239c60c7ad4SBarry Smith .  pc - the preconditioner context
240c60c7ad4SBarry Smith 
241c60c7ad4SBarry Smith    Output Parameter:
242c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
243c60c7ad4SBarry Smith 
244c60c7ad4SBarry Smith    Options Database Key:
245c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
246c60c7ad4SBarry Smith 
247c60c7ad4SBarry Smith    Level: intermediate
248c60c7ad4SBarry Smith 
249c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
250c60c7ad4SBarry Smith 
251c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
252c60c7ad4SBarry Smith @*/
253c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
254c60c7ad4SBarry Smith {
255c60c7ad4SBarry Smith   PetscErrorCode ierr;
256c60c7ad4SBarry Smith 
257c60c7ad4SBarry Smith   PetscFunctionBegin;
258c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
259c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
260c60c7ad4SBarry Smith   PetscFunctionReturn(0);
261c60c7ad4SBarry Smith }
262c60c7ad4SBarry Smith 
263c60c7ad4SBarry Smith /*@
264c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
265c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
266c60c7ad4SBarry Smith 
267c60c7ad4SBarry Smith    Logically Collective on PC
268c60c7ad4SBarry Smith 
269c60c7ad4SBarry Smith    Input Parameter:
270c60c7ad4SBarry Smith .  pc - the preconditioner context
271c60c7ad4SBarry Smith 
272c60c7ad4SBarry Smith    Output Parameter:
273c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
274c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
275c60c7ad4SBarry Smith 
276c60c7ad4SBarry Smith    Options Database Key:
277c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
278c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
279c60c7ad4SBarry Smith 
280c60c7ad4SBarry Smith    Level: intermediate
281c60c7ad4SBarry Smith 
282c60c7ad4SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
283c60c7ad4SBarry Smith 
284c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
285c60c7ad4SBarry Smith 
286c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
287c60c7ad4SBarry Smith @*/
288c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
289c60c7ad4SBarry Smith {
290c60c7ad4SBarry Smith   PetscErrorCode ierr;
291c60c7ad4SBarry Smith 
292c60c7ad4SBarry Smith   PetscFunctionBegin;
293c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
294c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
295c60c7ad4SBarry Smith   PetscFunctionReturn(0);
296c60c7ad4SBarry Smith }
297c60c7ad4SBarry Smith 
2984b9ad928SBarry Smith /*@
2994b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
3004b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
3014b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
3024b9ad928SBarry Smith 
3033f9fe445SBarry Smith    Logically Collective on PC
3044b9ad928SBarry Smith 
3054b9ad928SBarry Smith    Input Parameters:
3064b9ad928SBarry Smith +  pc - the preconditioner context
3074b9ad928SBarry Smith -  flag - one of the following
3084b9ad928SBarry Smith .vb
3094b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3104b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3114b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3124b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3134b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3144b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3154b9ad928SBarry Smith .ve
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith    Options Database Keys:
3184b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3194b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3204b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3214b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3224b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith    Notes:
3254b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3264b9ad928SBarry Smith    which can be chosen with the option
3274b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith    Level: intermediate
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3344b9ad928SBarry Smith @*/
3357087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3364b9ad928SBarry Smith {
3374ac538c5SBarry Smith   PetscErrorCode ierr;
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith   PetscFunctionBegin;
3400700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
341c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3424ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3434b9ad928SBarry Smith   PetscFunctionReturn(0);
3444b9ad928SBarry Smith }
3454b9ad928SBarry Smith 
3464b9ad928SBarry Smith /*@
3474b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3484b9ad928SBarry Smith    (where omega = 1.0 by default).
3494b9ad928SBarry Smith 
3503f9fe445SBarry Smith    Logically Collective on PC
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith    Input Parameters:
3534b9ad928SBarry Smith +  pc - the preconditioner context
3544b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith    Options Database Key:
3574b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3584b9ad928SBarry Smith 
3594b9ad928SBarry Smith    Level: intermediate
3604b9ad928SBarry Smith 
3614b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
3644b9ad928SBarry Smith @*/
3657087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3664b9ad928SBarry Smith {
3674ac538c5SBarry Smith   PetscErrorCode ierr;
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith   PetscFunctionBegin;
370c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
371c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3724ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
3734b9ad928SBarry Smith   PetscFunctionReturn(0);
3744b9ad928SBarry Smith }
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith /*@
3774b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3784b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3794b9ad928SBarry Smith 
3803f9fe445SBarry Smith    Logically Collective on PC
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith    Input Parameters:
3834b9ad928SBarry Smith +  pc - the preconditioner context
3844b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
3854b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
3864b9ad928SBarry Smith 
3874b9ad928SBarry Smith    Options Database Key:
3884b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
3894b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
3904b9ad928SBarry Smith 
3914b9ad928SBarry Smith    Level: intermediate
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
3984b9ad928SBarry Smith @*/
3997087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
4004b9ad928SBarry Smith {
4014ac538c5SBarry Smith   PetscErrorCode ierr;
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith   PetscFunctionBegin;
4040700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
405c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
4064ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
4074b9ad928SBarry Smith   PetscFunctionReturn(0);
4084b9ad928SBarry Smith }
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith /*MC
4114b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4124b9ad928SBarry Smith 
4134b9ad928SBarry Smith    Options Database Keys:
4144b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4154b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
416a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4174b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
418a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4194b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4204b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
421422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
422a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
423a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith    Level: beginner
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
4284b9ad928SBarry Smith 
42937a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
4304b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4314b9ad928SBarry Smith           Jacobi with SOR on each block.
4324b9ad928SBarry Smith 
433422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
434422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
435422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
436422a814eSBarry Smith           zero pivot.
437422a814eSBarry Smith 
43837a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
43937a17b4dSBarry Smith 
440422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
441422a814eSBarry Smith           the computation is stopped with an error
442422a814eSBarry Smith 
44367991ba4SBarry Smith           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
44467991ba4SBarry Smith           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
44567991ba4SBarry Smith 
4464b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4474b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
4484b9ad928SBarry Smith M*/
4494b9ad928SBarry Smith 
4508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4514b9ad928SBarry Smith {
452dfbe8321SBarry Smith   PetscErrorCode ierr;
4534b9ad928SBarry Smith   PC_SOR         *jac;
4544b9ad928SBarry Smith 
4554b9ad928SBarry Smith   PetscFunctionBegin;
456b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4574b9ad928SBarry Smith 
4584b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4599d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4604b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4614b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4624b9ad928SBarry Smith   pc->ops->setup           = 0;
4634b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4644b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4654b9ad928SBarry Smith   pc->data                 = (void*)jac;
466d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4674b9ad928SBarry Smith   jac->omega               = 1.0;
46896fc60bcSBarry Smith   jac->fshift              = 0.0;
4694b9ad928SBarry Smith   jac->its                 = 1;
4704b9ad928SBarry Smith   jac->lits                = 1;
4714b9ad928SBarry Smith 
472bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
473bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
474bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
475c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
476c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
477c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
4784b9ad928SBarry Smith   PetscFunctionReturn(0);
4794b9ad928SBarry Smith }
4804b9ad928SBarry Smith 
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith 
4834b9ad928SBarry Smith 
4844b9ad928SBarry Smith 
485