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 */ 8b45d2f2cSJed Brown #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 174b9ad928SBarry Smith 184b9ad928SBarry Smith #undef __FUNCT__ 194b9ad928SBarry Smith #define __FUNCT__ "PCMult_Eisenstat" 206849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x) 214b9ad928SBarry Smith { 22dfbe8321SBarry Smith PetscErrorCode ierr; 234b9ad928SBarry Smith PC pc; 244b9ad928SBarry Smith PC_Eisenstat *eis; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 274b9ad928SBarry Smith ierr = MatShellGetContext(mat,(void**)&pc);CHKERRQ(ierr); 284b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2941f059aeSBarry Smith ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr); 304b9ad928SBarry Smith PetscFunctionReturn(0); 314b9ad928SBarry Smith } 324b9ad928SBarry Smith 334b9ad928SBarry Smith #undef __FUNCT__ 344b9ad928SBarry Smith #define __FUNCT__ "PCApply_Eisenstat" 356849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y) 364b9ad928SBarry Smith { 374b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 38dfbe8321SBarry Smith PetscErrorCode ierr; 39ace3abfcSBarry Smith PetscBool hasop; 404b9ad928SBarry Smith 414b9ad928SBarry Smith PetscFunctionBegin; 4289c6957cSBarry Smith if (eis->usediag) { 4389c6957cSBarry Smith ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 4489c6957cSBarry Smith if (hasop) { 4589c6957cSBarry Smith ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr); 4689c6957cSBarry Smith } else { 4789c6957cSBarry Smith ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr); 4889c6957cSBarry Smith } 4989c6957cSBarry Smith } else {ierr = VecCopy(x,y);CHKERRQ(ierr);} 504b9ad928SBarry Smith PetscFunctionReturn(0); 514b9ad928SBarry Smith } 524b9ad928SBarry Smith 534b9ad928SBarry Smith #undef __FUNCT__ 54dc231df0SBarry Smith #define __FUNCT__ "PCPreSolve_Eisenstat" 55946967b8SBarry Smith static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 564b9ad928SBarry Smith { 574b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 58ace3abfcSBarry Smith PetscBool nonzero; 59dfbe8321SBarry Smith PetscErrorCode ierr; 604b9ad928SBarry Smith 614b9ad928SBarry Smith PetscFunctionBegin; 6278c391d7SBarry Smith if (pc->presolvedone < 2) { 63ce94432eSBarry Smith if (pc->mat != pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot have different mat and pmat"); 644b9ad928SBarry Smith /* swap shell matrix and true matrix */ 654b9ad928SBarry Smith eis->A = pc->mat; 664b9ad928SBarry Smith pc->mat = eis->shell; 674b9ad928SBarry Smith } 684b9ad928SBarry Smith 6978c391d7SBarry Smith if (!eis->b[pc->presolvedone-1]) { 7078c391d7SBarry Smith ierr = VecDuplicate(b,&eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 71*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 7278c391d7SBarry Smith } 734b9ad928SBarry Smith 744b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 754b9ad928SBarry Smith ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 764b9ad928SBarry Smith if (nonzero) { 7778c391d7SBarry Smith ierr = VecCopy(x,eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 7878c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr); 794b9ad928SBarry Smith } 804b9ad928SBarry Smith 81121471adSBarry Smith /* save true b, other option is to swap pointers */ 8278c391d7SBarry Smith ierr = VecCopy(b,eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 83121471adSBarry Smith 845c99c7daSBarry Smith /* modify b by (L + D/omega)^{-1} */ 8578c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr); 864b9ad928SBarry Smith PetscFunctionReturn(0); 874b9ad928SBarry Smith } 884b9ad928SBarry Smith 894b9ad928SBarry Smith #undef __FUNCT__ 90dc231df0SBarry Smith #define __FUNCT__ "PCPostSolve_Eisenstat" 91946967b8SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 924b9ad928SBarry Smith { 934b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 94dfbe8321SBarry Smith PetscErrorCode ierr; 954b9ad928SBarry Smith 964b9ad928SBarry Smith PetscFunctionBegin; 974b9ad928SBarry Smith /* get back true b */ 9878c391d7SBarry Smith ierr = VecCopy(eis->b[pc->presolvedone],b);CHKERRQ(ierr); 99121471adSBarry Smith 100121471adSBarry Smith /* modify x by (U + D/omega)^{-1} */ 10178c391d7SBarry Smith ierr = VecCopy(x,eis->b[pc->presolvedone]);CHKERRQ(ierr); 10278c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr); 1032fa5cd67SKarl Rupp if (!pc->presolvedone) pc->mat = eis->A; 1044b9ad928SBarry Smith PetscFunctionReturn(0); 1054b9ad928SBarry Smith } 1064b9ad928SBarry Smith 1074b9ad928SBarry Smith #undef __FUNCT__ 10869d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Eisenstat" 10969d2c0f9SBarry Smith static PetscErrorCode PCReset_Eisenstat(PC pc) 1104b9ad928SBarry Smith { 1114b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 112dfbe8321SBarry Smith PetscErrorCode ierr; 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith PetscFunctionBegin; 11578c391d7SBarry Smith ierr = VecDestroy(&eis->b[0]);CHKERRQ(ierr); 11678c391d7SBarry Smith ierr = VecDestroy(&eis->b[1]);CHKERRQ(ierr); 1176bf464f9SBarry Smith ierr = MatDestroy(&eis->shell);CHKERRQ(ierr); 1186bf464f9SBarry Smith ierr = VecDestroy(&eis->diag);CHKERRQ(ierr); 11969d2c0f9SBarry Smith PetscFunctionReturn(0); 12069d2c0f9SBarry Smith } 12169d2c0f9SBarry Smith 12269d2c0f9SBarry Smith #undef __FUNCT__ 12369d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Eisenstat" 12469d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc) 12569d2c0f9SBarry Smith { 12669d2c0f9SBarry Smith PetscErrorCode ierr; 12769d2c0f9SBarry Smith 12869d2c0f9SBarry Smith PetscFunctionBegin; 1293f34a243SBarry Smith ierr = PCReset_Eisenstat(pc);CHKERRQ(ierr); 130c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1314b9ad928SBarry Smith PetscFunctionReturn(0); 1324b9ad928SBarry Smith } 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith #undef __FUNCT__ 1354b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Eisenstat" 1366849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc) 1374b9ad928SBarry Smith { 1384b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 139dfbe8321SBarry Smith PetscErrorCode ierr; 140ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith PetscFunctionBegin; 1434b9ad928SBarry Smith ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr); 1440298fd71SBarry Smith ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,NULL);CHKERRQ(ierr); 1450298fd71SBarry Smith ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,NULL);CHKERRQ(ierr); 1464b9ad928SBarry Smith if (flg) { 1474b9ad928SBarry Smith ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr); 1484b9ad928SBarry Smith } 1494b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1504b9ad928SBarry Smith PetscFunctionReturn(0); 1514b9ad928SBarry Smith } 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith #undef __FUNCT__ 1544b9ad928SBarry Smith #define __FUNCT__ "PCView_Eisenstat" 1556849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 1564b9ad928SBarry Smith { 1574b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 158dfbe8321SBarry Smith PetscErrorCode ierr; 159ace3abfcSBarry Smith PetscBool iascii; 1604b9ad928SBarry Smith 1614b9ad928SBarry Smith PetscFunctionBegin; 162251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 16332077d6dSBarry Smith if (iascii) { 164a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr); 1654b9ad928SBarry Smith if (eis->usediag) { 1664b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 1674b9ad928SBarry Smith } else { 1684b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 1694b9ad928SBarry Smith } 1704b9ad928SBarry Smith } 1714b9ad928SBarry Smith PetscFunctionReturn(0); 1724b9ad928SBarry Smith } 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith #undef __FUNCT__ 1754b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Eisenstat" 1766849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc) 1774b9ad928SBarry Smith { 178dfbe8321SBarry Smith PetscErrorCode ierr; 17913f74950SBarry Smith PetscInt M,N,m,n; 1804b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith PetscFunctionBegin; 1834b9ad928SBarry Smith if (!pc->setupcalled) { 1844b9ad928SBarry Smith ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 1854b9ad928SBarry Smith ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 186ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)pc),&eis->shell);CHKERRQ(ierr); 18764aae45aSBarry Smith ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr); 188f204ca49SKris Buschelman ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 1897ae8954aSSatish Balay ierr = MatSetUp(eis->shell);CHKERRQ(ierr); 190f204ca49SKris Buschelman ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 191*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->shell);CHKERRQ(ierr); 1924b9ad928SBarry Smith ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void (*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 1934b9ad928SBarry Smith } 1944b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 1954b9ad928SBarry Smith if (!pc->setupcalled) { 19623ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 197*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->diag);CHKERRQ(ierr); 1984b9ad928SBarry Smith } 1994b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 2004b9ad928SBarry Smith PetscFunctionReturn(0); 2014b9ad928SBarry Smith } 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 2044b9ad928SBarry Smith 2054b9ad928SBarry Smith #undef __FUNCT__ 2064b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat" 2071e6b0712SBarry Smith static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 2084b9ad928SBarry Smith { 2094b9ad928SBarry Smith PC_Eisenstat *eis; 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith PetscFunctionBegin; 212ce94432eSBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 2134b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2144b9ad928SBarry Smith eis->omega = omega; 2154b9ad928SBarry Smith PetscFunctionReturn(0); 2164b9ad928SBarry Smith } 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith #undef __FUNCT__ 2194b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat" 2201e6b0712SBarry Smith static PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc) 2214b9ad928SBarry Smith { 2224b9ad928SBarry Smith PC_Eisenstat *eis; 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith PetscFunctionBegin; 2254b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2264b9ad928SBarry Smith eis->usediag = PETSC_FALSE; 2274b9ad928SBarry Smith PetscFunctionReturn(0); 2284b9ad928SBarry Smith } 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith #undef __FUNCT__ 2314b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega" 2324b9ad928SBarry Smith /*@ 2334b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 2344b9ad928SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 2354b9ad928SBarry Smith 236ad4df100SBarry Smith Logically Collective on PC 2374b9ad928SBarry Smith 2384b9ad928SBarry Smith Input Parameters: 2394b9ad928SBarry Smith + pc - the preconditioner context 2404b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith Options Database Key: 2434b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2444b9ad928SBarry Smith 2454b9ad928SBarry Smith Notes: 2464b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2474b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2484b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2494b9ad928SBarry Smith and right preconditioning. 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 2524b9ad928SBarry Smith which can be chosen with the database options 2534b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2544b9ad928SBarry Smith 2554b9ad928SBarry Smith Level: intermediate 2564b9ad928SBarry Smith 2574b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith .seealso: PCSORSetOmega() 2604b9ad928SBarry Smith @*/ 2617087cfbeSBarry Smith PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 2624b9ad928SBarry Smith { 2634ac538c5SBarry Smith PetscErrorCode ierr; 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith PetscFunctionBegin; 2660700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 267c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 2684ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 2694b9ad928SBarry Smith PetscFunctionReturn(0); 2704b9ad928SBarry Smith } 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith #undef __FUNCT__ 2734b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling" 2744b9ad928SBarry Smith /*@ 2754b9ad928SBarry Smith PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner 2764b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2774b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2784b9ad928SBarry Smith 279ad4df100SBarry Smith Logically Collective on PC 2804b9ad928SBarry Smith 2814b9ad928SBarry Smith Input Parameter: 2824b9ad928SBarry Smith . pc - the preconditioner context 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith Options Database Key: 2854b9ad928SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith Level: intermediate 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith Note: 2904b9ad928SBarry Smith If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 2914b9ad928SBarry Smith likley want to use this routine since it will save you some unneeded flops. 2924b9ad928SBarry Smith 2934b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 2944b9ad928SBarry Smith 2954b9ad928SBarry Smith .seealso: PCEisenstatSetOmega() 2964b9ad928SBarry Smith @*/ 2977087cfbeSBarry Smith PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc) 2984b9ad928SBarry Smith { 2994ac538c5SBarry Smith PetscErrorCode ierr; 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith PetscFunctionBegin; 3020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3034ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));CHKERRQ(ierr); 3044b9ad928SBarry Smith PetscFunctionReturn(0); 3054b9ad928SBarry Smith } 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 3084b9ad928SBarry Smith 3094b9ad928SBarry Smith /*MC 3104b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 3114b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith Options Database Keys: 3144b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 3154b9ad928SBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 3164b9ad928SBarry Smith 3174b9ad928SBarry Smith Level: beginner 3184b9ad928SBarry Smith 3194b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 3204b9ad928SBarry Smith 3214b9ad928SBarry Smith Notes: Only implemented for the SeqAIJ matrix format. 3224b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3234b9ad928SBarry Smith Jacobi with SOR on each block. 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3264b9ad928SBarry Smith PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 3274b9ad928SBarry Smith M*/ 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith #undef __FUNCT__ 3304b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat" 3318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) 3324b9ad928SBarry Smith { 333dfbe8321SBarry Smith PetscErrorCode ierr; 3344b9ad928SBarry Smith PC_Eisenstat *eis; 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith PetscFunctionBegin; 33738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr); 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 340dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 341dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 3424b9ad928SBarry Smith pc->ops->applyrichardson = 0; 3434b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 3444b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 34569d2c0f9SBarry Smith pc->ops->reset = PCReset_Eisenstat; 3464b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 3474b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 3484b9ad928SBarry Smith 3494b9ad928SBarry Smith pc->data = (void*)eis; 3504b9ad928SBarry Smith eis->omega = 1.0; 35178c391d7SBarry Smith eis->b[0] = 0; 35278c391d7SBarry Smith eis->b[1] = 0; 3534b9ad928SBarry Smith eis->diag = 0; 3544b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 3554b9ad928SBarry Smith 356bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 357bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 3584b9ad928SBarry Smith PetscFunctionReturn(0); 3594b9ad928SBarry Smith } 360