xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    Defines a  Eisenstat trick SSOR  preconditioner. This uses about
34b9ad928SBarry Smith  %50 of the usual amount of floating point ops used for SSOR + Krylov
44b9ad928SBarry Smith  method. But it requires actually solving the preconditioned problem
54b9ad928SBarry Smith  with both left and right preconditioning.
64b9ad928SBarry Smith */
7af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
84b9ad928SBarry Smith 
94b9ad928SBarry Smith typedef struct {
104b9ad928SBarry Smith   Mat       shell, A;
1178c391d7SBarry Smith   Vec       b[2], diag; /* temporary storage for true right hand side */
124b9ad928SBarry Smith   PetscReal omega;
13ace3abfcSBarry Smith   PetscBool usediag; /* indicates preconditioner should include diagonal scaling*/
144b9ad928SBarry Smith } PC_Eisenstat;
154b9ad928SBarry Smith 
16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMult_Eisenstat(Mat mat, Vec b, Vec x)
17d71ae5a4SJacob Faibussowitsch {
184b9ad928SBarry Smith   PC            pc;
194b9ad928SBarry Smith   PC_Eisenstat *eis;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &pc));
234b9ad928SBarry Smith   eis = (PC_Eisenstat *)pc->data;
249566063dSJacob Faibussowitsch   PetscCall(MatSOR(eis->A, b, eis->omega, SOR_EISENSTAT, 0.0, 1, 1, x));
2554e23b33SJose E. Roman   PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
274b9ad928SBarry Smith }
284b9ad928SBarry Smith 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCNorm_Eisenstat(Mat mat, NormType type, PetscReal *nrm)
30d71ae5a4SJacob Faibussowitsch {
31a9e2c90aSJose E. Roman   PC            pc;
32a9e2c90aSJose E. Roman   PC_Eisenstat *eis;
33a9e2c90aSJose E. Roman 
34a9e2c90aSJose E. Roman   PetscFunctionBegin;
35a9e2c90aSJose E. Roman   PetscCall(MatShellGetContext(mat, &pc));
36a9e2c90aSJose E. Roman   eis = (PC_Eisenstat *)pc->data;
37a9e2c90aSJose E. Roman   PetscCall(MatNorm(eis->A, type, nrm));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39a9e2c90aSJose E. Roman }
40a9e2c90aSJose E. Roman 
41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Eisenstat(PC pc, Vec x, Vec y)
42d71ae5a4SJacob Faibussowitsch {
434b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
44ace3abfcSBarry Smith   PetscBool     hasop;
454b9ad928SBarry Smith 
464b9ad928SBarry Smith   PetscFunctionBegin;
4789c6957cSBarry Smith   if (eis->usediag) {
489566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop));
4989c6957cSBarry Smith     if (hasop) {
509566063dSJacob Faibussowitsch       PetscCall(MatMultDiagonalBlock(pc->pmat, x, y));
5189c6957cSBarry Smith     } else {
529566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(y, x, eis->diag));
5389c6957cSBarry Smith     }
549566063dSJacob Faibussowitsch   } else PetscCall(VecCopy(x, y));
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
564b9ad928SBarry Smith }
574b9ad928SBarry Smith 
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Eisenstat(PC pc, Vec x, Vec y)
59d71ae5a4SJacob Faibussowitsch {
6052d511d9SJose E. Roman   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
6152d511d9SJose E. Roman   PetscBool     hasop, set, sym;
6252d511d9SJose E. Roman 
6352d511d9SJose E. Roman   PetscFunctionBegin;
6452d511d9SJose E. Roman   PetscCall(MatIsSymmetricKnown(eis->A, &set, &sym));
6552d511d9SJose E. Roman   PetscCheck(set && sym, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Can only apply transpose of Eisenstat if matrix is symmetric");
6652d511d9SJose E. Roman   if (eis->usediag) {
6752d511d9SJose E. Roman     PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop));
6852d511d9SJose E. Roman     if (hasop) {
6952d511d9SJose E. Roman       PetscCall(MatMultDiagonalBlock(pc->pmat, x, y));
7052d511d9SJose E. Roman     } else {
7152d511d9SJose E. Roman       PetscCall(VecPointwiseMult(y, x, eis->diag));
7252d511d9SJose E. Roman     }
7352d511d9SJose E. Roman   } else PetscCall(VecCopy(x, y));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7552d511d9SJose E. Roman }
7652d511d9SJose E. Roman 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x)
78d71ae5a4SJacob Faibussowitsch {
794b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
80ace3abfcSBarry Smith   PetscBool     nonzero;
814b9ad928SBarry Smith 
824b9ad928SBarry Smith   PetscFunctionBegin;
8378c391d7SBarry Smith   if (pc->presolvedone < 2) {
8408401ef6SPierre Jolivet     PetscCheck(pc->mat == pc->pmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot have different mat and pmat");
854b9ad928SBarry Smith     /* swap shell matrix and true matrix */
864b9ad928SBarry Smith     eis->A  = pc->mat;
874b9ad928SBarry Smith     pc->mat = eis->shell;
884b9ad928SBarry Smith   }
894b9ad928SBarry Smith 
904dfa11a4SJacob Faibussowitsch   if (!eis->b[pc->presolvedone - 1]) { PetscCall(VecDuplicate(b, &eis->b[pc->presolvedone - 1])); }
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   /* if nonzero initial guess, modify x */
939566063dSJacob Faibussowitsch   PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero));
944b9ad928SBarry Smith   if (nonzero) {
959566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, eis->b[pc->presolvedone - 1]));
969566063dSJacob Faibussowitsch     PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, SOR_APPLY_UPPER, 0.0, 1, 1, x));
9754e23b33SJose E. Roman     PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
984b9ad928SBarry Smith   }
994b9ad928SBarry Smith 
100121471adSBarry Smith   /* save true b, other option is to swap pointers */
1019566063dSJacob Faibussowitsch   PetscCall(VecCopy(b, eis->b[pc->presolvedone - 1]));
102121471adSBarry Smith 
1035c99c7daSBarry Smith   /* modify b by (L + D/omega)^{-1} */
1049566063dSJacob Faibussowitsch   PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP), 0.0, 1, 1, b));
10554e23b33SJose E. Roman   PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1074b9ad928SBarry Smith }
1084b9ad928SBarry Smith 
109d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPostSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x)
110d71ae5a4SJacob Faibussowitsch {
1114b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith   PetscFunctionBegin;
1144b9ad928SBarry Smith   /* get back true b */
1159566063dSJacob Faibussowitsch   PetscCall(VecCopy(eis->b[pc->presolvedone], b));
116121471adSBarry Smith 
117121471adSBarry Smith   /* modify x by (U + D/omega)^{-1} */
1189566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, eis->b[pc->presolvedone]));
1199566063dSJacob Faibussowitsch   PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone], eis->omega, (MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP), 0.0, 1, 1, x));
12054e23b33SJose E. Roman   PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
1212fa5cd67SKarl Rupp   if (!pc->presolvedone) pc->mat = eis->A;
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1234b9ad928SBarry Smith }
1244b9ad928SBarry Smith 
125d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Eisenstat(PC pc)
126d71ae5a4SJacob Faibussowitsch {
1274b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->b[0]));
1319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->b[1]));
1329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&eis->shell));
1339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->diag));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13569d2c0f9SBarry Smith }
13669d2c0f9SBarry Smith 
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Eisenstat(PC pc)
138d71ae5a4SJacob Faibussowitsch {
13969d2c0f9SBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PCReset_Eisenstat(pc));
1412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", NULL));
1422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", NULL));
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", NULL));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", NULL));
1452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1484b9ad928SBarry Smith }
1494b9ad928SBarry Smith 
150d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc, PetscOptionItems *PetscOptionsObject)
151d71ae5a4SJacob Faibussowitsch {
1524b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1538afaa268SBarry Smith   PetscBool     set, flg;
154d8e4f26eSJose E. Roman   PetscReal     omega;
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith   PetscFunctionBegin;
157d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Eisenstat SSOR options");
158d8e4f26eSJose E. Roman   PetscCall(PetscOptionsReal("-pc_eisenstat_omega", "Relaxation factor 0 < omega < 2", "PCEisenstatSetOmega", eis->omega, &omega, &flg));
159d8e4f26eSJose E. Roman   if (flg) PetscCall(PCEisenstatSetOmega(pc, omega));
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling", "Do not use standard diagonal scaling", "PCEisenstatSetNoDiagonalScaling", eis->usediag ? PETSC_FALSE : PETSC_TRUE, &flg, &set));
1611baa6e33SBarry Smith   if (set) PetscCall(PCEisenstatSetNoDiagonalScaling(pc, flg));
162d0609cedSBarry Smith   PetscOptionsHeadEnd();
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1644b9ad928SBarry Smith }
1654b9ad928SBarry Smith 
166d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Eisenstat(PC pc, PetscViewer viewer)
167d71ae5a4SJacob Faibussowitsch {
1684b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
169ace3abfcSBarry Smith   PetscBool     iascii;
1704b9ad928SBarry Smith 
1714b9ad928SBarry Smith   PetscFunctionBegin;
1729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
17332077d6dSBarry Smith   if (iascii) {
1749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  omega = %g\n", (double)eis->omega));
1754b9ad928SBarry Smith     if (eis->usediag) {
1769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Using diagonal scaling (default)\n"));
1774b9ad928SBarry Smith     } else {
1789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Not using diagonal scaling\n"));
1794b9ad928SBarry Smith     }
1804b9ad928SBarry Smith   }
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1824b9ad928SBarry Smith }
1834b9ad928SBarry Smith 
184d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Eisenstat(PC pc)
185d71ae5a4SJacob Faibussowitsch {
18613f74950SBarry Smith   PetscInt      M, N, m, n;
18752d511d9SJose E. Roman   PetscBool     set, sym;
1884b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   PetscFunctionBegin;
1914b9ad928SBarry Smith   if (!pc->setupcalled) {
1929566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pc->mat, &M, &N));
1939566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m, &n));
19452d511d9SJose E. Roman     PetscCall(MatIsSymmetricKnown(pc->mat, &set, &sym));
1959566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &eis->shell));
1969566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(eis->shell, m, n, M, N));
1979566063dSJacob Faibussowitsch     PetscCall(MatSetType(eis->shell, MATSHELL));
1989566063dSJacob Faibussowitsch     PetscCall(MatSetUp(eis->shell));
1999566063dSJacob Faibussowitsch     PetscCall(MatShellSetContext(eis->shell, pc));
2009566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT, (void (*)(void))PCMult_Eisenstat));
20152d511d9SJose E. Roman     if (set && sym) PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT_TRANSPOSE, (void (*)(void))PCMult_Eisenstat));
202a9e2c90aSJose E. Roman     PetscCall(MatShellSetOperation(eis->shell, MATOP_NORM, (void (*)(void))PCNorm_Eisenstat));
2034b9ad928SBarry Smith   }
2043ba16761SJacob Faibussowitsch   if (!eis->usediag) PetscFunctionReturn(PETSC_SUCCESS);
2054dfa11a4SJacob Faibussowitsch   if (!pc->setupcalled) { PetscCall(MatCreateVecs(pc->pmat, &eis->diag, NULL)); }
2069566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(pc->pmat, eis->diag));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2084b9ad928SBarry Smith }
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith /* --------------------------------------------------------------------*/
2114b9ad928SBarry Smith 
212d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc, PetscReal omega)
213d71ae5a4SJacob Faibussowitsch {
214c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith   PetscFunctionBegin;
2172472a847SBarry Smith   PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range");
2184b9ad928SBarry Smith   eis->omega = omega;
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2204b9ad928SBarry Smith }
2214b9ad928SBarry Smith 
222d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc, PetscBool flg)
223d71ae5a4SJacob Faibussowitsch {
224c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
2254b9ad928SBarry Smith 
2264b9ad928SBarry Smith   PetscFunctionBegin;
227c60c7ad4SBarry Smith   eis->usediag = flg;
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229c60c7ad4SBarry Smith }
230c60c7ad4SBarry Smith 
231d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc, PetscReal *omega)
232d71ae5a4SJacob Faibussowitsch {
233c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
234c60c7ad4SBarry Smith 
235c60c7ad4SBarry Smith   PetscFunctionBegin;
236c60c7ad4SBarry Smith   *omega = eis->omega;
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
238c60c7ad4SBarry Smith }
239c60c7ad4SBarry Smith 
240d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc, PetscBool *flg)
241d71ae5a4SJacob Faibussowitsch {
242c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
243c60c7ad4SBarry Smith 
244c60c7ad4SBarry Smith   PetscFunctionBegin;
245c60c7ad4SBarry Smith   *flg = eis->usediag;
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2474b9ad928SBarry Smith }
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith /*@
2504b9ad928SBarry Smith   PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
251f1580f4eSBarry Smith   to use with Eisenstat's trick (where omega = 1.0 by default)
2524b9ad928SBarry Smith 
253c3339decSBarry Smith   Logically Collective
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith   Input Parameters:
2564b9ad928SBarry Smith + pc    - the preconditioner context
2574b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2)
2584b9ad928SBarry Smith 
2594b9ad928SBarry Smith   Options Database Key:
2604b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega
2614b9ad928SBarry Smith 
2622fe279fdSBarry Smith   Level: intermediate
2632fe279fdSBarry Smith 
2644b9ad928SBarry Smith   Notes:
2654b9ad928SBarry Smith   The Eisenstat trick implementation of SSOR requires about 50% of the
2664b9ad928SBarry Smith   usual amount of floating point operations used for SSOR + Krylov method;
2674b9ad928SBarry Smith   however, the preconditioned problem must be solved with both left
2684b9ad928SBarry Smith   and right preconditioning.
2694b9ad928SBarry Smith 
270f1580f4eSBarry Smith   To use SSOR without the Eisenstat trick, employ the `PCSOR` preconditioner,
27137fdd005SBarry Smith   which can be chosen with the database options `-pc_type sor -pc_sor_symmetric`
2724b9ad928SBarry Smith 
273*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSORSetOmega()`, `PCEISENSTAT`
2744b9ad928SBarry Smith @*/
275d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatSetOmega(PC pc, PetscReal omega)
276d71ae5a4SJacob Faibussowitsch {
2774b9ad928SBarry Smith   PetscFunctionBegin;
2780700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
279c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, omega, 2);
280cac4c232SBarry Smith   PetscTryMethod(pc, "PCEisenstatSetOmega_C", (PC, PetscReal), (pc, omega));
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2824b9ad928SBarry Smith }
2834b9ad928SBarry Smith 
2844b9ad928SBarry Smith /*@
285f1580f4eSBarry Smith   PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner, `PCEISENSTAT`
2864b9ad928SBarry Smith   not to do additional diagonal preconditioning. For matrices with a constant
2874b9ad928SBarry Smith   along the diagonal, this may save a small amount of work.
2884b9ad928SBarry Smith 
289c3339decSBarry Smith   Logically Collective
2904b9ad928SBarry Smith 
291c60c7ad4SBarry Smith   Input Parameters:
292c60c7ad4SBarry Smith + pc  - the preconditioner context
293f1580f4eSBarry Smith - flg - `PETSC_TRUE` turns off diagonal scaling inside the algorithm
2944b9ad928SBarry Smith 
2954b9ad928SBarry Smith   Options Database Key:
296f1580f4eSBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith   Level: intermediate
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith   Note:
301f1580f4eSBarry Smith   If you use the `KSPSetDiagonalScaling()` or -ksp_diagonal_scale option then you will
3026aad120cSJose E. Roman   likely want to use this routine since it will save you some unneeded flops.
3034b9ad928SBarry Smith 
304*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCEisenstatSetOmega()`, `PCEISENSTAT`
3054b9ad928SBarry Smith @*/
306d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc, PetscBool flg)
307d71ae5a4SJacob Faibussowitsch {
3084b9ad928SBarry Smith   PetscFunctionBegin;
3090700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
310cac4c232SBarry Smith   PetscTryMethod(pc, "PCEisenstatSetNoDiagonalScaling_C", (PC, PetscBool), (pc, flg));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312c60c7ad4SBarry Smith }
313c60c7ad4SBarry Smith 
314c60c7ad4SBarry Smith /*@
315c60c7ad4SBarry Smith   PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega,
316c60c7ad4SBarry Smith   to use with Eisenstat's trick (where omega = 1.0 by default).
317c60c7ad4SBarry Smith 
318c3339decSBarry Smith   Logically Collective
319c60c7ad4SBarry Smith 
320c60c7ad4SBarry Smith   Input Parameter:
321c60c7ad4SBarry Smith . pc - the preconditioner context
322c60c7ad4SBarry Smith 
323c60c7ad4SBarry Smith   Output Parameter:
324c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2)
325c60c7ad4SBarry Smith 
326c60c7ad4SBarry Smith   Options Database Key:
327c60c7ad4SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega
328c60c7ad4SBarry Smith 
329c60c7ad4SBarry Smith   Notes:
330c60c7ad4SBarry Smith   The Eisenstat trick implementation of SSOR requires about 50% of the
331c60c7ad4SBarry Smith   usual amount of floating point operations used for SSOR + Krylov method;
332c60c7ad4SBarry Smith   however, the preconditioned problem must be solved with both left
333c60c7ad4SBarry Smith   and right preconditioning.
334c60c7ad4SBarry Smith 
335c60c7ad4SBarry Smith   To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
33637fdd005SBarry Smith   which can be chosen with the database options `-pc_type sor -pc_sor_symmetric`
337c60c7ad4SBarry Smith 
338c60c7ad4SBarry Smith   Level: intermediate
339c60c7ad4SBarry Smith 
340*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCEISENSTAT`, `PCSORGetOmega()`, `PCEisenstatSetOmega()`
341c60c7ad4SBarry Smith @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatGetOmega(PC pc, PetscReal *omega)
343d71ae5a4SJacob Faibussowitsch {
344c60c7ad4SBarry Smith   PetscFunctionBegin;
345c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
346cac4c232SBarry Smith   PetscUseMethod(pc, "PCEisenstatGetOmega_C", (PC, PetscReal *), (pc, omega));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348c60c7ad4SBarry Smith }
349c60c7ad4SBarry Smith 
350c60c7ad4SBarry Smith /*@
351163d334eSBarry Smith   PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner
352c60c7ad4SBarry Smith   not to do additional diagonal preconditioning. For matrices with a constant
353c60c7ad4SBarry Smith   along the diagonal, this may save a small amount of work.
354c60c7ad4SBarry Smith 
355c3339decSBarry Smith   Logically Collective
356c60c7ad4SBarry Smith 
357c60c7ad4SBarry Smith   Input Parameter:
358c60c7ad4SBarry Smith . pc - the preconditioner context
359c60c7ad4SBarry Smith 
360c60c7ad4SBarry Smith   Output Parameter:
361f1580f4eSBarry Smith . flg - `PETSC_TRUE` means there is no diagonal scaling applied
362c60c7ad4SBarry Smith 
363c60c7ad4SBarry Smith   Options Database Key:
364f1580f4eSBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
365c60c7ad4SBarry Smith 
366c60c7ad4SBarry Smith   Level: intermediate
367c60c7ad4SBarry Smith 
368c60c7ad4SBarry Smith   Note:
369f9ff08acSPierre Jolivet   If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will
3706aad120cSJose E. Roman   likely want to use this routine since it will save you some unneeded flops.
371c60c7ad4SBarry Smith 
372f1580f4eSBarry Smith .seealso: , `PCEISENSTAT`, `PCEisenstatGetOmega()`
373c60c7ad4SBarry Smith @*/
374d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc, PetscBool *flg)
375d71ae5a4SJacob Faibussowitsch {
376c60c7ad4SBarry Smith   PetscFunctionBegin;
377c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
378cac4c232SBarry Smith   PetscUseMethod(pc, "PCEisenstatGetNoDiagonalScaling_C", (PC, PetscBool *), (pc, flg));
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3804b9ad928SBarry Smith }
3814b9ad928SBarry Smith 
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool *change)
383d71ae5a4SJacob Faibussowitsch {
3848066bbecSBarry Smith   PetscFunctionBegin;
3858066bbecSBarry Smith   *change = PETSC_TRUE;
3863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3878066bbecSBarry Smith }
3888066bbecSBarry Smith 
3894b9ad928SBarry Smith /*MC
3904b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
3914b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Options Database Keys:
3944b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
395f1580f4eSBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Level: beginner
3984b9ad928SBarry Smith 
39995452b02SPatrick Sanan    Notes:
400f1580f4eSBarry Smith    Only implemented for the `MATAIJ` matrix format.
4014b9ad928SBarry Smith 
402f1580f4eSBarry Smith    Not a true parallel SOR, in parallel this implementation corresponds to block Jacobi with SOR on each block.
403f1580f4eSBarry Smith 
404f1580f4eSBarry Smith    Developer Note:
405f1580f4eSBarry Smith    Since this algorithm runs the Krylov method on a transformed linear system the implementation provides `PCPreSolve()` and `PCPostSolve()`
406f1580f4eSBarry Smith    routines that `KSP` uses to set up the transformed linear system.
407f1580f4eSBarry Smith 
408*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCEisenstatGetOmega()`,
409db781477SPatrick Sanan           `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR`
4104b9ad928SBarry Smith M*/
4114b9ad928SBarry Smith 
412d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc)
413d71ae5a4SJacob Faibussowitsch {
4144b9ad928SBarry Smith   PC_Eisenstat *eis;
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith   PetscFunctionBegin;
4174dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&eis));
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
42052d511d9SJose E. Roman   pc->ops->applytranspose  = PCApplyTranspose_Eisenstat;
421dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
422dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
4230a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
4244b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
4254b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
42669d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Eisenstat;
4274b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
4284b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
4294b9ad928SBarry Smith 
4303ec1f749SStefano Zampini   pc->data     = eis;
4314b9ad928SBarry Smith   eis->omega   = 1.0;
4320a545947SLisandro Dalcin   eis->b[0]    = NULL;
4330a545947SLisandro Dalcin   eis->b[1]    = NULL;
4340a545947SLisandro Dalcin   eis->diag    = NULL;
4354b9ad928SBarry Smith   eis->usediag = PETSC_TRUE;
4364b9ad928SBarry Smith 
4379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", PCEisenstatSetOmega_Eisenstat));
4389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", PCEisenstatSetNoDiagonalScaling_Eisenstat));
4399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", PCEisenstatGetOmega_Eisenstat));
4409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", PCEisenstatGetNoDiagonalScaling_Eisenstat));
4419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Eisenstat));
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4434b9ad928SBarry Smith }
444