xref: /petsc/src/ksp/pc/impls/sor/sor.c (revision 7827d75ba736e00c19d105c058b1c2ddcca945c7)
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 {
164b9ad928SBarry Smith   PetscFunctionBegin;
179566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
184b9ad928SBarry Smith   PetscFunctionReturn(0);
194b9ad928SBarry Smith }
204b9ad928SBarry Smith 
216849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y)
224b9ad928SBarry Smith {
234b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
24c1ac3661SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
254b9ad928SBarry Smith 
264b9ad928SBarry Smith   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y));
289566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason));
294b9ad928SBarry Smith   PetscFunctionReturn(0);
304b9ad928SBarry Smith }
314b9ad928SBarry Smith 
329d2471e0SBarry Smith static PetscErrorCode PCApplyTranspose_SOR(PC pc,Vec x,Vec y)
339d2471e0SBarry Smith {
349d2471e0SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
359d2471e0SBarry Smith   PetscInt       flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
369d2471e0SBarry Smith   PetscBool      set,sym;
379d2471e0SBarry Smith 
389d2471e0SBarry Smith   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetricKnown(pc->pmat,&set,&sym));
40*7827d75bSBarry Smith   PetscCheck(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");
419566063dSJacob Faibussowitsch   PetscCall(MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y));
429566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason));
439d2471e0SBarry Smith   PetscFunctionReturn(0);
449d2471e0SBarry Smith }
459d2471e0SBarry Smith 
46ace3abfcSBarry 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)
474b9ad928SBarry Smith {
484b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
497319c654SBarry Smith   MatSORType     stype = jac->sym;
504b9ad928SBarry Smith 
514b9ad928SBarry Smith   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"Warning, convergence critera ignored, using %D iterations\n",its));
532fa5cd67SKarl Rupp   if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS);
549566063dSJacob Faibussowitsch   PetscCall(MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y));
559566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason));
564d0a8057SBarry Smith   *outits = its;
574d0a8057SBarry Smith   *reason = PCRICHARDSON_CONVERGED_ITS;
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc)
624b9ad928SBarry Smith {
634b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
64ace3abfcSBarry Smith   PetscBool      flg;
654b9ad928SBarry Smith 
664b9ad928SBarry Smith   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"(S)SOR options"));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg));
739566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg));
759566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg));
779566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP));
789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg));
799566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP));
809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg));
819566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP));
829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg));
839566063dSJacob Faibussowitsch   if (flg) PetscCall(PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP));
849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
854b9ad928SBarry Smith   PetscFunctionReturn(0);
864b9ad928SBarry Smith }
874b9ad928SBarry Smith 
88dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer)
894b9ad928SBarry Smith {
904b9ad928SBarry Smith   PC_SOR         *jac = (PC_SOR*)pc->data;
914b9ad928SBarry Smith   MatSORType     sym  = jac->sym;
922fc52814SBarry Smith   const char     *sortype;
93ace3abfcSBarry Smith   PetscBool      iascii;
944b9ad928SBarry Smith 
954b9ad928SBarry Smith   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
9732077d6dSBarry Smith   if (iascii) {
989566063dSJacob Faibussowitsch     if (sym & SOR_ZERO_INITIAL_GUESS) PetscCall(PetscViewerASCIIPrintf(viewer,"  zero initial guess\n"));
994b9ad928SBarry Smith     if (sym == SOR_APPLY_UPPER)                                              sortype = "apply_upper";
1004b9ad928SBarry Smith     else if (sym == SOR_APPLY_LOWER)                                         sortype = "apply_lower";
1014b9ad928SBarry Smith     else if (sym & SOR_EISENSTAT)                                            sortype = "Eisenstat";
102db4deed7SKarl Rupp     else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP)             sortype = "symmetric";
1034b9ad928SBarry Smith     else if (sym & SOR_BACKWARD_SWEEP)                                       sortype = "backward";
1044b9ad928SBarry Smith     else if (sym & SOR_FORWARD_SWEEP)                                        sortype = "forward";
105db4deed7SKarl Rupp     else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric";
1064b9ad928SBarry Smith     else if (sym & SOR_LOCAL_FORWARD_SWEEP)                                  sortype = "local_forward";
1074b9ad928SBarry Smith     else if (sym & SOR_LOCAL_BACKWARD_SWEEP)                                 sortype = "local_backward";
1084b9ad928SBarry Smith     else                                                                     sortype = "unknown";
1099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega));
1104b9ad928SBarry Smith   }
1114b9ad928SBarry Smith   PetscFunctionReturn(0);
1124b9ad928SBarry Smith }
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
115f7a08781SBarry Smith static PetscErrorCode  PCSORSetSymmetric_SOR(PC pc,MatSORType flag)
1164b9ad928SBarry Smith {
117c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1184b9ad928SBarry Smith 
1194b9ad928SBarry Smith   PetscFunctionBegin;
1204b9ad928SBarry Smith   jac->sym = flag;
1214b9ad928SBarry Smith   PetscFunctionReturn(0);
1224b9ad928SBarry Smith }
1234b9ad928SBarry Smith 
124f7a08781SBarry Smith static PetscErrorCode  PCSORSetOmega_SOR(PC pc,PetscReal omega)
1254b9ad928SBarry Smith {
126c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   PetscFunctionBegin;
1292c71b3e2SJacob Faibussowitsch   PetscCheckFalse(omega >= 2.0 || omega <= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1304b9ad928SBarry Smith   jac->omega = omega;
1314b9ad928SBarry Smith   PetscFunctionReturn(0);
1324b9ad928SBarry Smith }
1334b9ad928SBarry Smith 
134f7a08781SBarry Smith static PetscErrorCode  PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)
1354b9ad928SBarry Smith {
136c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
1374b9ad928SBarry Smith 
1384b9ad928SBarry Smith   PetscFunctionBegin;
1394b9ad928SBarry Smith   jac->its  = its;
1404b9ad928SBarry Smith   jac->lits = lits;
1414b9ad928SBarry Smith   PetscFunctionReturn(0);
1424b9ad928SBarry Smith }
1434b9ad928SBarry Smith 
144c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetSymmetric_SOR(PC pc,MatSORType *flag)
145c60c7ad4SBarry Smith {
146c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
147c60c7ad4SBarry Smith 
148c60c7ad4SBarry Smith   PetscFunctionBegin;
149c60c7ad4SBarry Smith   *flag = jac->sym;
150c60c7ad4SBarry Smith   PetscFunctionReturn(0);
151c60c7ad4SBarry Smith }
152c60c7ad4SBarry Smith 
153c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetOmega_SOR(PC pc,PetscReal *omega)
154c60c7ad4SBarry Smith {
155c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
156c60c7ad4SBarry Smith 
157c60c7ad4SBarry Smith   PetscFunctionBegin;
158c60c7ad4SBarry Smith   *omega = jac->omega;
159c60c7ad4SBarry Smith   PetscFunctionReturn(0);
160c60c7ad4SBarry Smith }
161c60c7ad4SBarry Smith 
162c60c7ad4SBarry Smith static PetscErrorCode  PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits)
163c60c7ad4SBarry Smith {
164c60c7ad4SBarry Smith   PC_SOR *jac = (PC_SOR*)pc->data;
165c60c7ad4SBarry Smith 
166c60c7ad4SBarry Smith   PetscFunctionBegin;
167c60c7ad4SBarry Smith   if (its)  *its = jac->its;
168c60c7ad4SBarry Smith   if (lits) *lits = jac->lits;
169c60c7ad4SBarry Smith   PetscFunctionReturn(0);
170c60c7ad4SBarry Smith }
171c60c7ad4SBarry Smith 
1724b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
173c60c7ad4SBarry Smith /*@
1741ff2113eSBarry Smith    PCSORGetSymmetric - Gets the form the SOR preconditioner is using;   backward, or forward relaxation.  The local variants perform SOR on
175c60c7ad4SBarry Smith    each processor.  By default forward relaxation is used.
176c60c7ad4SBarry Smith 
177c60c7ad4SBarry Smith    Logically Collective on PC
178c60c7ad4SBarry Smith 
179c60c7ad4SBarry Smith    Input Parameter:
180c60c7ad4SBarry Smith .  pc - the preconditioner context
181c60c7ad4SBarry Smith 
182c60c7ad4SBarry Smith    Output Parameter:
183c60c7ad4SBarry Smith .  flag - one of the following
184c60c7ad4SBarry Smith .vb
185c60c7ad4SBarry Smith     SOR_FORWARD_SWEEP
186c60c7ad4SBarry Smith     SOR_BACKWARD_SWEEP
187c60c7ad4SBarry Smith     SOR_SYMMETRIC_SWEEP
188c60c7ad4SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
189c60c7ad4SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
190c60c7ad4SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
191c60c7ad4SBarry Smith .ve
192c60c7ad4SBarry Smith 
193c60c7ad4SBarry Smith    Options Database Keys:
194c60c7ad4SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
195c60c7ad4SBarry Smith .  -pc_sor_backward - Activates backward version
196c60c7ad4SBarry Smith .  -pc_sor_local_forward - Activates local forward version
197c60c7ad4SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
198c60c7ad4SBarry Smith -  -pc_sor_local_backward - Activates local backward version
199c60c7ad4SBarry Smith 
200c60c7ad4SBarry Smith    Notes:
201c60c7ad4SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
202c60c7ad4SBarry Smith    which can be chosen with the option
203c60c7ad4SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
204c60c7ad4SBarry Smith 
205c60c7ad4SBarry Smith    Level: intermediate
206c60c7ad4SBarry Smith 
207c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric()
208c60c7ad4SBarry Smith @*/
209c60c7ad4SBarry Smith PetscErrorCode  PCSORGetSymmetric(PC pc,MatSORType *flag)
210c60c7ad4SBarry Smith {
211c60c7ad4SBarry Smith   PetscFunctionBegin;
212c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2139566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag)));
214c60c7ad4SBarry Smith   PetscFunctionReturn(0);
215c60c7ad4SBarry Smith }
216c60c7ad4SBarry Smith 
217c60c7ad4SBarry Smith /*@
218c60c7ad4SBarry Smith    PCSORGetOmega - Gets the SOR relaxation coefficient, omega
219c60c7ad4SBarry Smith    (where omega = 1.0 by default).
220c60c7ad4SBarry Smith 
221c60c7ad4SBarry Smith    Logically Collective on PC
222c60c7ad4SBarry Smith 
223c60c7ad4SBarry Smith    Input Parameter:
224c60c7ad4SBarry Smith .  pc - the preconditioner context
225c60c7ad4SBarry Smith 
226c60c7ad4SBarry Smith    Output Parameter:
227c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2).
228c60c7ad4SBarry Smith 
229c60c7ad4SBarry Smith    Options Database Key:
230c60c7ad4SBarry Smith .  -pc_sor_omega <omega> - Sets omega
231c60c7ad4SBarry Smith 
232c60c7ad4SBarry Smith    Level: intermediate
233c60c7ad4SBarry Smith 
234c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega()
235c60c7ad4SBarry Smith @*/
236c60c7ad4SBarry Smith PetscErrorCode  PCSORGetOmega(PC pc,PetscReal *omega)
237c60c7ad4SBarry Smith {
238c60c7ad4SBarry Smith   PetscFunctionBegin;
239c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2409566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega)));
241c60c7ad4SBarry Smith   PetscFunctionReturn(0);
242c60c7ad4SBarry Smith }
243c60c7ad4SBarry Smith 
244c60c7ad4SBarry Smith /*@
245c60c7ad4SBarry Smith    PCSORGetIterations - Gets the number of inner iterations to
246c60c7ad4SBarry Smith    be used by the SOR preconditioner. The default is 1.
247c60c7ad4SBarry Smith 
248c60c7ad4SBarry Smith    Logically Collective on PC
249c60c7ad4SBarry Smith 
250c60c7ad4SBarry Smith    Input Parameter:
251c60c7ad4SBarry Smith .  pc - the preconditioner context
252c60c7ad4SBarry Smith 
253d8d19677SJose E. Roman    Output Parameters:
254c60c7ad4SBarry Smith +  lits - number of local iterations, smoothings over just variables on processor
255c60c7ad4SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
256c60c7ad4SBarry Smith 
257c60c7ad4SBarry Smith    Options Database Key:
258c60c7ad4SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
259c60c7ad4SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
260c60c7ad4SBarry Smith 
261c60c7ad4SBarry Smith    Level: intermediate
262c60c7ad4SBarry Smith 
26395452b02SPatrick Sanan    Notes:
26495452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
265c60c7ad4SBarry Smith 
266c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations()
267c60c7ad4SBarry Smith @*/
268c60c7ad4SBarry Smith PetscErrorCode  PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits)
269c60c7ad4SBarry Smith {
270c60c7ad4SBarry Smith   PetscFunctionBegin;
271c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2729566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits)));
273c60c7ad4SBarry Smith   PetscFunctionReturn(0);
274c60c7ad4SBarry Smith }
275c60c7ad4SBarry Smith 
2764b9ad928SBarry Smith /*@
2774b9ad928SBarry Smith    PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
2784b9ad928SBarry Smith    backward, or forward relaxation.  The local variants perform SOR on
2794b9ad928SBarry Smith    each processor.  By default forward relaxation is used.
2804b9ad928SBarry Smith 
2813f9fe445SBarry Smith    Logically Collective on PC
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith    Input Parameters:
2844b9ad928SBarry Smith +  pc - the preconditioner context
2854b9ad928SBarry Smith -  flag - one of the following
2864b9ad928SBarry Smith .vb
2874b9ad928SBarry Smith     SOR_FORWARD_SWEEP
2884b9ad928SBarry Smith     SOR_BACKWARD_SWEEP
2894b9ad928SBarry Smith     SOR_SYMMETRIC_SWEEP
2904b9ad928SBarry Smith     SOR_LOCAL_FORWARD_SWEEP
2914b9ad928SBarry Smith     SOR_LOCAL_BACKWARD_SWEEP
2924b9ad928SBarry Smith     SOR_LOCAL_SYMMETRIC_SWEEP
2934b9ad928SBarry Smith .ve
2944b9ad928SBarry Smith 
2954b9ad928SBarry Smith    Options Database Keys:
2964b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
2974b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
2984b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
2994b9ad928SBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version
3004b9ad928SBarry Smith -  -pc_sor_local_backward - Activates local backward version
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith    Notes:
3034b9ad928SBarry Smith    To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3044b9ad928SBarry Smith    which can be chosen with the option
3054b9ad928SBarry Smith .  -pc_type eisenstat - Activates Eisenstat trick
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith    Level: intermediate
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega()
3104b9ad928SBarry Smith @*/
3117087cfbeSBarry Smith PetscErrorCode  PCSORSetSymmetric(PC pc,MatSORType flag)
3124b9ad928SBarry Smith {
3134b9ad928SBarry Smith   PetscFunctionBegin;
3140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
315c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,flag,2);
3169566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag)));
3174b9ad928SBarry Smith   PetscFunctionReturn(0);
3184b9ad928SBarry Smith }
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith /*@
3214b9ad928SBarry Smith    PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3224b9ad928SBarry Smith    (where omega = 1.0 by default).
3234b9ad928SBarry Smith 
3243f9fe445SBarry Smith    Logically Collective on PC
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith    Input Parameters:
3274b9ad928SBarry Smith +  pc - the preconditioner context
3284b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2).
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Options Database Key:
3314b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith    Level: intermediate
3344b9ad928SBarry Smith 
335485168efSMatthew Knepley    Note:
336485168efSMatthew Knepley    If omega != 1, you will need to set the MAT_USE_INODES option to PETSC_FALSE on the matrix.
337485168efSMatthew Knepley 
338485168efSMatthew Knepley .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), MatSetOption()
3394b9ad928SBarry Smith @*/
3407087cfbeSBarry Smith PetscErrorCode  PCSORSetOmega(PC pc,PetscReal omega)
3414b9ad928SBarry Smith {
3424b9ad928SBarry Smith   PetscFunctionBegin;
343c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
344c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
3459566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega)));
3464b9ad928SBarry Smith   PetscFunctionReturn(0);
3474b9ad928SBarry Smith }
3484b9ad928SBarry Smith 
3494b9ad928SBarry Smith /*@
3504b9ad928SBarry Smith    PCSORSetIterations - Sets the number of inner iterations to
3514b9ad928SBarry Smith    be used by the SOR preconditioner. The default is 1.
3524b9ad928SBarry Smith 
3533f9fe445SBarry Smith    Logically Collective on PC
3544b9ad928SBarry Smith 
3554b9ad928SBarry Smith    Input Parameters:
3564b9ad928SBarry Smith +  pc - the preconditioner context
3574b9ad928SBarry Smith .  lits - number of local iterations, smoothings over just variables on processor
3584b9ad928SBarry Smith -  its - number of parallel iterations to use; each parallel iteration has lits local iterations
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith    Options Database Key:
3614b9ad928SBarry Smith +  -pc_sor_its <its> - Sets number of iterations
3624b9ad928SBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations
3634b9ad928SBarry Smith 
3644b9ad928SBarry Smith    Level: intermediate
3654b9ad928SBarry Smith 
36695452b02SPatrick Sanan    Notes:
36795452b02SPatrick Sanan     When run on one processor the number of smoothings is lits*its
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric()
3704b9ad928SBarry Smith @*/
3717087cfbeSBarry Smith PetscErrorCode  PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)
3724b9ad928SBarry Smith {
3734b9ad928SBarry Smith   PetscFunctionBegin;
3740700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
375c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,its,2);
3769566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits)));
3774b9ad928SBarry Smith   PetscFunctionReturn(0);
3784b9ad928SBarry Smith }
3794b9ad928SBarry Smith 
3804b9ad928SBarry Smith /*MC
3814b9ad928SBarry Smith      PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith    Options Database Keys:
3844b9ad928SBarry Smith +  -pc_sor_symmetric - Activates symmetric version
3854b9ad928SBarry Smith .  -pc_sor_backward - Activates backward version
386a9510f2eSBarry Smith .  -pc_sor_forward - Activates forward version
3874b9ad928SBarry Smith .  -pc_sor_local_forward - Activates local forward version
388a9510f2eSBarry Smith .  -pc_sor_local_symmetric - Activates local symmetric version  (default version)
3894b9ad928SBarry Smith .  -pc_sor_local_backward - Activates local backward version
3904b9ad928SBarry Smith .  -pc_sor_omega <omega> - Sets omega
391422a814eSBarry Smith .  -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
392a9510f2eSBarry Smith .  -pc_sor_its <its> - Sets number of iterations   (default 1)
393a9510f2eSBarry Smith -  -pc_sor_lits <lits> - Sets number of local iterations  (default 1)
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith    Level: beginner
3964b9ad928SBarry Smith 
39795452b02SPatrick Sanan    Notes:
39895452b02SPatrick Sanan     Only implemented for the AIJ  and SeqBAIJ matrix formats.
3994b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4004b9ad928SBarry Smith           Jacobi with SOR on each block.
4014b9ad928SBarry Smith 
402422a814eSBarry Smith           For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
403422a814eSBarry Smith           zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option
404422a814eSBarry Smith           KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the
405422a814eSBarry Smith           zero pivot.
406422a814eSBarry Smith 
40737a17b4dSBarry Smith           For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported.
40837a17b4dSBarry Smith 
409422a814eSBarry Smith           For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
410422a814eSBarry Smith           the computation is stopped with an error
411422a814eSBarry Smith 
41267991ba4SBarry Smith           If used with KSPRICHARDSON and no monitors the convergence test is skipped to improve speed, thus it always iterates
41367991ba4SBarry Smith           the maximum number of iterations you've selected for KSP. It is usually used in this mode as a smoother for multigrid.
41467991ba4SBarry Smith 
415485168efSMatthew Knepley           If omega != 1, you will need to set the MAT_USE_INODES option to PETSC_FALSE on the matrix.
416485168efSMatthew Knepley 
4174b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
418485168efSMatthew Knepley            PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT, MatSetOption()
4194b9ad928SBarry Smith M*/
4204b9ad928SBarry Smith 
4218cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
4224b9ad928SBarry Smith {
4234b9ad928SBarry Smith   PC_SOR         *jac;
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith   PetscFunctionBegin;
4269566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&jac));
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith   pc->ops->apply           = PCApply_SOR;
4299d2471e0SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_SOR;
4304b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_SOR;
4314b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SOR;
4320a545947SLisandro Dalcin   pc->ops->setup           = NULL;
4334b9ad928SBarry Smith   pc->ops->view            = PCView_SOR;
4344b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_SOR;
4354b9ad928SBarry Smith   pc->data                 = (void*)jac;
436d9bc8e36SBarry Smith   jac->sym                 = SOR_LOCAL_SYMMETRIC_SWEEP;
4374b9ad928SBarry Smith   jac->omega               = 1.0;
43896fc60bcSBarry Smith   jac->fshift              = 0.0;
4394b9ad928SBarry Smith   jac->its                 = 1;
4404b9ad928SBarry Smith   jac->lits                = 1;
4414b9ad928SBarry Smith 
4429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR));
4439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR));
4449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR));
4459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR));
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR));
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR));
4484b9ad928SBarry Smith   PetscFunctionReturn(0);
4494b9ad928SBarry Smith }
450