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 86*4dfa11a4SJacob Faibussowitsch if (!eis->b[pc->presolvedone - 1]) { PetscCall(VecDuplicate(b, &eis->b[pc->presolvedone - 1])); } 874b9ad928SBarry Smith 884b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 899566063dSJacob Faibussowitsch PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero)); 904b9ad928SBarry Smith if (nonzero) { 919566063dSJacob Faibussowitsch PetscCall(VecCopy(x, eis->b[pc->presolvedone - 1])); 929566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, SOR_APPLY_UPPER, 0.0, 1, 1, x)); 9354e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 944b9ad928SBarry Smith } 954b9ad928SBarry Smith 96121471adSBarry Smith /* save true b, other option is to swap pointers */ 979566063dSJacob Faibussowitsch PetscCall(VecCopy(b, eis->b[pc->presolvedone - 1])); 98121471adSBarry Smith 995c99c7daSBarry Smith /* modify b by (L + D/omega)^{-1} */ 1009566063dSJacob 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)); 10154e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 1024b9ad928SBarry Smith PetscFunctionReturn(0); 1034b9ad928SBarry Smith } 1044b9ad928SBarry Smith 1059371c9d4SSatish Balay static PetscErrorCode PCPostSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) { 1064b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith PetscFunctionBegin; 1094b9ad928SBarry Smith /* get back true b */ 1109566063dSJacob Faibussowitsch PetscCall(VecCopy(eis->b[pc->presolvedone], b)); 111121471adSBarry Smith 112121471adSBarry Smith /* modify x by (U + D/omega)^{-1} */ 1139566063dSJacob Faibussowitsch PetscCall(VecCopy(x, eis->b[pc->presolvedone])); 1149566063dSJacob 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)); 11554e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 1162fa5cd67SKarl Rupp if (!pc->presolvedone) pc->mat = eis->A; 1174b9ad928SBarry Smith PetscFunctionReturn(0); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 1209371c9d4SSatish Balay static PetscErrorCode PCReset_Eisenstat(PC pc) { 1214b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith PetscFunctionBegin; 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[0])); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[1])); 1269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&eis->shell)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->diag)); 12869d2c0f9SBarry Smith PetscFunctionReturn(0); 12969d2c0f9SBarry Smith } 13069d2c0f9SBarry Smith 1319371c9d4SSatish Balay static PetscErrorCode PCDestroy_Eisenstat(PC pc) { 13269d2c0f9SBarry Smith PetscFunctionBegin; 1339566063dSJacob Faibussowitsch PetscCall(PCReset_Eisenstat(pc)); 1342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", NULL)); 1352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", NULL)); 1362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", NULL)); 1372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", NULL)); 1382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1404b9ad928SBarry Smith PetscFunctionReturn(0); 1414b9ad928SBarry Smith } 1424b9ad928SBarry Smith 1439371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc, PetscOptionItems *PetscOptionsObject) { 1444b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1458afaa268SBarry Smith PetscBool set, flg; 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith PetscFunctionBegin; 148d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Eisenstat SSOR options"); 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_eisenstat_omega", "Relaxation factor 0 < omega < 2", "PCEisenstatSetOmega", eis->omega, &eis->omega, NULL)); 1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling", "Do not use standard diagonal scaling", "PCEisenstatSetNoDiagonalScaling", eis->usediag ? PETSC_FALSE : PETSC_TRUE, &flg, &set)); 1511baa6e33SBarry Smith if (set) PetscCall(PCEisenstatSetNoDiagonalScaling(pc, flg)); 152d0609cedSBarry Smith PetscOptionsHeadEnd(); 1534b9ad928SBarry Smith PetscFunctionReturn(0); 1544b9ad928SBarry Smith } 1554b9ad928SBarry Smith 1569371c9d4SSatish Balay static PetscErrorCode PCView_Eisenstat(PC pc, PetscViewer viewer) { 1574b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 158ace3abfcSBarry Smith PetscBool iascii; 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith PetscFunctionBegin; 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 16232077d6dSBarry Smith if (iascii) { 1639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " omega = %g\n", (double)eis->omega)); 1644b9ad928SBarry Smith if (eis->usediag) { 1659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using diagonal scaling (default)\n")); 1664b9ad928SBarry Smith } else { 1679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Not using diagonal scaling\n")); 1684b9ad928SBarry Smith } 1694b9ad928SBarry Smith } 1704b9ad928SBarry Smith PetscFunctionReturn(0); 1714b9ad928SBarry Smith } 1724b9ad928SBarry Smith 1739371c9d4SSatish Balay static PetscErrorCode PCSetUp_Eisenstat(PC pc) { 17413f74950SBarry Smith PetscInt M, N, m, n; 17552d511d9SJose E. Roman PetscBool set, sym; 1764b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1774b9ad928SBarry Smith 1784b9ad928SBarry Smith PetscFunctionBegin; 1794b9ad928SBarry Smith if (!pc->setupcalled) { 1809566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &M, &N)); 1819566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m, &n)); 18252d511d9SJose E. Roman PetscCall(MatIsSymmetricKnown(pc->mat, &set, &sym)); 1839566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &eis->shell)); 1849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(eis->shell, m, n, M, N)); 1859566063dSJacob Faibussowitsch PetscCall(MatSetType(eis->shell, MATSHELL)); 1869566063dSJacob Faibussowitsch PetscCall(MatSetUp(eis->shell)); 1879566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(eis->shell, pc)); 1889566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT, (void (*)(void))PCMult_Eisenstat)); 18952d511d9SJose E. Roman if (set && sym) PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT_TRANSPOSE, (void (*)(void))PCMult_Eisenstat)); 190a9e2c90aSJose E. Roman PetscCall(MatShellSetOperation(eis->shell, MATOP_NORM, (void (*)(void))PCNorm_Eisenstat)); 1914b9ad928SBarry Smith } 1924b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 193*4dfa11a4SJacob Faibussowitsch if (!pc->setupcalled) { PetscCall(MatCreateVecs(pc->pmat, &eis->diag, NULL)); } 1949566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, eis->diag)); 1954b9ad928SBarry Smith PetscFunctionReturn(0); 1964b9ad928SBarry Smith } 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 1994b9ad928SBarry Smith 2009371c9d4SSatish Balay static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc, PetscReal omega) { 201c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith PetscFunctionBegin; 2042472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range"); 2054b9ad928SBarry Smith eis->omega = omega; 2064b9ad928SBarry Smith PetscFunctionReturn(0); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith 2099371c9d4SSatish Balay static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc, PetscBool flg) { 210c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 2114b9ad928SBarry Smith 2124b9ad928SBarry Smith PetscFunctionBegin; 213c60c7ad4SBarry Smith eis->usediag = flg; 214c60c7ad4SBarry Smith PetscFunctionReturn(0); 215c60c7ad4SBarry Smith } 216c60c7ad4SBarry Smith 2179371c9d4SSatish Balay static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc, PetscReal *omega) { 218c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 219c60c7ad4SBarry Smith 220c60c7ad4SBarry Smith PetscFunctionBegin; 221c60c7ad4SBarry Smith *omega = eis->omega; 222c60c7ad4SBarry Smith PetscFunctionReturn(0); 223c60c7ad4SBarry Smith } 224c60c7ad4SBarry Smith 2259371c9d4SSatish Balay static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc, PetscBool *flg) { 226c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 227c60c7ad4SBarry Smith 228c60c7ad4SBarry Smith PetscFunctionBegin; 229c60c7ad4SBarry Smith *flg = eis->usediag; 2304b9ad928SBarry Smith PetscFunctionReturn(0); 2314b9ad928SBarry Smith } 2324b9ad928SBarry Smith 2334b9ad928SBarry Smith /*@ 2344b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 235f1580f4eSBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default) 2364b9ad928SBarry Smith 237f1580f4eSBarry Smith Logically Collective on pc 2384b9ad928SBarry Smith 2394b9ad928SBarry Smith Input Parameters: 2404b9ad928SBarry Smith + pc - the preconditioner context 2414b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2424b9ad928SBarry Smith 2434b9ad928SBarry Smith Options Database Key: 2444b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith Notes: 2474b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2484b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2494b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2504b9ad928SBarry Smith and right preconditioning. 2514b9ad928SBarry Smith 252f1580f4eSBarry Smith To use SSOR without the Eisenstat trick, employ the `PCSOR` preconditioner, 2534b9ad928SBarry Smith which can be chosen with the database options 2544b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith Level: intermediate 2574b9ad928SBarry Smith 258f1580f4eSBarry Smith .seealso: `PCSORSetOmega()`, `PCEISENSTAT` 2594b9ad928SBarry Smith @*/ 2609371c9d4SSatish Balay PetscErrorCode PCEisenstatSetOmega(PC pc, PetscReal omega) { 2614b9ad928SBarry Smith PetscFunctionBegin; 2620700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 263c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, omega, 2); 264cac4c232SBarry Smith PetscTryMethod(pc, "PCEisenstatSetOmega_C", (PC, PetscReal), (pc, omega)); 2654b9ad928SBarry Smith PetscFunctionReturn(0); 2664b9ad928SBarry Smith } 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith /*@ 269f1580f4eSBarry Smith PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner, `PCEISENSTAT` 2704b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2714b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2724b9ad928SBarry Smith 273f1580f4eSBarry Smith Logically Collective on pc 2744b9ad928SBarry Smith 275c60c7ad4SBarry Smith Input Parameters: 276c60c7ad4SBarry Smith + pc - the preconditioner context 277f1580f4eSBarry Smith - flg - `PETSC_TRUE` turns off diagonal scaling inside the algorithm 2784b9ad928SBarry Smith 2794b9ad928SBarry Smith Options Database Key: 280f1580f4eSBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 2814b9ad928SBarry Smith 2824b9ad928SBarry Smith Level: intermediate 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith Note: 285f1580f4eSBarry Smith If you use the `KSPSetDiagonalScaling()` or -ksp_diagonal_scale option then you will 2866aad120cSJose E. Roman likely want to use this routine since it will save you some unneeded flops. 2874b9ad928SBarry Smith 288f1580f4eSBarry Smith .seealso: `PCEisenstatSetOmega()`, `PCEISENSTAT` 2894b9ad928SBarry Smith @*/ 2909371c9d4SSatish Balay PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc, PetscBool flg) { 2914b9ad928SBarry Smith PetscFunctionBegin; 2920700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 293cac4c232SBarry Smith PetscTryMethod(pc, "PCEisenstatSetNoDiagonalScaling_C", (PC, PetscBool), (pc, flg)); 294c60c7ad4SBarry Smith PetscFunctionReturn(0); 295c60c7ad4SBarry Smith } 296c60c7ad4SBarry Smith 297c60c7ad4SBarry Smith /*@ 298c60c7ad4SBarry Smith PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega, 299c60c7ad4SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 300c60c7ad4SBarry Smith 301f1580f4eSBarry Smith Logically Collective on pc 302c60c7ad4SBarry Smith 303c60c7ad4SBarry Smith Input Parameter: 304c60c7ad4SBarry Smith . pc - the preconditioner context 305c60c7ad4SBarry Smith 306c60c7ad4SBarry Smith Output Parameter: 307c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2) 308c60c7ad4SBarry Smith 309c60c7ad4SBarry Smith Options Database Key: 310c60c7ad4SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 311c60c7ad4SBarry Smith 312c60c7ad4SBarry Smith Notes: 313c60c7ad4SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 314c60c7ad4SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 315c60c7ad4SBarry Smith however, the preconditioned problem must be solved with both left 316c60c7ad4SBarry Smith and right preconditioning. 317c60c7ad4SBarry Smith 318c60c7ad4SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 319c60c7ad4SBarry Smith which can be chosen with the database options 320c60c7ad4SBarry Smith $ -pc_type sor -pc_sor_symmetric 321c60c7ad4SBarry Smith 322c60c7ad4SBarry Smith Level: intermediate 323c60c7ad4SBarry Smith 324f1580f4eSBarry Smith .seealso: `PCEISENSTAT`, `PCSORGetOmega()`, `PCEisenstatSetOmega()` 325c60c7ad4SBarry Smith @*/ 3269371c9d4SSatish Balay PetscErrorCode PCEisenstatGetOmega(PC pc, PetscReal *omega) { 327c60c7ad4SBarry Smith PetscFunctionBegin; 328c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 329cac4c232SBarry Smith PetscUseMethod(pc, "PCEisenstatGetOmega_C", (PC, PetscReal *), (pc, omega)); 330c60c7ad4SBarry Smith PetscFunctionReturn(0); 331c60c7ad4SBarry Smith } 332c60c7ad4SBarry Smith 333c60c7ad4SBarry Smith /*@ 334163d334eSBarry Smith PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner 335c60c7ad4SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 336c60c7ad4SBarry Smith along the diagonal, this may save a small amount of work. 337c60c7ad4SBarry Smith 338f1580f4eSBarry Smith Logically Collective on pc 339c60c7ad4SBarry Smith 340c60c7ad4SBarry Smith Input Parameter: 341c60c7ad4SBarry Smith . pc - the preconditioner context 342c60c7ad4SBarry Smith 343c60c7ad4SBarry Smith Output Parameter: 344f1580f4eSBarry Smith . flg - `PETSC_TRUE` means there is no diagonal scaling applied 345c60c7ad4SBarry Smith 346c60c7ad4SBarry Smith Options Database Key: 347f1580f4eSBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 348c60c7ad4SBarry Smith 349c60c7ad4SBarry Smith Level: intermediate 350c60c7ad4SBarry Smith 351c60c7ad4SBarry Smith Note: 352f9ff08acSPierre Jolivet If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will 3536aad120cSJose E. Roman likely want to use this routine since it will save you some unneeded flops. 354c60c7ad4SBarry Smith 355f1580f4eSBarry Smith .seealso: , `PCEISENSTAT`, `PCEisenstatGetOmega()` 356c60c7ad4SBarry Smith @*/ 3579371c9d4SSatish Balay PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc, PetscBool *flg) { 358c60c7ad4SBarry Smith PetscFunctionBegin; 359c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 360cac4c232SBarry Smith PetscUseMethod(pc, "PCEisenstatGetNoDiagonalScaling_C", (PC, PetscBool *), (pc, flg)); 3614b9ad928SBarry Smith PetscFunctionReturn(0); 3624b9ad928SBarry Smith } 3634b9ad928SBarry Smith 3649371c9d4SSatish Balay static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool *change) { 3658066bbecSBarry Smith PetscFunctionBegin; 3668066bbecSBarry Smith *change = PETSC_TRUE; 3678066bbecSBarry Smith PetscFunctionReturn(0); 3688066bbecSBarry Smith } 3698066bbecSBarry Smith 3704b9ad928SBarry Smith /*MC 3714b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 3724b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 3734b9ad928SBarry Smith 3744b9ad928SBarry Smith Options Database Keys: 3754b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 376f1580f4eSBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith Level: beginner 3794b9ad928SBarry Smith 38095452b02SPatrick Sanan Notes: 381f1580f4eSBarry Smith Only implemented for the `MATAIJ` matrix format. 3824b9ad928SBarry Smith 383f1580f4eSBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block Jacobi with SOR on each block. 384f1580f4eSBarry Smith 385f1580f4eSBarry Smith Developer Note: 386f1580f4eSBarry Smith Since this algorithm runs the Krylov method on a transformed linear system the implementation provides `PCPreSolve()` and `PCPostSolve()` 387f1580f4eSBarry Smith routines that `KSP` uses to set up the transformed linear system. 388f1580f4eSBarry Smith 389f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCEisenstatGetOmega()`, 390db781477SPatrick Sanan `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR` 3914b9ad928SBarry Smith M*/ 3924b9ad928SBarry Smith 3939371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) { 3944b9ad928SBarry Smith PC_Eisenstat *eis; 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith PetscFunctionBegin; 397*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&eis)); 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 40052d511d9SJose E. Roman pc->ops->applytranspose = PCApplyTranspose_Eisenstat; 401dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 402dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 4030a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 4044b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 4054b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 40669d2c0f9SBarry Smith pc->ops->reset = PCReset_Eisenstat; 4074b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 4084b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 4094b9ad928SBarry Smith 4103ec1f749SStefano Zampini pc->data = eis; 4114b9ad928SBarry Smith eis->omega = 1.0; 4120a545947SLisandro Dalcin eis->b[0] = NULL; 4130a545947SLisandro Dalcin eis->b[1] = NULL; 4140a545947SLisandro Dalcin eis->diag = NULL; 4154b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 4164b9ad928SBarry Smith 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", PCEisenstatSetOmega_Eisenstat)); 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", PCEisenstatSetNoDiagonalScaling_Eisenstat)); 4199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", PCEisenstatGetOmega_Eisenstat)); 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", PCEisenstatGetNoDiagonalScaling_Eisenstat)); 4219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Eisenstat)); 4224b9ad928SBarry Smith PetscFunctionReturn(0); 4234b9ad928SBarry Smith } 424