xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision be1e21885f3a7f32653323ef855a5a8d07798310)
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;
311230317dSBarry 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");
461230317dSBarry 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;
564b9ad928SBarry Smith 
574b9ad928SBarry Smith   PetscFunctionBegin;
58ae15b995SBarry Smith   ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
592fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
60*be1e2188SBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
61539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
624d0a8057SBarry Smith   *outits = its;
634d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith 
674416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
684b9ad928SBarry Smith {
694b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
70dfbe8321SBarry Smith   PetscErrorCode ierr;
71ace3abfcSBarry Smith   PetscBool      flg;
724b9ad928SBarry Smith 
734b9ad928SBarry Smith   PetscFunctionBegin;
74e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
7594ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
7694ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
7794ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
7894ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
79acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
804b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
81acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
824b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
83acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
84a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
85acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
864b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
87acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
884b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
89acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
904b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
914b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
924b9ad928SBarry Smith   PetscFunctionReturn(0);
934b9ad928SBarry Smith }
944b9ad928SBarry Smith 
95dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
964b9ad928SBarry Smith {
974b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
984b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
992fc52814SBarry Smith   const char     *sortype;
100dfbe8321SBarry Smith   PetscErrorCode ierr;
101ace3abfcSBarry Smith   PetscBool      iascii;
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith   PetscFunctionBegin;
104251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10532077d6dSBarry Smith   if (iascii) {
106efd4aadfSBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  zero initial guess\n");CHKERRQ(ierr);}
1074b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1084b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1094b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
110db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1114b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1124b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
113db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1144b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1154b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1164b9ad928SBarry Smith     else                                                                     sortype = "unknown";
117efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
1184b9ad928SBarry Smith   }
1194b9ad928SBarry Smith   PetscFunctionReturn(0);
1204b9ad928SBarry Smith }
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
124f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1254b9ad928SBarry Smith {
126c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   PetscFunctionBegin;
1294b9ad928SBarry Smith   jac->sym = flag;
1304b9ad928SBarry Smith   PetscFunctionReturn(0);
1314b9ad928SBarry Smith }
1324b9ad928SBarry Smith 
133f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1344b9ad928SBarry Smith {
135c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1364b9ad928SBarry Smith 
1374b9ad928SBarry Smith   PetscFunctionBegin;
138ce94432eSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1394b9ad928SBarry Smith   jac->omega = omega;
1404b9ad928SBarry Smith   PetscFunctionReturn(0);
1414b9ad928SBarry Smith }
1424b9ad928SBarry Smith 
143f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1444b9ad928SBarry Smith {
145c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1464b9ad928SBarry Smith 
1474b9ad928SBarry Smith   PetscFunctionBegin;
1484b9ad928SBarry Smith   jac->its  = its;
1494b9ad928SBarry Smith   jac->lits = lits;
1504b9ad928SBarry Smith   PetscFunctionReturn(0);
1514b9ad928SBarry Smith }
1524b9ad928SBarry Smith 
153c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
154c60c7ad4SBarry Smith {
155c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
156c60c7ad4SBarry Smith 
157c60c7ad4SBarry Smith   PetscFunctionBegin;
158c60c7ad4SBarry Smith   *flag = jac->sym;
159c60c7ad4SBarry Smith   PetscFunctionReturn(0);
160c60c7ad4SBarry Smith }
161c60c7ad4SBarry Smith 
162c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
163c60c7ad4SBarry Smith {
164c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
165c60c7ad4SBarry Smith 
166c60c7ad4SBarry Smith   PetscFunctionBegin;
167c60c7ad4SBarry Smith   *omega = jac->omega;
168c60c7ad4SBarry Smith   PetscFunctionReturn(0);
169c60c7ad4SBarry Smith }
170c60c7ad4SBarry Smith 
171c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
172c60c7ad4SBarry Smith {
173c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
174c60c7ad4SBarry Smith 
175c60c7ad4SBarry Smith   PetscFunctionBegin;
176c60c7ad4SBarry Smith   if (its)  *its = jac->its;
177c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
178c60c7ad4SBarry Smith   PetscFunctionReturn(0);
179c60c7ad4SBarry Smith }
180c60c7ad4SBarry Smith 
1814b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
182c60c7ad4SBarry Smith /*@
1831ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
184c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
185c60c7ad4SBarry Smith 
186c60c7ad4SBarry Smith    Logically Collective on PC
187c60c7ad4SBarry Smith 
188c60c7ad4SBarry Smith    Input Parameter:
189c60c7ad4SBarry Smith .  pc - the preconditioner context
190c60c7ad4SBarry Smith 
191c60c7ad4SBarry Smith    Output Parameter:
192c60c7ad4SBarry Smith .  flag - one of the following
193c60c7ad4SBarry Smith .vb
194c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
195c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
196c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
197c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
198c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
199c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
200c60c7ad4SBarry Smith .ve
201c60c7ad4SBarry Smith 
202c60c7ad4SBarry Smith    Options Database Keys:
203c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
204c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
205c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
206c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
207c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
208c60c7ad4SBarry Smith 
209c60c7ad4SBarry Smith    Notes:
210c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
211c60c7ad4SBarry Smith    which can be chosen with the option
212c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
213c60c7ad4SBarry Smith 
214c60c7ad4SBarry Smith    Level: intermediate
215c60c7ad4SBarry Smith 
216c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
217c60c7ad4SBarry Smith 
218c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
219c60c7ad4SBarry Smith @*/
220c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
221c60c7ad4SBarry Smith {
222c60c7ad4SBarry Smith   PetscErrorCode ierr;
223c60c7ad4SBarry Smith 
224c60c7ad4SBarry Smith   PetscFunctionBegin;
225c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
226c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
227c60c7ad4SBarry Smith   PetscFunctionReturn(0);
228c60c7ad4SBarry Smith }
229c60c7ad4SBarry Smith 
230c60c7ad4SBarry Smith /*@
231c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
232c60c7ad4SBarry Smith    (where omega = 1.0 by default).
233c60c7ad4SBarry Smith 
234c60c7ad4SBarry Smith    Logically Collective on PC
235c60c7ad4SBarry Smith 
236c60c7ad4SBarry Smith    Input Parameter:
237c60c7ad4SBarry Smith .  pc - the preconditioner context
238c60c7ad4SBarry Smith 
239c60c7ad4SBarry Smith    Output Parameter:
240c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
241c60c7ad4SBarry Smith 
242c60c7ad4SBarry Smith    Options Database Key:
243c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
244c60c7ad4SBarry Smith 
245c60c7ad4SBarry Smith    Level: intermediate
246c60c7ad4SBarry Smith 
247c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
248c60c7ad4SBarry Smith 
249c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
250c60c7ad4SBarry Smith @*/
251c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
252c60c7ad4SBarry Smith {
253c60c7ad4SBarry Smith   PetscErrorCode ierr;
254c60c7ad4SBarry Smith 
255c60c7ad4SBarry Smith   PetscFunctionBegin;
256c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
257c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
258c60c7ad4SBarry Smith   PetscFunctionReturn(0);
259c60c7ad4SBarry Smith }
260c60c7ad4SBarry Smith 
261c60c7ad4SBarry Smith /*@
262c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
263c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
264c60c7ad4SBarry Smith 
265c60c7ad4SBarry Smith    Logically Collective on PC
266c60c7ad4SBarry Smith 
267c60c7ad4SBarry Smith    Input Parameter:
268c60c7ad4SBarry Smith .  pc - the preconditioner context
269c60c7ad4SBarry Smith 
270c60c7ad4SBarry Smith    Output Parameter:
271c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
272c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
273c60c7ad4SBarry Smith 
274c60c7ad4SBarry Smith    Options Database Key:
275c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
276c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
277c60c7ad4SBarry Smith 
278c60c7ad4SBarry Smith    Level: intermediate
279c60c7ad4SBarry Smith 
280c60c7ad4SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
281c60c7ad4SBarry Smith 
282c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
283c60c7ad4SBarry Smith 
284c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
285c60c7ad4SBarry Smith @*/
286c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
287c60c7ad4SBarry Smith {
288c60c7ad4SBarry Smith   PetscErrorCode ierr;
289c60c7ad4SBarry Smith 
290c60c7ad4SBarry Smith   PetscFunctionBegin;
291c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
292c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
293c60c7ad4SBarry Smith   PetscFunctionReturn(0);
294c60c7ad4SBarry Smith }
295c60c7ad4SBarry Smith 
2964b9ad928SBarry Smith /*@
2974b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
2984b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
2994b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
3004b9ad928SBarry Smith 
3013f9fe445SBarry Smith    Logically Collective on PC
3024b9ad928SBarry Smith 
3034b9ad928SBarry Smith    Input Parameters:
3044b9ad928SBarry Smith +  pc - the preconditioner context
3054b9ad928SBarry Smith -  flag - one of the following
3064b9ad928SBarry Smith .vb
3074b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3084b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3094b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3104b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3114b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3124b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3134b9ad928SBarry Smith .ve
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith    Options Database Keys:
3164b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3174b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3184b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3194b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3204b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Notes:
3234b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3244b9ad928SBarry Smith    which can be chosen with the option
3254b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3264b9ad928SBarry Smith 
3274b9ad928SBarry Smith    Level: intermediate
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3324b9ad928SBarry Smith @*/
3337087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3344b9ad928SBarry Smith {
3354ac538c5SBarry Smith   PetscErrorCode ierr;
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith   PetscFunctionBegin;
3380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
339c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3404ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3414b9ad928SBarry Smith   PetscFunctionReturn(0);
3424b9ad928SBarry Smith }
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith /*@
3454b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3464b9ad928SBarry Smith    (where omega = 1.0 by default).
3474b9ad928SBarry Smith 
3483f9fe445SBarry Smith    Logically Collective on PC
3494b9ad928SBarry Smith 
3504b9ad928SBarry Smith    Input Parameters:
3514b9ad928SBarry Smith +  pc - the preconditioner context
3524b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith    Options Database Key:
3554b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith    Level: intermediate
3584b9ad928SBarry Smith 
3594b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega
3604b9ad928SBarry Smith 
3614b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega()
3624b9ad928SBarry Smith @*/
3637087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3644b9ad928SBarry Smith {
3654ac538c5SBarry Smith   PetscErrorCode ierr;
3664b9ad928SBarry Smith 
3674b9ad928SBarry Smith   PetscFunctionBegin;
368c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
369c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3704ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
3714b9ad928SBarry Smith   PetscFunctionReturn(0);
3724b9ad928SBarry Smith }
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith /*@
3754b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3764b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3774b9ad928SBarry Smith 
3783f9fe445SBarry Smith    Logically Collective on PC
3794b9ad928SBarry Smith 
3804b9ad928SBarry Smith    Input Parameters:
3814b9ad928SBarry Smith +  pc - the preconditioner context
3824b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
3834b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
3844b9ad928SBarry Smith 
3854b9ad928SBarry Smith    Options Database Key:
3864b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
3874b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
3884b9ad928SBarry Smith 
3894b9ad928SBarry Smith    Level: intermediate
3904b9ad928SBarry Smith 
3914b9ad928SBarry Smith    Notes: When run on one processor the number of smoothings is lits*its
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
3964b9ad928SBarry Smith @*/
3977087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
3984b9ad928SBarry Smith {
3994ac538c5SBarry Smith   PetscErrorCode ierr;
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith   PetscFunctionBegin;
4020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
403c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
4044ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
4054b9ad928SBarry Smith   PetscFunctionReturn(0);
4064b9ad928SBarry Smith }
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith /*MC
4094b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith    Options Database Keys:
4124b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4134b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
414a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4154b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
416a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4174b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4184b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
419422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
420a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
421a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith    Level: beginner
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel
4264b9ad928SBarry Smith 
42737a17b4dSBarry Smith    Notes: Only implemented for the AIJ  and SeqBAIJ matrix formats.
4284b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4294b9ad928SBarry Smith           Jacobi with SOR on each block.
4304b9ad928SBarry Smith 
431422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
432422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
433422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
434422a814eSBarry Smith           zero pivot.
435422a814eSBarry Smith 
43637a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
43737a17b4dSBarry Smith 
438422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
439422a814eSBarry Smith           the computation is stopped with an error
440422a814eSBarry Smith 
44167991ba4SBarry Smith           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
44267991ba4SBarry Smith           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
44367991ba4SBarry Smith 
4444b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4454b9ad928SBarry Smith            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT
4464b9ad928SBarry Smith M*/
4474b9ad928SBarry Smith 
4488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4494b9ad928SBarry Smith {
450dfbe8321SBarry Smith   PetscErrorCode ierr;
4514b9ad928SBarry Smith   PC_SOR         *jac;
4524b9ad928SBarry Smith 
4534b9ad928SBarry Smith   PetscFunctionBegin;
454b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4579d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4584b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4594b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4604b9ad928SBarry Smith   pc->ops->setup           = 0;
4614b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4624b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4634b9ad928SBarry Smith   pc->data                 = (void*)jac;
464d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4654b9ad928SBarry Smith   jac->omega               = 1.0;
46696fc60bcSBarry Smith   jac->fshift              = 0.0;
4674b9ad928SBarry Smith   jac->its                 = 1;
4684b9ad928SBarry Smith   jac->lits                = 1;
4694b9ad928SBarry Smith 
470bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
471bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
472bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
473c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
474c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
475c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
4764b9ad928SBarry Smith   PetscFunctionReturn(0);
4774b9ad928SBarry Smith }
4784b9ad928SBarry Smith 
4794b9ad928SBarry Smith 
4804b9ad928SBarry Smith 
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith 
483