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
PCDestroy_SOR(PC pc)14d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SOR(PC pc)
15d71ae5a4SJacob Faibussowitsch {
164b9ad928SBarry Smith PetscFunctionBegin;
172e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetSymmetric_C", NULL));
182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetOmega_C", NULL));
192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetIterations_C", NULL));
202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetSymmetric_C", NULL));
212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetOmega_C", NULL));
222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetIterations_C", NULL));
239566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
254b9ad928SBarry Smith }
264b9ad928SBarry Smith
PCApply_SOR(PC pc,Vec x,Vec y)27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SOR(PC pc, Vec x, Vec y)
28d71ae5a4SJacob Faibussowitsch {
294b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
30c1ac3661SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
314b9ad928SBarry Smith
324b9ad928SBarry Smith PetscFunctionBegin;
339566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, x, jac->omega, (MatSORType)flag, jac->fshift, jac->its, jac->lits, y));
349566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason));
353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
364b9ad928SBarry Smith }
374b9ad928SBarry Smith
PCApplyTranspose_SOR(PC pc,Vec x,Vec y)38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_SOR(PC pc, Vec x, Vec y)
39d71ae5a4SJacob Faibussowitsch {
409d2471e0SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
419d2471e0SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS;
429d2471e0SBarry Smith PetscBool set, sym;
439d2471e0SBarry Smith
449d2471e0SBarry Smith PetscFunctionBegin;
459566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &sym));
467827d75bSBarry 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");
479566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, x, jac->omega, (MatSORType)flag, jac->fshift, jac->its, jac->lits, y));
489566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason));
493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
509d2471e0SBarry Smith }
519d2471e0SBarry Smith
PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt * outits,PCRichardsonConvergedReason * reason)52d71ae5a4SJacob Faibussowitsch 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)
53d71ae5a4SJacob Faibussowitsch {
544b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
557319c654SBarry Smith MatSORType stype = jac->sym;
564b9ad928SBarry Smith
574b9ad928SBarry Smith PetscFunctionBegin;
582fa5cd67SKarl Rupp if (guesszero) stype = (MatSORType)(stype | SOR_ZERO_INITIAL_GUESS);
599566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, b, jac->omega, stype, jac->fshift, its * jac->its, jac->lits, y));
609566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason));
614d0a8057SBarry Smith *outits = its;
624d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS;
633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith
PCSetFromOptions_SOR(PC pc,PetscOptionItems PetscOptionsObject)66ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_SOR(PC pc, PetscOptionItems PetscOptionsObject)
67d71ae5a4SJacob Faibussowitsch {
684b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
69ace3abfcSBarry Smith PetscBool flg;
70d8e4f26eSJose E. Roman PetscReal omega;
714b9ad928SBarry Smith
724b9ad928SBarry Smith PetscFunctionBegin;
73d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "(S)SOR options");
74d8e4f26eSJose E. Roman PetscCall(PetscOptionsReal("-pc_sor_omega", "relaxation factor (0 < omega < 2)", "PCSORSetOmega", jac->omega, &omega, &flg));
75d8e4f26eSJose E. Roman if (flg) PetscCall(PCSORSetOmega(pc, omega));
769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_sor_diagonal_shift", "Add to the diagonal entries", "", jac->fshift, &jac->fshift, NULL));
779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_its", "number of inner SOR iterations", "PCSORSetIterations", jac->its, &jac->its, NULL));
789566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_lits", "number of local inner SOR iterations", "PCSORSetIterations", jac->lits, &jac->lits, NULL));
799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupBegin("-pc_sor_symmetric", "SSOR, not SOR", "PCSORSetSymmetric", &flg));
809566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_SYMMETRIC_SWEEP));
819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_backward", "use backward sweep instead of forward", "PCSORSetSymmetric", &flg));
829566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_BACKWARD_SWEEP));
839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_forward", "use forward sweep", "PCSORSetSymmetric", &flg));
849566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_FORWARD_SWEEP));
859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_symmetric", "use SSOR separately on each processor", "PCSORSetSymmetric", &flg));
869566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_SYMMETRIC_SWEEP));
879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_backward", "use backward sweep locally", "PCSORSetSymmetric", &flg));
889566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_BACKWARD_SWEEP));
899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupEnd("-pc_sor_local_forward", "use forward sweep locally", "PCSORSetSymmetric", &flg));
909566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_FORWARD_SWEEP));
91d0609cedSBarry Smith PetscOptionsHeadEnd();
923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
934b9ad928SBarry Smith }
944b9ad928SBarry Smith
PCView_SOR(PC pc,PetscViewer viewer)9566976f2fSJacob Faibussowitsch static PetscErrorCode PCView_SOR(PC pc, PetscViewer viewer)
96d71ae5a4SJacob Faibussowitsch {
974b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
984b9ad928SBarry Smith MatSORType sym = jac->sym;
992fc52814SBarry Smith const char *sortype;
100*9f196a02SMartin Diehl PetscBool isascii;
1014b9ad928SBarry Smith
1024b9ad928SBarry Smith PetscFunctionBegin;
103*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
104*9f196a02SMartin Diehl if (isascii) {
1059566063dSJacob Faibussowitsch if (sym & SOR_ZERO_INITIAL_GUESS) PetscCall(PetscViewerASCIIPrintf(viewer, " zero initial guess\n"));
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";
11663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type = %s, iterations = %" PetscInt_FMT ", local iterations = %" PetscInt_FMT ", omega = %g\n", sortype, jac->its, jac->lits, (double)jac->omega));
1174b9ad928SBarry Smith }
1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith
PCSORSetSymmetric_SOR(PC pc,MatSORType flag)121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetSymmetric_SOR(PC pc, MatSORType flag)
122d71ae5a4SJacob Faibussowitsch {
123c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
1244b9ad928SBarry Smith
1254b9ad928SBarry Smith PetscFunctionBegin;
1264b9ad928SBarry Smith jac->sym = flag;
1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1284b9ad928SBarry Smith }
1294b9ad928SBarry Smith
PCSORSetOmega_SOR(PC pc,PetscReal omega)130d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetOmega_SOR(PC pc, PetscReal omega)
131d71ae5a4SJacob Faibussowitsch {
132c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
1334b9ad928SBarry Smith
1344b9ad928SBarry Smith PetscFunctionBegin;
1352472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range");
1364b9ad928SBarry Smith jac->omega = omega;
1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith
PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits)140d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetIterations_SOR(PC pc, PetscInt its, PetscInt lits)
141d71ae5a4SJacob Faibussowitsch {
142c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
1434b9ad928SBarry Smith
1444b9ad928SBarry Smith PetscFunctionBegin;
1454b9ad928SBarry Smith jac->its = its;
1464b9ad928SBarry Smith jac->lits = lits;
1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1484b9ad928SBarry Smith }
1494b9ad928SBarry Smith
PCSORGetSymmetric_SOR(PC pc,MatSORType * flag)150d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetSymmetric_SOR(PC pc, MatSORType *flag)
151d71ae5a4SJacob Faibussowitsch {
152c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
153c60c7ad4SBarry Smith
154c60c7ad4SBarry Smith PetscFunctionBegin;
155c60c7ad4SBarry Smith *flag = jac->sym;
1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
157c60c7ad4SBarry Smith }
158c60c7ad4SBarry Smith
PCSORGetOmega_SOR(PC pc,PetscReal * omega)159d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetOmega_SOR(PC pc, PetscReal *omega)
160d71ae5a4SJacob Faibussowitsch {
161c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
162c60c7ad4SBarry Smith
163c60c7ad4SBarry Smith PetscFunctionBegin;
164c60c7ad4SBarry Smith *omega = jac->omega;
1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
166c60c7ad4SBarry Smith }
167c60c7ad4SBarry Smith
PCSORGetIterations_SOR(PC pc,PetscInt * its,PetscInt * lits)168d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetIterations_SOR(PC pc, PetscInt *its, PetscInt *lits)
169d71ae5a4SJacob Faibussowitsch {
170c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data;
171c60c7ad4SBarry Smith
172c60c7ad4SBarry Smith PetscFunctionBegin;
173c60c7ad4SBarry Smith if (its) *its = jac->its;
174c60c7ad4SBarry Smith if (lits) *lits = jac->lits;
1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
176c60c7ad4SBarry Smith }
177c60c7ad4SBarry Smith
178c60c7ad4SBarry Smith /*@
1791ff2113eSBarry Smith PCSORGetSymmetric - Gets the form the SOR preconditioner is using; backward, or forward relaxation. The local variants perform SOR on
180c60c7ad4SBarry Smith each processor. By default forward relaxation is used.
181c60c7ad4SBarry Smith
182c3339decSBarry Smith Logically Collective
183c60c7ad4SBarry Smith
184c60c7ad4SBarry Smith Input Parameter:
185c60c7ad4SBarry Smith . pc - the preconditioner context
186c60c7ad4SBarry Smith
187c60c7ad4SBarry Smith Output Parameter:
188c60c7ad4SBarry Smith . flag - one of the following
189c60c7ad4SBarry Smith .vb
190c60c7ad4SBarry Smith SOR_FORWARD_SWEEP
191c60c7ad4SBarry Smith SOR_BACKWARD_SWEEP
192c60c7ad4SBarry Smith SOR_SYMMETRIC_SWEEP
193c60c7ad4SBarry Smith SOR_LOCAL_FORWARD_SWEEP
194c60c7ad4SBarry Smith SOR_LOCAL_BACKWARD_SWEEP
195c60c7ad4SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP
196c60c7ad4SBarry Smith .ve
197c60c7ad4SBarry Smith
198c60c7ad4SBarry Smith Options Database Keys:
199c60c7ad4SBarry Smith + -pc_sor_symmetric - Activates symmetric version
200c60c7ad4SBarry Smith . -pc_sor_backward - Activates backward version
201c60c7ad4SBarry Smith . -pc_sor_local_forward - Activates local forward version
202c60c7ad4SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version
203c60c7ad4SBarry Smith - -pc_sor_local_backward - Activates local backward version
204c60c7ad4SBarry Smith
205f1580f4eSBarry Smith Note:
206f1580f4eSBarry Smith To use the Eisenstat trick with SSOR, employ the `PCEISENSTAT` preconditioner,
207c60c7ad4SBarry Smith which can be chosen with the option
208c60c7ad4SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick
209c60c7ad4SBarry Smith
210c60c7ad4SBarry Smith Level: intermediate
211c60c7ad4SBarry Smith
212562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()`, `PCSORSetSymmetric()`
213c60c7ad4SBarry Smith @*/
PCSORGetSymmetric(PC pc,MatSORType * flag)214d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetSymmetric(PC pc, MatSORType *flag)
215d71ae5a4SJacob Faibussowitsch {
216c60c7ad4SBarry Smith PetscFunctionBegin;
217c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
218cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetSymmetric_C", (PC, MatSORType *), (pc, flag));
2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
220c60c7ad4SBarry Smith }
221c60c7ad4SBarry Smith
222c60c7ad4SBarry Smith /*@
223c60c7ad4SBarry Smith PCSORGetOmega - Gets the SOR relaxation coefficient, omega
224c60c7ad4SBarry Smith (where omega = 1.0 by default).
225c60c7ad4SBarry Smith
226c3339decSBarry Smith Logically Collective
227c60c7ad4SBarry Smith
228c60c7ad4SBarry Smith Input Parameter:
229c60c7ad4SBarry Smith . pc - the preconditioner context
230c60c7ad4SBarry Smith
231c60c7ad4SBarry Smith Output Parameter:
232c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2).
233c60c7ad4SBarry Smith
234c60c7ad4SBarry Smith Options Database Key:
235c60c7ad4SBarry Smith . -pc_sor_omega <omega> - Sets omega
236c60c7ad4SBarry Smith
237c60c7ad4SBarry Smith Level: intermediate
238c60c7ad4SBarry Smith
239562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `PCSORSetOmega()`
240c60c7ad4SBarry Smith @*/
PCSORGetOmega(PC pc,PetscReal * omega)241d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetOmega(PC pc, PetscReal *omega)
242d71ae5a4SJacob Faibussowitsch {
243c60c7ad4SBarry Smith PetscFunctionBegin;
244c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
245cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetOmega_C", (PC, PetscReal *), (pc, omega));
2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
247c60c7ad4SBarry Smith }
248c60c7ad4SBarry Smith
249c60c7ad4SBarry Smith /*@
250c60c7ad4SBarry Smith PCSORGetIterations - Gets the number of inner iterations to
251c60c7ad4SBarry Smith be used by the SOR preconditioner. The default is 1.
252c60c7ad4SBarry Smith
253c3339decSBarry Smith Logically Collective
254c60c7ad4SBarry Smith
255c60c7ad4SBarry Smith Input Parameter:
256c60c7ad4SBarry Smith . pc - the preconditioner context
257c60c7ad4SBarry Smith
258d8d19677SJose E. Roman Output Parameters:
259c60c7ad4SBarry Smith + lits - number of local iterations, smoothings over just variables on processor
260c60c7ad4SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations
261c60c7ad4SBarry Smith
262f1580f4eSBarry Smith Options Database Keys:
263c60c7ad4SBarry Smith + -pc_sor_its <its> - Sets number of iterations
264c60c7ad4SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations
265c60c7ad4SBarry Smith
266c60c7ad4SBarry Smith Level: intermediate
267c60c7ad4SBarry Smith
268f1580f4eSBarry Smith Note:
26995452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its
270c60c7ad4SBarry Smith
271562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()`, `PCSORSetIterations()`
272c60c7ad4SBarry Smith @*/
PCSORGetIterations(PC pc,PetscInt * its,PetscInt * lits)273d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetIterations(PC pc, PetscInt *its, PetscInt *lits)
274d71ae5a4SJacob Faibussowitsch {
275c60c7ad4SBarry Smith PetscFunctionBegin;
276c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
277cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetIterations_C", (PC, PetscInt *, PetscInt *), (pc, its, lits));
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
279c60c7ad4SBarry Smith }
280c60c7ad4SBarry Smith
2814b9ad928SBarry Smith /*@
2824b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR),
2834b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on
2844b9ad928SBarry Smith each processor. By default forward relaxation is used.
2854b9ad928SBarry Smith
286c3339decSBarry Smith Logically Collective
2874b9ad928SBarry Smith
2884b9ad928SBarry Smith Input Parameters:
2894b9ad928SBarry Smith + pc - the preconditioner context
2904b9ad928SBarry Smith - flag - one of the following
2914b9ad928SBarry Smith .vb
2924b9ad928SBarry Smith SOR_FORWARD_SWEEP
2934b9ad928SBarry Smith SOR_BACKWARD_SWEEP
2944b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP
2954b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP
2964b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP
2974b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP
2984b9ad928SBarry Smith .ve
2994b9ad928SBarry Smith
3004b9ad928SBarry Smith Options Database Keys:
3014b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version
3024b9ad928SBarry Smith . -pc_sor_backward - Activates backward version
3034b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version
3044b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version
3054b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version
3064b9ad928SBarry Smith
307f1580f4eSBarry Smith Note:
3084b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner,
3094b9ad928SBarry Smith which can be chosen with the option
3104b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick
3114b9ad928SBarry Smith
3124b9ad928SBarry Smith Level: intermediate
3134b9ad928SBarry Smith
314562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()`
3154b9ad928SBarry Smith @*/
PCSORSetSymmetric(PC pc,MatSORType flag)316d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetSymmetric(PC pc, MatSORType flag)
317d71ae5a4SJacob Faibussowitsch {
3184b9ad928SBarry Smith PetscFunctionBegin;
3190700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
320c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, flag, 2);
321cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetSymmetric_C", (PC, MatSORType), (pc, flag));
3223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3234b9ad928SBarry Smith }
3244b9ad928SBarry Smith
3254b9ad928SBarry Smith /*@
3264b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega
3274b9ad928SBarry Smith (where omega = 1.0 by default).
3284b9ad928SBarry Smith
329c3339decSBarry Smith Logically Collective
3304b9ad928SBarry Smith
3314b9ad928SBarry Smith Input Parameters:
3324b9ad928SBarry Smith + pc - the preconditioner context
3334b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2).
3344b9ad928SBarry Smith
3354b9ad928SBarry Smith Options Database Key:
3364b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega
3374b9ad928SBarry Smith
3384b9ad928SBarry Smith Level: intermediate
3394b9ad928SBarry Smith
340485168efSMatthew Knepley Note:
341f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODE`S option to `PETSC_FALSE` on the matrix.
342485168efSMatthew Knepley
343562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `MatSetOption()`
3444b9ad928SBarry Smith @*/
PCSORSetOmega(PC pc,PetscReal omega)345d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetOmega(PC pc, PetscReal omega)
346d71ae5a4SJacob Faibussowitsch {
3474b9ad928SBarry Smith PetscFunctionBegin;
348c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
349c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, omega, 2);
350cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetOmega_C", (PC, PetscReal), (pc, omega));
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3524b9ad928SBarry Smith }
3534b9ad928SBarry Smith
3544b9ad928SBarry Smith /*@
3554b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to
3564b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1.
3574b9ad928SBarry Smith
358c3339decSBarry Smith Logically Collective
3594b9ad928SBarry Smith
3604b9ad928SBarry Smith Input Parameters:
3614b9ad928SBarry Smith + pc - the preconditioner context
3624b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor
3634b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations
3644b9ad928SBarry Smith
365f1580f4eSBarry Smith Options Database Keys:
3664b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations
3674b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations
3684b9ad928SBarry Smith
3694b9ad928SBarry Smith Level: intermediate
3704b9ad928SBarry Smith
371f1580f4eSBarry Smith Note:
37295452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its
3734b9ad928SBarry Smith
374562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()`
3754b9ad928SBarry Smith @*/
PCSORSetIterations(PC pc,PetscInt its,PetscInt lits)376d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetIterations(PC pc, PetscInt its, PetscInt lits)
377d71ae5a4SJacob Faibussowitsch {
3784b9ad928SBarry Smith PetscFunctionBegin;
3790700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
380c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, its, 2);
381cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetIterations_C", (PC, PetscInt, PetscInt), (pc, its, lits));
3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3834b9ad928SBarry Smith }
3844b9ad928SBarry Smith
3854b9ad928SBarry Smith /*MC
3864b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning
3874b9ad928SBarry Smith
3884b9ad928SBarry Smith Options Database Keys:
3894b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version
3904b9ad928SBarry Smith . -pc_sor_backward - Activates backward version
391a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version
3924b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version
393a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version)
3944b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version
3954b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega
396422a814eSBarry Smith . -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal
397a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1)
398a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1)
3994b9ad928SBarry Smith
4004b9ad928SBarry Smith Level: beginner
4014b9ad928SBarry Smith
40295452b02SPatrick Sanan Notes:
403f1580f4eSBarry Smith Only implemented for the `MATAIJ` and `MATSEQBAIJ` matrix formats.
404f1580f4eSBarry Smith
4054b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block
4064b9ad928SBarry Smith Jacobi with SOR on each block.
4074b9ad928SBarry Smith
408f1580f4eSBarry Smith For `MATAIJ` matrices if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that
40948773899SPierre Jolivet zero will be used and hence the `KSPSolve()` will terminate with `KSP_DIVERGED_NANORINF`. If the option
410f1580f4eSBarry Smith `KSPSetErrorIfNotConverged()` or -ksp_error_if_not_converged the code will terminate as soon as it detects the
411422a814eSBarry Smith zero pivot.
412422a814eSBarry Smith
413f1580f4eSBarry Smith For `MATSEQBAIJ` matrices this implements point-block SOR, but the omega, its, lits options are not supported.
41437a17b4dSBarry Smith
415f1580f4eSBarry Smith For `MATSEQBAIJ` the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected
416422a814eSBarry Smith the computation is stopped with an error
417422a814eSBarry Smith
418f1580f4eSBarry Smith If used with `KSPRICHARDSON` and no monitors the convergence test is skipped to improve speed, thus it always iterates
419f1580f4eSBarry Smith the maximum number of iterations you've selected for `KSP`. It is usually used in this mode as a smoother for multigrid.
42067991ba4SBarry Smith
421f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODES` option to `PETSC_FALSE` on the matrix.
422485168efSMatthew Knepley
423562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`,
424db781477SPatrick Sanan `PCSORSetIterations()`, `PCSORSetSymmetric()`, `PCSORSetOmega()`, `PCEISENSTAT`, `MatSetOption()`
4254b9ad928SBarry Smith M*/
4264b9ad928SBarry Smith
PCCreate_SOR(PC pc)427d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc)
428d71ae5a4SJacob Faibussowitsch {
4294b9ad928SBarry Smith PC_SOR *jac;
4304b9ad928SBarry Smith
4314b9ad928SBarry Smith PetscFunctionBegin;
4324dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac));
4334b9ad928SBarry Smith
4344b9ad928SBarry Smith pc->ops->apply = PCApply_SOR;
4359d2471e0SBarry Smith pc->ops->applytranspose = PCApplyTranspose_SOR;
4364b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR;
4374b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR;
4380a545947SLisandro Dalcin pc->ops->setup = NULL;
4394b9ad928SBarry Smith pc->ops->view = PCView_SOR;
4404b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR;
4414b9ad928SBarry Smith pc->data = (void *)jac;
442d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP;
4434b9ad928SBarry Smith jac->omega = 1.0;
44496fc60bcSBarry Smith jac->fshift = 0.0;
4454b9ad928SBarry Smith jac->its = 1;
4464b9ad928SBarry Smith jac->lits = 1;
4474b9ad928SBarry Smith
4489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetSymmetric_C", PCSORSetSymmetric_SOR));
4499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetOmega_C", PCSORSetOmega_SOR));
4509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetIterations_C", PCSORSetIterations_SOR));
4519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetSymmetric_C", PCSORGetSymmetric_SOR));
4529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetOmega_C", PCSORGetOmega_SOR));
4539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetIterations_C", PCSORGetIterations_SOR));
4543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4554b9ad928SBarry Smith }
456