xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Defines a  Eisenstat trick SSOR  preconditioner. This uses about
44b9ad928SBarry Smith  %50 of the usual amount of floating point ops used for SSOR + Krylov
54b9ad928SBarry Smith  method. But it requires actually solving the preconditioned problem
64b9ad928SBarry Smith  with both left and right preconditioning.
74b9ad928SBarry Smith */
8af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
94b9ad928SBarry Smith 
104b9ad928SBarry Smith typedef struct {
114b9ad928SBarry Smith   Mat       shell, A;
1278c391d7SBarry Smith   Vec       b[2], diag; /* temporary storage for true right hand side */
134b9ad928SBarry Smith   PetscReal omega;
14ace3abfcSBarry Smith   PetscBool usediag; /* indicates preconditioner should include diagonal scaling*/
154b9ad928SBarry Smith } PC_Eisenstat;
164b9ad928SBarry Smith 
179371c9d4SSatish Balay static PetscErrorCode PCMult_Eisenstat(Mat mat, Vec b, Vec x) {
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));
264b9ad928SBarry Smith   PetscFunctionReturn(0);
274b9ad928SBarry Smith }
284b9ad928SBarry Smith 
29a9e2c90aSJose E. Roman static PetscErrorCode PCNorm_Eisenstat(Mat mat, NormType type, PetscReal *nrm) {
30a9e2c90aSJose E. Roman   PC            pc;
31a9e2c90aSJose E. Roman   PC_Eisenstat *eis;
32a9e2c90aSJose E. Roman 
33a9e2c90aSJose E. Roman   PetscFunctionBegin;
34a9e2c90aSJose E. Roman   PetscCall(MatShellGetContext(mat, &pc));
35a9e2c90aSJose E. Roman   eis = (PC_Eisenstat *)pc->data;
36a9e2c90aSJose E. Roman   PetscCall(MatNorm(eis->A, type, nrm));
37a9e2c90aSJose E. Roman   PetscFunctionReturn(0);
38a9e2c90aSJose E. Roman }
39a9e2c90aSJose E. Roman 
409371c9d4SSatish Balay static PetscErrorCode PCApply_Eisenstat(PC pc, Vec x, Vec y) {
414b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
42ace3abfcSBarry Smith   PetscBool     hasop;
434b9ad928SBarry Smith 
444b9ad928SBarry Smith   PetscFunctionBegin;
4589c6957cSBarry Smith   if (eis->usediag) {
469566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop));
4789c6957cSBarry Smith     if (hasop) {
489566063dSJacob Faibussowitsch       PetscCall(MatMultDiagonalBlock(pc->pmat, x, y));
4989c6957cSBarry Smith     } else {
509566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(y, x, eis->diag));
5189c6957cSBarry Smith     }
529566063dSJacob Faibussowitsch   } else PetscCall(VecCopy(x, y));
534b9ad928SBarry Smith   PetscFunctionReturn(0);
544b9ad928SBarry Smith }
554b9ad928SBarry Smith 
5652d511d9SJose E. Roman static PetscErrorCode PCApplyTranspose_Eisenstat(PC pc, Vec x, Vec y) {
5752d511d9SJose E. Roman   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
5852d511d9SJose E. Roman   PetscBool     hasop, set, sym;
5952d511d9SJose E. Roman 
6052d511d9SJose E. Roman   PetscFunctionBegin;
6152d511d9SJose E. Roman   PetscCall(MatIsSymmetricKnown(eis->A, &set, &sym));
6252d511d9SJose E. Roman   PetscCheck(set && sym, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Can only apply transpose of Eisenstat if matrix is symmetric");
6352d511d9SJose E. Roman   if (eis->usediag) {
6452d511d9SJose E. Roman     PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop));
6552d511d9SJose E. Roman     if (hasop) {
6652d511d9SJose E. Roman       PetscCall(MatMultDiagonalBlock(pc->pmat, x, y));
6752d511d9SJose E. Roman     } else {
6852d511d9SJose E. Roman       PetscCall(VecPointwiseMult(y, x, eis->diag));
6952d511d9SJose E. Roman     }
7052d511d9SJose E. Roman   } else PetscCall(VecCopy(x, y));
7152d511d9SJose E. Roman   PetscFunctionReturn(0);
7252d511d9SJose E. Roman }
7352d511d9SJose E. Roman 
749371c9d4SSatish Balay static PetscErrorCode PCPreSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) {
754b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
76ace3abfcSBarry Smith   PetscBool     nonzero;
774b9ad928SBarry Smith 
784b9ad928SBarry Smith   PetscFunctionBegin;
7978c391d7SBarry Smith   if (pc->presolvedone < 2) {
8008401ef6SPierre Jolivet     PetscCheck(pc->mat == pc->pmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot have different mat and pmat");
814b9ad928SBarry Smith     /* swap shell matrix and true matrix */
824b9ad928SBarry Smith     eis->A  = pc->mat;
834b9ad928SBarry Smith     pc->mat = eis->shell;
844b9ad928SBarry Smith   }
854b9ad928SBarry Smith 
8678c391d7SBarry Smith   if (!eis->b[pc->presolvedone - 1]) {
879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(b, &eis->b[pc->presolvedone - 1]));
889566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)eis->b[pc->presolvedone - 1]));
8978c391d7SBarry Smith   }
904b9ad928SBarry Smith 
914b9ad928SBarry Smith   /* if nonzero initial guess, modify x */
929566063dSJacob Faibussowitsch   PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero));
934b9ad928SBarry Smith   if (nonzero) {
949566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, eis->b[pc->presolvedone - 1]));
959566063dSJacob Faibussowitsch     PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, SOR_APPLY_UPPER, 0.0, 1, 1, x));
9654e23b33SJose E. Roman     PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
974b9ad928SBarry Smith   }
984b9ad928SBarry Smith 
99121471adSBarry Smith   /* save true b, other option is to swap pointers */
1009566063dSJacob Faibussowitsch   PetscCall(VecCopy(b, eis->b[pc->presolvedone - 1]));
101121471adSBarry Smith 
1025c99c7daSBarry Smith   /* modify b by (L + D/omega)^{-1} */
1039566063dSJacob 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));
10454e23b33SJose E. Roman   PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
1054b9ad928SBarry Smith   PetscFunctionReturn(0);
1064b9ad928SBarry Smith }
1074b9ad928SBarry Smith 
1089371c9d4SSatish Balay static PetscErrorCode PCPostSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) {
1094b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith   PetscFunctionBegin;
1124b9ad928SBarry Smith   /* get back true b */
1139566063dSJacob Faibussowitsch   PetscCall(VecCopy(eis->b[pc->presolvedone], b));
114121471adSBarry Smith 
115121471adSBarry Smith   /* modify x by (U + D/omega)^{-1} */
1169566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, eis->b[pc->presolvedone]));
1179566063dSJacob 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));
11854e23b33SJose E. Roman   PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason));
1192fa5cd67SKarl Rupp   if (!pc->presolvedone) pc->mat = eis->A;
1204b9ad928SBarry Smith   PetscFunctionReturn(0);
1214b9ad928SBarry Smith }
1224b9ad928SBarry Smith 
1239371c9d4SSatish Balay static PetscErrorCode PCReset_Eisenstat(PC pc) {
1244b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith   PetscFunctionBegin;
1279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->b[0]));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->b[1]));
1299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&eis->shell));
1309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&eis->diag));
13169d2c0f9SBarry Smith   PetscFunctionReturn(0);
13269d2c0f9SBarry Smith }
13369d2c0f9SBarry Smith 
1349371c9d4SSatish Balay static PetscErrorCode PCDestroy_Eisenstat(PC pc) {
13569d2c0f9SBarry Smith   PetscFunctionBegin;
1369566063dSJacob Faibussowitsch   PetscCall(PCReset_Eisenstat(pc));
1372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", NULL));
1382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", NULL));
1392e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", NULL));
1402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", NULL));
1412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1434b9ad928SBarry Smith   PetscFunctionReturn(0);
1444b9ad928SBarry Smith }
1454b9ad928SBarry Smith 
1469371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc, PetscOptionItems *PetscOptionsObject) {
1474b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1488afaa268SBarry Smith   PetscBool     set, flg;
1494b9ad928SBarry Smith 
1504b9ad928SBarry Smith   PetscFunctionBegin;
151d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Eisenstat SSOR options");
1529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_eisenstat_omega", "Relaxation factor 0 < omega < 2", "PCEisenstatSetOmega", eis->omega, &eis->omega, NULL));
1539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling", "Do not use standard diagonal scaling", "PCEisenstatSetNoDiagonalScaling", eis->usediag ? PETSC_FALSE : PETSC_TRUE, &flg, &set));
1541baa6e33SBarry Smith   if (set) PetscCall(PCEisenstatSetNoDiagonalScaling(pc, flg));
155d0609cedSBarry Smith   PetscOptionsHeadEnd();
1564b9ad928SBarry Smith   PetscFunctionReturn(0);
1574b9ad928SBarry Smith }
1584b9ad928SBarry Smith 
1599371c9d4SSatish Balay static PetscErrorCode PCView_Eisenstat(PC pc, PetscViewer viewer) {
1604b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
161ace3abfcSBarry Smith   PetscBool     iascii;
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
16532077d6dSBarry Smith   if (iascii) {
1669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  omega = %g\n", (double)eis->omega));
1674b9ad928SBarry Smith     if (eis->usediag) {
1689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Using diagonal scaling (default)\n"));
1694b9ad928SBarry Smith     } else {
1709566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Not using diagonal scaling\n"));
1714b9ad928SBarry Smith     }
1724b9ad928SBarry Smith   }
1734b9ad928SBarry Smith   PetscFunctionReturn(0);
1744b9ad928SBarry Smith }
1754b9ad928SBarry Smith 
1769371c9d4SSatish Balay static PetscErrorCode PCSetUp_Eisenstat(PC pc) {
17713f74950SBarry Smith   PetscInt      M, N, m, n;
17852d511d9SJose E. Roman   PetscBool     set, sym;
1794b9ad928SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
1804b9ad928SBarry Smith 
1814b9ad928SBarry Smith   PetscFunctionBegin;
1824b9ad928SBarry Smith   if (!pc->setupcalled) {
1839566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pc->mat, &M, &N));
1849566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(pc->mat, &m, &n));
18552d511d9SJose E. Roman     PetscCall(MatIsSymmetricKnown(pc->mat, &set, &sym));
1869566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &eis->shell));
1879566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(eis->shell, m, n, M, N));
1889566063dSJacob Faibussowitsch     PetscCall(MatSetType(eis->shell, MATSHELL));
1899566063dSJacob Faibussowitsch     PetscCall(MatSetUp(eis->shell));
1909566063dSJacob Faibussowitsch     PetscCall(MatShellSetContext(eis->shell, pc));
1919566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)eis->shell));
1929566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT, (void (*)(void))PCMult_Eisenstat));
19352d511d9SJose E. Roman     if (set && sym) PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT_TRANSPOSE, (void (*)(void))PCMult_Eisenstat));
194a9e2c90aSJose E. Roman     PetscCall(MatShellSetOperation(eis->shell, MATOP_NORM, (void (*)(void))PCNorm_Eisenstat));
1954b9ad928SBarry Smith   }
1964b9ad928SBarry Smith   if (!eis->usediag) PetscFunctionReturn(0);
1974b9ad928SBarry Smith   if (!pc->setupcalled) {
1989566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &eis->diag, NULL));
1999566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)eis->diag));
2004b9ad928SBarry Smith   }
2019566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(pc->pmat, eis->diag));
2024b9ad928SBarry Smith   PetscFunctionReturn(0);
2034b9ad928SBarry Smith }
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith /* --------------------------------------------------------------------*/
2064b9ad928SBarry Smith 
2079371c9d4SSatish Balay static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc, PetscReal omega) {
208c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith   PetscFunctionBegin;
2112472a847SBarry Smith   PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range");
2124b9ad928SBarry Smith   eis->omega = omega;
2134b9ad928SBarry Smith   PetscFunctionReturn(0);
2144b9ad928SBarry Smith }
2154b9ad928SBarry Smith 
2169371c9d4SSatish Balay static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc, PetscBool flg) {
217c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
2184b9ad928SBarry Smith 
2194b9ad928SBarry Smith   PetscFunctionBegin;
220c60c7ad4SBarry Smith   eis->usediag = flg;
221c60c7ad4SBarry Smith   PetscFunctionReturn(0);
222c60c7ad4SBarry Smith }
223c60c7ad4SBarry Smith 
2249371c9d4SSatish Balay static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc, PetscReal *omega) {
225c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
226c60c7ad4SBarry Smith 
227c60c7ad4SBarry Smith   PetscFunctionBegin;
228c60c7ad4SBarry Smith   *omega = eis->omega;
229c60c7ad4SBarry Smith   PetscFunctionReturn(0);
230c60c7ad4SBarry Smith }
231c60c7ad4SBarry Smith 
2329371c9d4SSatish Balay static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc, PetscBool *flg) {
233c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
234c60c7ad4SBarry Smith 
235c60c7ad4SBarry Smith   PetscFunctionBegin;
236c60c7ad4SBarry Smith   *flg = eis->usediag;
2374b9ad928SBarry Smith   PetscFunctionReturn(0);
2384b9ad928SBarry Smith }
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith /*@
2414b9ad928SBarry Smith    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
242*f1580f4eSBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default)
2434b9ad928SBarry Smith 
244*f1580f4eSBarry Smith    Logically Collective on pc
2454b9ad928SBarry Smith 
2464b9ad928SBarry Smith    Input Parameters:
2474b9ad928SBarry Smith +  pc - the preconditioner context
2484b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2)
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith    Options Database Key:
2514b9ad928SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
2524b9ad928SBarry Smith 
2534b9ad928SBarry Smith    Notes:
2544b9ad928SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
2554b9ad928SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
2564b9ad928SBarry Smith    however, the preconditioned problem must be solved with both left
2574b9ad928SBarry Smith    and right preconditioning.
2584b9ad928SBarry Smith 
259*f1580f4eSBarry Smith    To use SSOR without the Eisenstat trick, employ the `PCSOR` preconditioner,
2604b9ad928SBarry Smith    which can be chosen with the database options
2614b9ad928SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
2624b9ad928SBarry Smith 
2634b9ad928SBarry Smith    Level: intermediate
2644b9ad928SBarry Smith 
265*f1580f4eSBarry Smith .seealso: `PCSORSetOmega()`, `PCEISENSTAT`
2664b9ad928SBarry Smith @*/
2679371c9d4SSatish Balay PetscErrorCode PCEisenstatSetOmega(PC pc, PetscReal omega) {
2684b9ad928SBarry Smith   PetscFunctionBegin;
2690700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
270c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc, omega, 2);
271cac4c232SBarry Smith   PetscTryMethod(pc, "PCEisenstatSetOmega_C", (PC, PetscReal), (pc, omega));
2724b9ad928SBarry Smith   PetscFunctionReturn(0);
2734b9ad928SBarry Smith }
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith /*@
276*f1580f4eSBarry Smith    PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner, `PCEISENSTAT`
2774b9ad928SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
2784b9ad928SBarry Smith    along the diagonal, this may save a small amount of work.
2794b9ad928SBarry Smith 
280*f1580f4eSBarry Smith    Logically Collective on pc
2814b9ad928SBarry Smith 
282c60c7ad4SBarry Smith    Input Parameters:
283c60c7ad4SBarry Smith +  pc - the preconditioner context
284*f1580f4eSBarry Smith -  flg - `PETSC_TRUE` turns off diagonal scaling inside the algorithm
2854b9ad928SBarry Smith 
2864b9ad928SBarry Smith    Options Database Key:
287*f1580f4eSBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
2884b9ad928SBarry Smith 
2894b9ad928SBarry Smith    Level: intermediate
2904b9ad928SBarry Smith 
2914b9ad928SBarry Smith    Note:
292*f1580f4eSBarry Smith      If you use the `KSPSetDiagonalScaling()` or -ksp_diagonal_scale option then you will
2936aad120cSJose E. Roman    likely want to use this routine since it will save you some unneeded flops.
2944b9ad928SBarry Smith 
295*f1580f4eSBarry Smith .seealso: `PCEisenstatSetOmega()`, `PCEISENSTAT`
2964b9ad928SBarry Smith @*/
2979371c9d4SSatish Balay PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc, PetscBool flg) {
2984b9ad928SBarry Smith   PetscFunctionBegin;
2990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
300cac4c232SBarry Smith   PetscTryMethod(pc, "PCEisenstatSetNoDiagonalScaling_C", (PC, PetscBool), (pc, flg));
301c60c7ad4SBarry Smith   PetscFunctionReturn(0);
302c60c7ad4SBarry Smith }
303c60c7ad4SBarry Smith 
304c60c7ad4SBarry Smith /*@
305c60c7ad4SBarry Smith    PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega,
306c60c7ad4SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
307c60c7ad4SBarry Smith 
308*f1580f4eSBarry Smith    Logically Collective on pc
309c60c7ad4SBarry Smith 
310c60c7ad4SBarry Smith    Input Parameter:
311c60c7ad4SBarry Smith .  pc - the preconditioner context
312c60c7ad4SBarry Smith 
313c60c7ad4SBarry Smith    Output Parameter:
314c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2)
315c60c7ad4SBarry Smith 
316c60c7ad4SBarry Smith    Options Database Key:
317c60c7ad4SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
318c60c7ad4SBarry Smith 
319c60c7ad4SBarry Smith    Notes:
320c60c7ad4SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
321c60c7ad4SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
322c60c7ad4SBarry Smith    however, the preconditioned problem must be solved with both left
323c60c7ad4SBarry Smith    and right preconditioning.
324c60c7ad4SBarry Smith 
325c60c7ad4SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
326c60c7ad4SBarry Smith    which can be chosen with the database options
327c60c7ad4SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
328c60c7ad4SBarry Smith 
329c60c7ad4SBarry Smith    Level: intermediate
330c60c7ad4SBarry Smith 
331*f1580f4eSBarry Smith .seealso: `PCEISENSTAT`, `PCSORGetOmega()`, `PCEisenstatSetOmega()`
332c60c7ad4SBarry Smith @*/
3339371c9d4SSatish Balay PetscErrorCode PCEisenstatGetOmega(PC pc, PetscReal *omega) {
334c60c7ad4SBarry Smith   PetscFunctionBegin;
335c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
336cac4c232SBarry Smith   PetscUseMethod(pc, "PCEisenstatGetOmega_C", (PC, PetscReal *), (pc, omega));
337c60c7ad4SBarry Smith   PetscFunctionReturn(0);
338c60c7ad4SBarry Smith }
339c60c7ad4SBarry Smith 
340c60c7ad4SBarry Smith /*@
341163d334eSBarry Smith    PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner
342c60c7ad4SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
343c60c7ad4SBarry Smith    along the diagonal, this may save a small amount of work.
344c60c7ad4SBarry Smith 
345*f1580f4eSBarry Smith    Logically Collective on pc
346c60c7ad4SBarry Smith 
347c60c7ad4SBarry Smith    Input Parameter:
348c60c7ad4SBarry Smith .  pc - the preconditioner context
349c60c7ad4SBarry Smith 
350c60c7ad4SBarry Smith    Output Parameter:
351*f1580f4eSBarry Smith .  flg - `PETSC_TRUE` means there is no diagonal scaling applied
352c60c7ad4SBarry Smith 
353c60c7ad4SBarry Smith    Options Database Key:
354*f1580f4eSBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
355c60c7ad4SBarry Smith 
356c60c7ad4SBarry Smith    Level: intermediate
357c60c7ad4SBarry Smith 
358c60c7ad4SBarry Smith    Note:
359f9ff08acSPierre Jolivet      If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will
3606aad120cSJose E. Roman    likely want to use this routine since it will save you some unneeded flops.
361c60c7ad4SBarry Smith 
362*f1580f4eSBarry Smith .seealso: , `PCEISENSTAT`, `PCEisenstatGetOmega()`
363c60c7ad4SBarry Smith @*/
3649371c9d4SSatish Balay PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc, PetscBool *flg) {
365c60c7ad4SBarry Smith   PetscFunctionBegin;
366c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
367cac4c232SBarry Smith   PetscUseMethod(pc, "PCEisenstatGetNoDiagonalScaling_C", (PC, PetscBool *), (pc, flg));
3684b9ad928SBarry Smith   PetscFunctionReturn(0);
3694b9ad928SBarry Smith }
3704b9ad928SBarry Smith 
3719371c9d4SSatish Balay static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool *change) {
3728066bbecSBarry Smith   PetscFunctionBegin;
3738066bbecSBarry Smith   *change = PETSC_TRUE;
3748066bbecSBarry Smith   PetscFunctionReturn(0);
3758066bbecSBarry Smith }
3768066bbecSBarry Smith 
3774b9ad928SBarry Smith /*MC
3784b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
3794b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith    Options Database Keys:
3824b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
383*f1580f4eSBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()`
3844b9ad928SBarry Smith 
3854b9ad928SBarry Smith    Level: beginner
3864b9ad928SBarry Smith 
38795452b02SPatrick Sanan    Notes:
388*f1580f4eSBarry Smith    Only implemented for the `MATAIJ` matrix format.
3894b9ad928SBarry Smith 
390*f1580f4eSBarry Smith    Not a true parallel SOR, in parallel this implementation corresponds to block Jacobi with SOR on each block.
391*f1580f4eSBarry Smith 
392*f1580f4eSBarry Smith    Developer Note:
393*f1580f4eSBarry Smith    Since this algorithm runs the Krylov method on a transformed linear system the implementation provides `PCPreSolve()` and `PCPostSolve()`
394*f1580f4eSBarry Smith    routines that `KSP` uses to set up the transformed linear system.
395*f1580f4eSBarry Smith 
396*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCEisenstatGetOmega()`,
397db781477SPatrick Sanan           `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR`
3984b9ad928SBarry Smith M*/
3994b9ad928SBarry Smith 
4009371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) {
4014b9ad928SBarry Smith   PC_Eisenstat *eis;
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith   PetscFunctionBegin;
4049566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &eis));
4054b9ad928SBarry Smith 
4064b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
40752d511d9SJose E. Roman   pc->ops->applytranspose  = PCApplyTranspose_Eisenstat;
408dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
409dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
4100a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
4114b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
4124b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
41369d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Eisenstat;
4144b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
4154b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
4164b9ad928SBarry Smith 
4173ec1f749SStefano Zampini   pc->data     = eis;
4184b9ad928SBarry Smith   eis->omega   = 1.0;
4190a545947SLisandro Dalcin   eis->b[0]    = NULL;
4200a545947SLisandro Dalcin   eis->b[1]    = NULL;
4210a545947SLisandro Dalcin   eis->diag    = NULL;
4224b9ad928SBarry Smith   eis->usediag = PETSC_TRUE;
4234b9ad928SBarry Smith 
4249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", PCEisenstatSetOmega_Eisenstat));
4259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", PCEisenstatSetNoDiagonalScaling_Eisenstat));
4269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", PCEisenstatGetOmega_Eisenstat));
4279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", PCEisenstatGetNoDiagonalScaling_Eisenstat));
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Eisenstat));
4294b9ad928SBarry Smith   PetscFunctionReturn(0);
4304b9ad928SBarry Smith }
431