14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Defines a Eisenstat trick SSOR preconditioner. This uses about 34b9ad928SBarry Smith %50 of the usual amount of floating point ops used for SSOR + Krylov 44b9ad928SBarry Smith method. But it requires actually solving the preconditioned problem 54b9ad928SBarry Smith with both left and right preconditioning. 64b9ad928SBarry Smith */ 74b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 84b9ad928SBarry Smith 94b9ad928SBarry Smith typedef struct { 104b9ad928SBarry Smith Mat shell,A; 114b9ad928SBarry Smith Vec b,diag; /* temporary storage for true right hand side */ 124b9ad928SBarry Smith PetscReal omega; 134b9ad928SBarry Smith PetscTruth usediag; /* indicates preconditioner should include diagonal scaling*/ 144b9ad928SBarry Smith } PC_Eisenstat; 154b9ad928SBarry Smith 164b9ad928SBarry Smith 174b9ad928SBarry Smith #undef __FUNCT__ 184b9ad928SBarry Smith #define __FUNCT__ "PCMult_Eisenstat" 19*6849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x) 204b9ad928SBarry Smith { 21dfbe8321SBarry Smith PetscErrorCode ierr; 224b9ad928SBarry Smith PC pc; 234b9ad928SBarry Smith PC_Eisenstat *eis; 244b9ad928SBarry Smith 254b9ad928SBarry Smith PetscFunctionBegin; 264b9ad928SBarry Smith ierr = MatShellGetContext(mat,(void **)&pc);CHKERRQ(ierr); 274b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 284b9ad928SBarry Smith ierr = MatRelax(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr); 294b9ad928SBarry Smith PetscFunctionReturn(0); 304b9ad928SBarry Smith } 314b9ad928SBarry Smith 324b9ad928SBarry Smith #undef __FUNCT__ 334b9ad928SBarry Smith #define __FUNCT__ "PCApply_Eisenstat" 34*6849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y) 354b9ad928SBarry Smith { 364b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 37dfbe8321SBarry Smith PetscErrorCode ierr; 384b9ad928SBarry Smith 394b9ad928SBarry Smith PetscFunctionBegin; 404b9ad928SBarry Smith if (eis->usediag) {ierr = VecPointwiseMult(x,eis->diag,y);CHKERRQ(ierr);} 414b9ad928SBarry Smith else {ierr = VecCopy(x,y);CHKERRQ(ierr);} 424b9ad928SBarry Smith PetscFunctionReturn(0); 434b9ad928SBarry Smith } 444b9ad928SBarry Smith 454b9ad928SBarry Smith #undef __FUNCT__ 464b9ad928SBarry Smith #define __FUNCT__ "PCPre_Eisenstat" 47*6849ba73SBarry Smith static PetscErrorCode PCPre_Eisenstat(PC pc,KSP ksp,Vec x,Vec b) 484b9ad928SBarry Smith { 494b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 504b9ad928SBarry Smith PetscTruth nonzero; 51dfbe8321SBarry Smith PetscErrorCode ierr; 524b9ad928SBarry Smith 534b9ad928SBarry Smith PetscFunctionBegin; 544b9ad928SBarry Smith if (pc->mat != pc->pmat) SETERRQ(PETSC_ERR_SUP,"Cannot have different mat and pmat"); 554b9ad928SBarry Smith 564b9ad928SBarry Smith /* swap shell matrix and true matrix */ 574b9ad928SBarry Smith eis->A = pc->mat; 584b9ad928SBarry Smith pc->mat = eis->shell; 594b9ad928SBarry Smith 604b9ad928SBarry Smith if (!eis->b) { 614b9ad928SBarry Smith ierr = VecDuplicate(b,&eis->b);CHKERRQ(ierr); 624b9ad928SBarry Smith PetscLogObjectParent(pc,eis->b); 634b9ad928SBarry Smith } 644b9ad928SBarry Smith 654b9ad928SBarry Smith /* save true b, other option is to swap pointers */ 664b9ad928SBarry Smith ierr = VecCopy(b,eis->b);CHKERRQ(ierr); 674b9ad928SBarry Smith 684b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 694b9ad928SBarry Smith ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 704b9ad928SBarry Smith if (nonzero) { 714b9ad928SBarry Smith ierr = MatRelax(eis->A,x,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr); 724b9ad928SBarry Smith } 734b9ad928SBarry Smith 744b9ad928SBarry Smith /* modify b by (L + D)^{-1} */ 754b9ad928SBarry Smith ierr = MatRelax(eis->A,b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | 764b9ad928SBarry Smith SOR_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr); 774b9ad928SBarry Smith PetscFunctionReturn(0); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith 804b9ad928SBarry Smith #undef __FUNCT__ 814b9ad928SBarry Smith #define __FUNCT__ "PCPost_Eisenstat" 82*6849ba73SBarry Smith static PetscErrorCode PCPost_Eisenstat(PC pc,KSP ksp,Vec x,Vec b) 834b9ad928SBarry Smith { 844b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 85dfbe8321SBarry Smith PetscErrorCode ierr; 864b9ad928SBarry Smith 874b9ad928SBarry Smith PetscFunctionBegin; 884b9ad928SBarry Smith ierr = MatRelax(eis->A,x,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | 894b9ad928SBarry Smith SOR_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr); 904b9ad928SBarry Smith pc->mat = eis->A; 914b9ad928SBarry Smith /* get back true b */ 924b9ad928SBarry Smith ierr = VecCopy(eis->b,b);CHKERRQ(ierr); 934b9ad928SBarry Smith PetscFunctionReturn(0); 944b9ad928SBarry Smith } 954b9ad928SBarry Smith 964b9ad928SBarry Smith #undef __FUNCT__ 974b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Eisenstat" 98*6849ba73SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc) 994b9ad928SBarry Smith { 1004b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 101dfbe8321SBarry Smith PetscErrorCode ierr; 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith PetscFunctionBegin; 1044b9ad928SBarry Smith if (eis->b) {ierr = VecDestroy(eis->b);CHKERRQ(ierr);} 1054b9ad928SBarry Smith if (eis->shell) {ierr = MatDestroy(eis->shell);CHKERRQ(ierr);} 1064b9ad928SBarry Smith if (eis->diag) {ierr = VecDestroy(eis->diag);CHKERRQ(ierr);} 1074b9ad928SBarry Smith ierr = PetscFree(eis);CHKERRQ(ierr); 1084b9ad928SBarry Smith PetscFunctionReturn(0); 1094b9ad928SBarry Smith } 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith #undef __FUNCT__ 1124b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Eisenstat" 113*6849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc) 1144b9ad928SBarry Smith { 1154b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 116dfbe8321SBarry Smith PetscErrorCode ierr; 1174b9ad928SBarry Smith PetscTruth flg; 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith PetscFunctionBegin; 1204b9ad928SBarry Smith ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr); 1214b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,0);CHKERRQ(ierr); 1224b9ad928SBarry Smith ierr = PetscOptionsName("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",&flg);CHKERRQ(ierr); 1234b9ad928SBarry Smith if (flg) { 1244b9ad928SBarry Smith ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr); 1254b9ad928SBarry Smith } 1264b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1274b9ad928SBarry Smith PetscFunctionReturn(0); 1284b9ad928SBarry Smith } 1294b9ad928SBarry Smith 1304b9ad928SBarry Smith #undef __FUNCT__ 1314b9ad928SBarry Smith #define __FUNCT__ "PCView_Eisenstat" 132*6849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 1334b9ad928SBarry Smith { 1344b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 135dfbe8321SBarry Smith PetscErrorCode ierr; 13632077d6dSBarry Smith PetscTruth iascii; 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith PetscFunctionBegin; 13932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 14032077d6dSBarry Smith if (iascii) { 1414b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",eis->omega);CHKERRQ(ierr); 1424b9ad928SBarry Smith if (eis->usediag) { 1434b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 1444b9ad928SBarry Smith } else { 1454b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 1464b9ad928SBarry Smith } 1474b9ad928SBarry Smith } else { 1484b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith PetscFunctionReturn(0); 1514b9ad928SBarry Smith } 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith #undef __FUNCT__ 1544b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Eisenstat" 155*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc) 1564b9ad928SBarry Smith { 157dfbe8321SBarry Smith PetscErrorCode ierr; 158dfbe8321SBarry Smith int M,N,m,n; 1594b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1604b9ad928SBarry Smith 1614b9ad928SBarry Smith PetscFunctionBegin; 1624b9ad928SBarry Smith if (!pc->setupcalled) { 1634b9ad928SBarry Smith ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 1644b9ad928SBarry Smith ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 165f204ca49SKris Buschelman ierr = MatCreate(pc->comm,m,N,M,N,&eis->shell);CHKERRQ(ierr); 166f204ca49SKris Buschelman ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 167f204ca49SKris Buschelman ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 1684b9ad928SBarry Smith PetscLogObjectParent(pc,eis->shell); 1694b9ad928SBarry Smith ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 1704b9ad928SBarry Smith } 1714b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 1724b9ad928SBarry Smith if (!pc->setupcalled) { 17323ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 1744b9ad928SBarry Smith PetscLogObjectParent(pc,eis->diag); 1754b9ad928SBarry Smith } 1764b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 1774b9ad928SBarry Smith PetscFunctionReturn(0); 1784b9ad928SBarry Smith } 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith EXTERN_C_BEGIN 1834b9ad928SBarry Smith #undef __FUNCT__ 1844b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat" 185dfbe8321SBarry Smith PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 1864b9ad928SBarry Smith { 1874b9ad928SBarry Smith PC_Eisenstat *eis; 1884b9ad928SBarry Smith 1894b9ad928SBarry Smith PetscFunctionBegin; 1904b9ad928SBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1914b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 1924b9ad928SBarry Smith eis->omega = omega; 1934b9ad928SBarry Smith PetscFunctionReturn(0); 1944b9ad928SBarry Smith } 1954b9ad928SBarry Smith EXTERN_C_END 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith EXTERN_C_BEGIN 1984b9ad928SBarry Smith #undef __FUNCT__ 1994b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat" 200dfbe8321SBarry Smith PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc) 2014b9ad928SBarry Smith { 2024b9ad928SBarry Smith PC_Eisenstat *eis; 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith PetscFunctionBegin; 2054b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2064b9ad928SBarry Smith eis->usediag = PETSC_FALSE; 2074b9ad928SBarry Smith PetscFunctionReturn(0); 2084b9ad928SBarry Smith } 2094b9ad928SBarry Smith EXTERN_C_END 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith #undef __FUNCT__ 2124b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega" 2134b9ad928SBarry Smith /*@ 2144b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 2154b9ad928SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 2164b9ad928SBarry Smith 2174b9ad928SBarry Smith Collective on PC 2184b9ad928SBarry Smith 2194b9ad928SBarry Smith Input Parameters: 2204b9ad928SBarry Smith + pc - the preconditioner context 2214b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith Options Database Key: 2244b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Notes: 2274b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2284b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2294b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2304b9ad928SBarry Smith and right preconditioning. 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 2334b9ad928SBarry Smith which can be chosen with the database options 2344b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith Level: intermediate 2374b9ad928SBarry Smith 2384b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 2394b9ad928SBarry Smith 2404b9ad928SBarry Smith .seealso: PCSORSetOmega() 2414b9ad928SBarry Smith @*/ 242dfbe8321SBarry Smith PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 2434b9ad928SBarry Smith { 244dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith PetscFunctionBegin; 2474482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2484b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr); 2494b9ad928SBarry Smith if (f) { 2504b9ad928SBarry Smith ierr = (*f)(pc,omega);CHKERRQ(ierr); 2514b9ad928SBarry Smith } 2524b9ad928SBarry Smith PetscFunctionReturn(0); 2534b9ad928SBarry Smith } 2544b9ad928SBarry Smith 2554b9ad928SBarry Smith #undef __FUNCT__ 2564b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling" 2574b9ad928SBarry Smith /*@ 2584b9ad928SBarry Smith PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner 2594b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2604b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith Collective on PC 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith Input Parameter: 2654b9ad928SBarry Smith . pc - the preconditioner context 2664b9ad928SBarry Smith 2674b9ad928SBarry Smith Options Database Key: 2684b9ad928SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 2694b9ad928SBarry Smith 2704b9ad928SBarry Smith Level: intermediate 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith Note: 2734b9ad928SBarry Smith If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 2744b9ad928SBarry Smith likley want to use this routine since it will save you some unneeded flops. 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith .seealso: PCEisenstatSetOmega() 2794b9ad928SBarry Smith @*/ 280dfbe8321SBarry Smith PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc) 2814b9ad928SBarry Smith { 282dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith PetscFunctionBegin; 2854482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2864b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);CHKERRQ(ierr); 2874b9ad928SBarry Smith if (f) { 2884b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 2894b9ad928SBarry Smith } 2904b9ad928SBarry Smith PetscFunctionReturn(0); 2914b9ad928SBarry Smith } 2924b9ad928SBarry Smith 2934b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 2944b9ad928SBarry Smith 2954b9ad928SBarry Smith /*MC 2964b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 2974b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith Options Database Keys: 3004b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 3014b9ad928SBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Level: beginner 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith Notes: Only implemented for the SeqAIJ matrix format. 3084b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3094b9ad928SBarry Smith Jacobi with SOR on each block. 3104b9ad928SBarry Smith 3114b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3124b9ad928SBarry Smith PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 3134b9ad928SBarry Smith M*/ 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith EXTERN_C_BEGIN 3164b9ad928SBarry Smith #undef __FUNCT__ 3174b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat" 318dfbe8321SBarry Smith PetscErrorCode PCCreate_Eisenstat(PC pc) 3194b9ad928SBarry Smith { 320dfbe8321SBarry Smith PetscErrorCode ierr; 3214b9ad928SBarry Smith PC_Eisenstat *eis; 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith PetscFunctionBegin; 3244b9ad928SBarry Smith ierr = PetscNew(PC_Eisenstat,&eis);CHKERRQ(ierr); 3254b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Eisenstat)); 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 3284b9ad928SBarry Smith pc->ops->presolve = PCPre_Eisenstat; 3294b9ad928SBarry Smith pc->ops->postsolve = PCPost_Eisenstat; 3304b9ad928SBarry Smith pc->ops->applyrichardson = 0; 3314b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 3324b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 3334b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 3344b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith pc->data = (void*)eis; 3374b9ad928SBarry Smith eis->omega = 1.0; 3384b9ad928SBarry Smith eis->b = 0; 3394b9ad928SBarry Smith eis->diag = 0; 3404b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 3414b9ad928SBarry Smith 3424b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat", 3434b9ad928SBarry Smith PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 3444b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C", 3454b9ad928SBarry Smith "PCEisenstatNoDiagonalScaling_Eisenstat", 3464b9ad928SBarry Smith PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 3474b9ad928SBarry Smith PetscFunctionReturn(0); 3484b9ad928SBarry Smith } 3494b9ad928SBarry Smith EXTERN_C_END 350