xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    Defines a  (S)SOR  preconditioner for any Mat implementation
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcimpl.h>               /*I "petscpc.h" I*/
54b9ad928SBarry Smith 
64b9ad928SBarry Smith typedef struct {
7c1ac3661SBarry Smith   PetscInt   its;         /* inner iterations, number of sweeps */
8c1ac3661SBarry Smith   PetscInt   lits;        /* local inner iterations, number of sweeps applied by the local matrix mat->A */
94b9ad928SBarry Smith   MatSORType sym;         /* forward, reverse, symmetric etc. */
104b9ad928SBarry Smith   PetscReal  omega;
1129c1d7e0SHong Zhang   PetscReal  fshift;
124b9ad928SBarry Smith } PC_SOR;
134b9ad928SBarry Smith 
146849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc)
154b9ad928SBarry Smith {
16dfbe8321SBarry Smith   PetscErrorCode ierr;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
19c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
204b9ad928SBarry Smith   PetscFunctionReturn(0);
214b9ad928SBarry Smith }
224b9ad928SBarry Smith 
236849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
244b9ad928SBarry Smith {
254b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
26dfbe8321SBarry Smith   PetscErrorCode ierr;
27c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
284b9ad928SBarry Smith 
294b9ad928SBarry Smith   PetscFunctionBegin;
301230317dSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
31539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
324b9ad928SBarry Smith   PetscFunctionReturn(0);
334b9ad928SBarry Smith }
344b9ad928SBarry Smith 
359d2471e0SBarry Smith static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
369d2471e0SBarry Smith {
379d2471e0SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
389d2471e0SBarry Smith   PetscErrorCode ierr;
399d2471e0SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
409d2471e0SBarry Smith   PetscBool      set,sym;
419d2471e0SBarry Smith 
429d2471e0SBarry Smith   PetscFunctionBegin;
439d2471e0SBarry Smith   ierr = MatIsSymmetricKnown(pc->pmat,&set,&sym);CHKERRQ(ierr);
44*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!set || !sym || (jac->sym != SOR_SYMMETRIC_SWEEP && jac->sym != SOR_LOCAL_SYMMETRIC_SWEEP),PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Can only apply transpose of SOR if matrix is symmetric and sweep is symmetric");
451230317dSBarry Smith   ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr);
469d2471e0SBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
479d2471e0SBarry Smith   PetscFunctionReturn(0);
489d2471e0SBarry Smith }
499d2471e0SBarry Smith 
50ace3abfcSBarry 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)
514b9ad928SBarry Smith {
524b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
547319c654SBarry Smith   MatSORType     stype = jac->sym;
554b9ad928SBarry Smith 
564b9ad928SBarry Smith   PetscFunctionBegin;
577d3de750SJacob Faibussowitsch   ierr = PetscInfo(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr);
582fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
59be1e2188SBarry Smith   ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr);
60539c167fSBarry Smith   ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr);
614d0a8057SBarry Smith   *outits = its;
624d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith 
664416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
674b9ad928SBarry Smith {
684b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
70ace3abfcSBarry Smith   PetscBool      flg;
714b9ad928SBarry Smith 
724b9ad928SBarry Smith   PetscFunctionBegin;
73e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr);
7494ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr);
7594ae4db5SBarry Smith   ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr);
7694ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr);
7794ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr);
78acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
794b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
80acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
814b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);}
82acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
83a9510f2eSBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);}
84acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
854b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);}
86acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
874b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);}
88acfcf0e5SJed Brown   ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr);
894b9ad928SBarry Smith   if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);}
904b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
914b9ad928SBarry Smith   PetscFunctionReturn(0);
924b9ad928SBarry Smith }
934b9ad928SBarry Smith 
94dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
954b9ad928SBarry Smith {
964b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
974b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
982fc52814SBarry Smith   const char     *sortype;
99dfbe8321SBarry Smith   PetscErrorCode ierr;
100ace3abfcSBarry Smith   PetscBool      iascii;
1014b9ad928SBarry Smith 
1024b9ad928SBarry Smith   PetscFunctionBegin;
103251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10432077d6dSBarry Smith   if (iascii) {
105efd4aadfSBarry Smith     if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer,"  zero initial guess\n");CHKERRQ(ierr);}
1064b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1074b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1084b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
109db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1104b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1114b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
112db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1134b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1144b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1154b9ad928SBarry Smith     else                                                                     sortype = "unknown";
116efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr);
1174b9ad928SBarry Smith   }
1184b9ad928SBarry Smith   PetscFunctionReturn(0);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
122f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1234b9ad928SBarry Smith {
124c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith   PetscFunctionBegin;
1274b9ad928SBarry Smith   jac->sym = flag;
1284b9ad928SBarry Smith   PetscFunctionReturn(0);
1294b9ad928SBarry Smith }
1304b9ad928SBarry Smith 
131f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1324b9ad928SBarry Smith {
133c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith   PetscFunctionBegin;
136*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(omega >= 2.0 || omega <= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1374b9ad928SBarry Smith   jac->omega = omega;
1384b9ad928SBarry Smith   PetscFunctionReturn(0);
1394b9ad928SBarry Smith }
1404b9ad928SBarry Smith 
141f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1424b9ad928SBarry Smith {
143c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1444b9ad928SBarry Smith 
1454b9ad928SBarry Smith   PetscFunctionBegin;
1464b9ad928SBarry Smith   jac->its  = its;
1474b9ad928SBarry Smith   jac->lits = lits;
1484b9ad928SBarry Smith   PetscFunctionReturn(0);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith 
151c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
152c60c7ad4SBarry Smith {
153c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
154c60c7ad4SBarry Smith 
155c60c7ad4SBarry Smith   PetscFunctionBegin;
156c60c7ad4SBarry Smith   *flag = jac->sym;
157c60c7ad4SBarry Smith   PetscFunctionReturn(0);
158c60c7ad4SBarry Smith }
159c60c7ad4SBarry Smith 
160c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
161c60c7ad4SBarry Smith {
162c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
163c60c7ad4SBarry Smith 
164c60c7ad4SBarry Smith   PetscFunctionBegin;
165c60c7ad4SBarry Smith   *omega = jac->omega;
166c60c7ad4SBarry Smith   PetscFunctionReturn(0);
167c60c7ad4SBarry Smith }
168c60c7ad4SBarry Smith 
169c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
170c60c7ad4SBarry Smith {
171c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
172c60c7ad4SBarry Smith 
173c60c7ad4SBarry Smith   PetscFunctionBegin;
174c60c7ad4SBarry Smith   if (its)  *its = jac->its;
175c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
176c60c7ad4SBarry Smith   PetscFunctionReturn(0);
177c60c7ad4SBarry Smith }
178c60c7ad4SBarry Smith 
1794b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
180c60c7ad4SBarry Smith /*@
1811ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
182c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
183c60c7ad4SBarry Smith 
184c60c7ad4SBarry Smith    Logically Collective on PC
185c60c7ad4SBarry Smith 
186c60c7ad4SBarry Smith    Input Parameter:
187c60c7ad4SBarry Smith .  pc - the preconditioner context
188c60c7ad4SBarry Smith 
189c60c7ad4SBarry Smith    Output Parameter:
190c60c7ad4SBarry Smith .  flag - one of the following
191c60c7ad4SBarry Smith .vb
192c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
193c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
194c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
195c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
196c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
197c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
198c60c7ad4SBarry Smith .ve
199c60c7ad4SBarry Smith 
200c60c7ad4SBarry Smith    Options Database Keys:
201c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
202c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
203c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
204c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
205c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
206c60c7ad4SBarry Smith 
207c60c7ad4SBarry Smith    Notes:
208c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
209c60c7ad4SBarry Smith    which can be chosen with the option
210c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
211c60c7ad4SBarry Smith 
212c60c7ad4SBarry Smith    Level: intermediate
213c60c7ad4SBarry Smith 
214c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
215c60c7ad4SBarry Smith @*/
216c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
217c60c7ad4SBarry Smith {
218c60c7ad4SBarry Smith   PetscErrorCode ierr;
219c60c7ad4SBarry Smith 
220c60c7ad4SBarry Smith   PetscFunctionBegin;
221c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
222c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr);
223c60c7ad4SBarry Smith   PetscFunctionReturn(0);
224c60c7ad4SBarry Smith }
225c60c7ad4SBarry Smith 
226c60c7ad4SBarry Smith /*@
227c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
228c60c7ad4SBarry Smith    (where omega = 1.0 by default).
229c60c7ad4SBarry Smith 
230c60c7ad4SBarry Smith    Logically Collective on PC
231c60c7ad4SBarry Smith 
232c60c7ad4SBarry Smith    Input Parameter:
233c60c7ad4SBarry Smith .  pc - the preconditioner context
234c60c7ad4SBarry Smith 
235c60c7ad4SBarry Smith    Output Parameter:
236c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
237c60c7ad4SBarry Smith 
238c60c7ad4SBarry Smith    Options Database Key:
239c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
240c60c7ad4SBarry Smith 
241c60c7ad4SBarry Smith    Level: intermediate
242c60c7ad4SBarry Smith 
243c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
244c60c7ad4SBarry Smith @*/
245c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
246c60c7ad4SBarry Smith {
247c60c7ad4SBarry Smith   PetscErrorCode ierr;
248c60c7ad4SBarry Smith 
249c60c7ad4SBarry Smith   PetscFunctionBegin;
250c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
251c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
252c60c7ad4SBarry Smith   PetscFunctionReturn(0);
253c60c7ad4SBarry Smith }
254c60c7ad4SBarry Smith 
255c60c7ad4SBarry Smith /*@
256c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
257c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
258c60c7ad4SBarry Smith 
259c60c7ad4SBarry Smith    Logically Collective on PC
260c60c7ad4SBarry Smith 
261c60c7ad4SBarry Smith    Input Parameter:
262c60c7ad4SBarry Smith .  pc - the preconditioner context
263c60c7ad4SBarry Smith 
264d8d19677SJose E. Roman    Output Parameters:
265c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
266c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
267c60c7ad4SBarry Smith 
268c60c7ad4SBarry Smith    Options Database Key:
269c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
270c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
271c60c7ad4SBarry Smith 
272c60c7ad4SBarry Smith    Level: intermediate
273c60c7ad4SBarry Smith 
27495452b02SPatrick Sanan    Notes:
27595452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
276c60c7ad4SBarry Smith 
277c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
278c60c7ad4SBarry Smith @*/
279c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
280c60c7ad4SBarry Smith {
281c60c7ad4SBarry Smith   PetscErrorCode ierr;
282c60c7ad4SBarry Smith 
283c60c7ad4SBarry Smith   PetscFunctionBegin;
284c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
285c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr);
286c60c7ad4SBarry Smith   PetscFunctionReturn(0);
287c60c7ad4SBarry Smith }
288c60c7ad4SBarry Smith 
2894b9ad928SBarry Smith /*@
2904b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
2914b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
2924b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
2934b9ad928SBarry Smith 
2943f9fe445SBarry Smith    Logically Collective on PC
2954b9ad928SBarry Smith 
2964b9ad928SBarry Smith    Input Parameters:
2974b9ad928SBarry Smith +  pc - the preconditioner context
2984b9ad928SBarry Smith -  flag - one of the following
2994b9ad928SBarry Smith .vb
3004b9ad928SBarry Smith     SOR_FORWARD_SWEEP
3014b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
3024b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
3034b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
3044b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
3054b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
3064b9ad928SBarry Smith .ve
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith    Options Database Keys:
3094b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3104b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
3114b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
3124b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3134b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith    Notes:
3164b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3174b9ad928SBarry Smith    which can be chosen with the option
3184b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith    Level: intermediate
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3234b9ad928SBarry Smith @*/
3247087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3254b9ad928SBarry Smith {
3264ac538c5SBarry Smith   PetscErrorCode ierr;
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith   PetscFunctionBegin;
3290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
330c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3314ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr);
3324b9ad928SBarry Smith   PetscFunctionReturn(0);
3334b9ad928SBarry Smith }
3344b9ad928SBarry Smith 
3354b9ad928SBarry Smith /*@
3364b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3374b9ad928SBarry Smith    (where omega = 1.0 by default).
3384b9ad928SBarry Smith 
3393f9fe445SBarry Smith    Logically Collective on PC
3404b9ad928SBarry Smith 
3414b9ad928SBarry Smith    Input Parameters:
3424b9ad928SBarry Smith +  pc - the preconditioner context
3434b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith    Options Database Key:
3464b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith    Level: intermediate
3494b9ad928SBarry Smith 
350485168efSMatthew Knepley    Note:
351485168efSMatthew Knepley    If omega != 1, you will need to set the MAT_USE_INODES option to PETSC_FALSE on the matrix.
352485168efSMatthew Knepley 
353485168efSMatthew Knepley .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), MatSetOption()
3544b9ad928SBarry Smith @*/
3557087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3564b9ad928SBarry Smith {
3574ac538c5SBarry Smith   PetscErrorCode ierr;
3584b9ad928SBarry Smith 
3594b9ad928SBarry Smith   PetscFunctionBegin;
360c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
361c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3624ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
3634b9ad928SBarry Smith   PetscFunctionReturn(0);
3644b9ad928SBarry Smith }
3654b9ad928SBarry Smith 
3664b9ad928SBarry Smith /*@
3674b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3684b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3694b9ad928SBarry Smith 
3703f9fe445SBarry Smith    Logically Collective on PC
3714b9ad928SBarry Smith 
3724b9ad928SBarry Smith    Input Parameters:
3734b9ad928SBarry Smith +  pc - the preconditioner context
3744b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
3754b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith    Options Database Key:
3784b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
3794b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith    Level: intermediate
3824b9ad928SBarry Smith 
38395452b02SPatrick Sanan    Notes:
38495452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
3854b9ad928SBarry Smith 
3864b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
3874b9ad928SBarry Smith @*/
3887087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
3894b9ad928SBarry Smith {
3904ac538c5SBarry Smith   PetscErrorCode ierr;
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith   PetscFunctionBegin;
3930700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
394c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
3954ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr);
3964b9ad928SBarry Smith   PetscFunctionReturn(0);
3974b9ad928SBarry Smith }
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith /*MC
4004b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith    Options Database Keys:
4034b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
4044b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
405a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
4064b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
407a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
4084b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
4094b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
410422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
411a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
412a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith    Level: beginner
4154b9ad928SBarry Smith 
41695452b02SPatrick Sanan    Notes:
41795452b02SPatrick Sanan     Only implemented for the AIJ  and SeqBAIJ matrix formats.
4184b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4194b9ad928SBarry Smith           Jacobi with SOR on each block.
4204b9ad928SBarry Smith 
421422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
422422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
423422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
424422a814eSBarry Smith           zero pivot.
425422a814eSBarry Smith 
42637a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
42737a17b4dSBarry Smith 
428422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
429422a814eSBarry Smith           the computation is stopped with an error
430422a814eSBarry Smith 
43167991ba4SBarry Smith           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
43267991ba4SBarry Smith           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
43367991ba4SBarry Smith 
434485168efSMatthew Knepley           If omega != 1, you will need to set the MAT_USE_INODES option to PETSC_FALSE on the matrix.
435485168efSMatthew Knepley 
4364b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
437485168efSMatthew Knepley            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT, MatSetOption()
4384b9ad928SBarry Smith M*/
4394b9ad928SBarry Smith 
4408cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4414b9ad928SBarry Smith {
442dfbe8321SBarry Smith   PetscErrorCode ierr;
4434b9ad928SBarry Smith   PC_SOR         *jac;
4444b9ad928SBarry Smith 
4454b9ad928SBarry Smith   PetscFunctionBegin;
446b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
4474b9ad928SBarry Smith 
4484b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4499d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4504b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4514b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4520a545947SLisandro Dalcin   pc->ops->setup           = NULL;
4534b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4544b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4554b9ad928SBarry Smith   pc->data                 = (void*)jac;
456d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4574b9ad928SBarry Smith   jac->omega               = 1.0;
45896fc60bcSBarry Smith   jac->fshift              = 0.0;
4594b9ad928SBarry Smith   jac->its                 = 1;
4604b9ad928SBarry Smith   jac->lits                = 1;
4614b9ad928SBarry Smith 
462bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr);
463bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr);
464bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr);
465c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr);
466c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr);
467c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr);
4684b9ad928SBarry Smith   PetscFunctionReturn(0);
4694b9ad928SBarry Smith }
470