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 17d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMult_Eisenstat(Mat mat, Vec b, Vec x) 18d71ae5a4SJacob Faibussowitsch { 194b9ad928SBarry Smith PC pc; 204b9ad928SBarry Smith PC_Eisenstat *eis; 214b9ad928SBarry Smith 224b9ad928SBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &pc)); 244b9ad928SBarry Smith eis = (PC_Eisenstat *)pc->data; 259566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A, b, eis->omega, SOR_EISENSTAT, 0.0, 1, 1, x)); 2654e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284b9ad928SBarry Smith } 294b9ad928SBarry Smith 30d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCNorm_Eisenstat(Mat mat, NormType type, PetscReal *nrm) 31d71ae5a4SJacob Faibussowitsch { 32a9e2c90aSJose E. Roman PC pc; 33a9e2c90aSJose E. Roman PC_Eisenstat *eis; 34a9e2c90aSJose E. Roman 35a9e2c90aSJose E. Roman PetscFunctionBegin; 36a9e2c90aSJose E. Roman PetscCall(MatShellGetContext(mat, &pc)); 37a9e2c90aSJose E. Roman eis = (PC_Eisenstat *)pc->data; 38a9e2c90aSJose E. Roman PetscCall(MatNorm(eis->A, type, nrm)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40a9e2c90aSJose E. Roman } 41a9e2c90aSJose E. Roman 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Eisenstat(PC pc, Vec x, Vec y) 43d71ae5a4SJacob Faibussowitsch { 444b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 45ace3abfcSBarry Smith PetscBool hasop; 464b9ad928SBarry Smith 474b9ad928SBarry Smith PetscFunctionBegin; 4889c6957cSBarry Smith if (eis->usediag) { 499566063dSJacob Faibussowitsch PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop)); 5089c6957cSBarry Smith if (hasop) { 519566063dSJacob Faibussowitsch PetscCall(MatMultDiagonalBlock(pc->pmat, x, y)); 5289c6957cSBarry Smith } else { 539566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y, x, eis->diag)); 5489c6957cSBarry Smith } 559566063dSJacob Faibussowitsch } else PetscCall(VecCopy(x, y)); 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 574b9ad928SBarry Smith } 584b9ad928SBarry Smith 59d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Eisenstat(PC pc, Vec x, Vec y) 60d71ae5a4SJacob Faibussowitsch { 6152d511d9SJose E. Roman PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 6252d511d9SJose E. Roman PetscBool hasop, set, sym; 6352d511d9SJose E. Roman 6452d511d9SJose E. Roman PetscFunctionBegin; 6552d511d9SJose E. Roman PetscCall(MatIsSymmetricKnown(eis->A, &set, &sym)); 6652d511d9SJose E. Roman PetscCheck(set && sym, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Can only apply transpose of Eisenstat if matrix is symmetric"); 6752d511d9SJose E. Roman if (eis->usediag) { 6852d511d9SJose E. Roman PetscCall(MatHasOperation(pc->pmat, MATOP_MULT_DIAGONAL_BLOCK, &hasop)); 6952d511d9SJose E. Roman if (hasop) { 7052d511d9SJose E. Roman PetscCall(MatMultDiagonalBlock(pc->pmat, x, y)); 7152d511d9SJose E. Roman } else { 7252d511d9SJose E. Roman PetscCall(VecPointwiseMult(y, x, eis->diag)); 7352d511d9SJose E. Roman } 7452d511d9SJose E. Roman } else PetscCall(VecCopy(x, y)); 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7652d511d9SJose E. Roman } 7752d511d9SJose E. Roman 78d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) 79d71ae5a4SJacob Faibussowitsch { 804b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 81ace3abfcSBarry Smith PetscBool nonzero; 824b9ad928SBarry Smith 834b9ad928SBarry Smith PetscFunctionBegin; 8478c391d7SBarry Smith if (pc->presolvedone < 2) { 8508401ef6SPierre Jolivet PetscCheck(pc->mat == pc->pmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot have different mat and pmat"); 864b9ad928SBarry Smith /* swap shell matrix and true matrix */ 874b9ad928SBarry Smith eis->A = pc->mat; 884b9ad928SBarry Smith pc->mat = eis->shell; 894b9ad928SBarry Smith } 904b9ad928SBarry Smith 914dfa11a4SJacob Faibussowitsch if (!eis->b[pc->presolvedone - 1]) { PetscCall(VecDuplicate(b, &eis->b[pc->presolvedone - 1])); } 924b9ad928SBarry Smith 934b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 949566063dSJacob Faibussowitsch PetscCall(KSPGetInitialGuessNonzero(ksp, &nonzero)); 954b9ad928SBarry Smith if (nonzero) { 969566063dSJacob Faibussowitsch PetscCall(VecCopy(x, eis->b[pc->presolvedone - 1])); 979566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A, eis->b[pc->presolvedone - 1], eis->omega, SOR_APPLY_UPPER, 0.0, 1, 1, x)); 9854e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 994b9ad928SBarry Smith } 1004b9ad928SBarry Smith 101121471adSBarry Smith /* save true b, other option is to swap pointers */ 1029566063dSJacob Faibussowitsch PetscCall(VecCopy(b, eis->b[pc->presolvedone - 1])); 103121471adSBarry Smith 1045c99c7daSBarry Smith /* modify b by (L + D/omega)^{-1} */ 1059566063dSJacob 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)); 10654e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1084b9ad928SBarry Smith } 1094b9ad928SBarry Smith 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPostSolve_Eisenstat(PC pc, KSP ksp, Vec b, Vec x) 111d71ae5a4SJacob Faibussowitsch { 1124b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith PetscFunctionBegin; 1154b9ad928SBarry Smith /* get back true b */ 1169566063dSJacob Faibussowitsch PetscCall(VecCopy(eis->b[pc->presolvedone], b)); 117121471adSBarry Smith 118121471adSBarry Smith /* modify x by (U + D/omega)^{-1} */ 1199566063dSJacob Faibussowitsch PetscCall(VecCopy(x, eis->b[pc->presolvedone])); 1209566063dSJacob 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)); 12154e23b33SJose E. Roman PetscCall(MatFactorGetError(eis->A, (MatFactorError *)&pc->failedreason)); 1222fa5cd67SKarl Rupp if (!pc->presolvedone) pc->mat = eis->A; 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1244b9ad928SBarry Smith } 1254b9ad928SBarry Smith 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Eisenstat(PC pc) 127d71ae5a4SJacob Faibussowitsch { 1284b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1294b9ad928SBarry Smith 1304b9ad928SBarry Smith PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[0])); 1329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[1])); 1339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&eis->shell)); 1349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->diag)); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13669d2c0f9SBarry Smith } 13769d2c0f9SBarry Smith 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Eisenstat(PC pc) 139d71ae5a4SJacob Faibussowitsch { 14069d2c0f9SBarry Smith PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(PCReset_Eisenstat(pc)); 1422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", NULL)); 1432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", NULL)); 1442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", NULL)); 1452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", NULL)); 1462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith 151d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc, PetscOptionItems *PetscOptionsObject) 152d71ae5a4SJacob Faibussowitsch { 1534b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1548afaa268SBarry Smith PetscBool set, flg; 1554b9ad928SBarry Smith 1564b9ad928SBarry Smith PetscFunctionBegin; 157d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Eisenstat SSOR options"); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_eisenstat_omega", "Relaxation factor 0 < omega < 2", "PCEisenstatSetOmega", eis->omega, &eis->omega, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling", "Do not use standard diagonal scaling", "PCEisenstatSetNoDiagonalScaling", eis->usediag ? PETSC_FALSE : PETSC_TRUE, &flg, &set)); 1601baa6e33SBarry Smith if (set) PetscCall(PCEisenstatSetNoDiagonalScaling(pc, flg)); 161d0609cedSBarry Smith PetscOptionsHeadEnd(); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1634b9ad928SBarry Smith } 1644b9ad928SBarry Smith 165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Eisenstat(PC pc, PetscViewer viewer) 166d71ae5a4SJacob Faibussowitsch { 1674b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 168ace3abfcSBarry Smith PetscBool iascii; 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith PetscFunctionBegin; 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17232077d6dSBarry Smith if (iascii) { 1739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " omega = %g\n", (double)eis->omega)); 1744b9ad928SBarry Smith if (eis->usediag) { 1759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using diagonal scaling (default)\n")); 1764b9ad928SBarry Smith } else { 1779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Not using diagonal scaling\n")); 1784b9ad928SBarry Smith } 1794b9ad928SBarry Smith } 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith 183d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Eisenstat(PC pc) 184d71ae5a4SJacob Faibussowitsch { 18513f74950SBarry Smith PetscInt M, N, m, n; 18652d511d9SJose E. Roman PetscBool set, sym; 1874b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 1884b9ad928SBarry Smith 1894b9ad928SBarry Smith PetscFunctionBegin; 1904b9ad928SBarry Smith if (!pc->setupcalled) { 1919566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &M, &N)); 1929566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m, &n)); 19352d511d9SJose E. Roman PetscCall(MatIsSymmetricKnown(pc->mat, &set, &sym)); 1949566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &eis->shell)); 1959566063dSJacob Faibussowitsch PetscCall(MatSetSizes(eis->shell, m, n, M, N)); 1969566063dSJacob Faibussowitsch PetscCall(MatSetType(eis->shell, MATSHELL)); 1979566063dSJacob Faibussowitsch PetscCall(MatSetUp(eis->shell)); 1989566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(eis->shell, pc)); 1999566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT, (void (*)(void))PCMult_Eisenstat)); 20052d511d9SJose E. Roman if (set && sym) PetscCall(MatShellSetOperation(eis->shell, MATOP_MULT_TRANSPOSE, (void (*)(void))PCMult_Eisenstat)); 201a9e2c90aSJose E. Roman PetscCall(MatShellSetOperation(eis->shell, MATOP_NORM, (void (*)(void))PCNorm_Eisenstat)); 2024b9ad928SBarry Smith } 2033ba16761SJacob Faibussowitsch if (!eis->usediag) PetscFunctionReturn(PETSC_SUCCESS); 2044dfa11a4SJacob Faibussowitsch if (!pc->setupcalled) { PetscCall(MatCreateVecs(pc->pmat, &eis->diag, NULL)); } 2059566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat, eis->diag)); 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 2104b9ad928SBarry Smith 211d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc, PetscReal omega) 212d71ae5a4SJacob Faibussowitsch { 213c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith PetscFunctionBegin; 2162472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range"); 2174b9ad928SBarry Smith eis->omega = omega; 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith 221d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc, PetscBool flg) 222d71ae5a4SJacob Faibussowitsch { 223c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 2244b9ad928SBarry Smith 2254b9ad928SBarry Smith PetscFunctionBegin; 226c60c7ad4SBarry Smith eis->usediag = flg; 2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 228c60c7ad4SBarry Smith } 229c60c7ad4SBarry Smith 230d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc, PetscReal *omega) 231d71ae5a4SJacob Faibussowitsch { 232c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 233c60c7ad4SBarry Smith 234c60c7ad4SBarry Smith PetscFunctionBegin; 235c60c7ad4SBarry Smith *omega = eis->omega; 2363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 237c60c7ad4SBarry Smith } 238c60c7ad4SBarry Smith 239d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc, PetscBool *flg) 240d71ae5a4SJacob Faibussowitsch { 241c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 242c60c7ad4SBarry Smith 243c60c7ad4SBarry Smith PetscFunctionBegin; 244c60c7ad4SBarry Smith *flg = eis->usediag; 2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2464b9ad928SBarry Smith } 2474b9ad928SBarry Smith 2484b9ad928SBarry Smith /*@ 2494b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 250f1580f4eSBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default) 2514b9ad928SBarry Smith 252c3339decSBarry Smith Logically Collective 2534b9ad928SBarry Smith 2544b9ad928SBarry Smith Input Parameters: 2554b9ad928SBarry Smith + pc - the preconditioner context 2564b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith Options Database Key: 2594b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2604b9ad928SBarry Smith 261*2fe279fdSBarry Smith Level: intermediate 262*2fe279fdSBarry Smith 2634b9ad928SBarry Smith Notes: 2644b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2654b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2664b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2674b9ad928SBarry Smith and right preconditioning. 2684b9ad928SBarry Smith 269f1580f4eSBarry Smith To use SSOR without the Eisenstat trick, employ the `PCSOR` preconditioner, 2704b9ad928SBarry Smith which can be chosen with the database options 2714b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2724b9ad928SBarry Smith 273f1580f4eSBarry Smith .seealso: `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 304f1580f4eSBarry Smith .seealso: `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, 336c60c7ad4SBarry Smith which can be chosen with the database options 337c60c7ad4SBarry Smith $ -pc_type sor -pc_sor_symmetric 338c60c7ad4SBarry Smith 339c60c7ad4SBarry Smith Level: intermediate 340c60c7ad4SBarry Smith 341f1580f4eSBarry Smith .seealso: `PCEISENSTAT`, `PCSORGetOmega()`, `PCEisenstatSetOmega()` 342c60c7ad4SBarry Smith @*/ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatGetOmega(PC pc, PetscReal *omega) 344d71ae5a4SJacob Faibussowitsch { 345c60c7ad4SBarry Smith PetscFunctionBegin; 346c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 347cac4c232SBarry Smith PetscUseMethod(pc, "PCEisenstatGetOmega_C", (PC, PetscReal *), (pc, omega)); 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 349c60c7ad4SBarry Smith } 350c60c7ad4SBarry Smith 351c60c7ad4SBarry Smith /*@ 352163d334eSBarry Smith PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner 353c60c7ad4SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 354c60c7ad4SBarry Smith along the diagonal, this may save a small amount of work. 355c60c7ad4SBarry Smith 356c3339decSBarry Smith Logically Collective 357c60c7ad4SBarry Smith 358c60c7ad4SBarry Smith Input Parameter: 359c60c7ad4SBarry Smith . pc - the preconditioner context 360c60c7ad4SBarry Smith 361c60c7ad4SBarry Smith Output Parameter: 362f1580f4eSBarry Smith . flg - `PETSC_TRUE` means there is no diagonal scaling applied 363c60c7ad4SBarry Smith 364c60c7ad4SBarry Smith Options Database Key: 365f1580f4eSBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 366c60c7ad4SBarry Smith 367c60c7ad4SBarry Smith Level: intermediate 368c60c7ad4SBarry Smith 369c60c7ad4SBarry Smith Note: 370f9ff08acSPierre Jolivet If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will 3716aad120cSJose E. Roman likely want to use this routine since it will save you some unneeded flops. 372c60c7ad4SBarry Smith 373f1580f4eSBarry Smith .seealso: , `PCEISENSTAT`, `PCEisenstatGetOmega()` 374c60c7ad4SBarry Smith @*/ 375d71ae5a4SJacob Faibussowitsch PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc, PetscBool *flg) 376d71ae5a4SJacob Faibussowitsch { 377c60c7ad4SBarry Smith PetscFunctionBegin; 378c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 379cac4c232SBarry Smith PetscUseMethod(pc, "PCEisenstatGetNoDiagonalScaling_C", (PC, PetscBool *), (pc, flg)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3814b9ad928SBarry Smith } 3824b9ad928SBarry Smith 383d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool *change) 384d71ae5a4SJacob Faibussowitsch { 3858066bbecSBarry Smith PetscFunctionBegin; 3868066bbecSBarry Smith *change = PETSC_TRUE; 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3888066bbecSBarry Smith } 3898066bbecSBarry Smith 3904b9ad928SBarry Smith /*MC 3914b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 3924b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith Options Database Keys: 3954b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 396f1580f4eSBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates `PCEisenstatSetNoDiagonalScaling()` 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith Level: beginner 3994b9ad928SBarry Smith 40095452b02SPatrick Sanan Notes: 401f1580f4eSBarry Smith Only implemented for the `MATAIJ` matrix format. 4024b9ad928SBarry Smith 403f1580f4eSBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block Jacobi with SOR on each block. 404f1580f4eSBarry Smith 405f1580f4eSBarry Smith Developer Note: 406f1580f4eSBarry Smith Since this algorithm runs the Krylov method on a transformed linear system the implementation provides `PCPreSolve()` and `PCPostSolve()` 407f1580f4eSBarry Smith routines that `KSP` uses to set up the transformed linear system. 408f1580f4eSBarry Smith 409f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCEisenstatGetOmega()`, 410db781477SPatrick Sanan `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR` 4114b9ad928SBarry Smith M*/ 4124b9ad928SBarry Smith 413d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) 414d71ae5a4SJacob Faibussowitsch { 4154b9ad928SBarry Smith PC_Eisenstat *eis; 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith PetscFunctionBegin; 4184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&eis)); 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 42152d511d9SJose E. Roman pc->ops->applytranspose = PCApplyTranspose_Eisenstat; 422dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 423dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 4240a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 4254b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 4264b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 42769d2c0f9SBarry Smith pc->ops->reset = PCReset_Eisenstat; 4284b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 4294b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 4304b9ad928SBarry Smith 4313ec1f749SStefano Zampini pc->data = eis; 4324b9ad928SBarry Smith eis->omega = 1.0; 4330a545947SLisandro Dalcin eis->b[0] = NULL; 4340a545947SLisandro Dalcin eis->b[1] = NULL; 4350a545947SLisandro Dalcin eis->diag = NULL; 4364b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 4374b9ad928SBarry Smith 4389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetOmega_C", PCEisenstatSetOmega_Eisenstat)); 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatSetNoDiagonalScaling_C", PCEisenstatSetNoDiagonalScaling_Eisenstat)); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetOmega_C", PCEisenstatGetOmega_Eisenstat)); 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCEisenstatGetNoDiagonalScaling_C", PCEisenstatGetNoDiagonalScaling_Eisenstat)); 4429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_Eisenstat)); 4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4444b9ad928SBarry Smith } 445