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" 196849ba73SBarry 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" 346849ba73SBarry 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__ 46dc231df0SBarry Smith #define __FUNCT__ "PCPreSolve_Eisenstat" 47dc231df0SBarry Smith static PetscErrorCode PCPreSolve_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); 62*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,eis->b);CHKERRQ(ierr); 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} */ 75dc231df0SBarry Smith ierr = MatRelax(eis->A,b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr); 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 794b9ad928SBarry Smith #undef __FUNCT__ 80dc231df0SBarry Smith #define __FUNCT__ "PCPostSolve_Eisenstat" 81dc231df0SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec x,Vec b) 824b9ad928SBarry Smith { 834b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 84dfbe8321SBarry Smith PetscErrorCode ierr; 854b9ad928SBarry Smith 864b9ad928SBarry Smith PetscFunctionBegin; 87dc231df0SBarry Smith ierr = MatRelax(eis->A,x,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr); 884b9ad928SBarry Smith pc->mat = eis->A; 894b9ad928SBarry Smith /* get back true b */ 904b9ad928SBarry Smith ierr = VecCopy(eis->b,b);CHKERRQ(ierr); 914b9ad928SBarry Smith PetscFunctionReturn(0); 924b9ad928SBarry Smith } 934b9ad928SBarry Smith 944b9ad928SBarry Smith #undef __FUNCT__ 954b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Eisenstat" 966849ba73SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc) 974b9ad928SBarry Smith { 984b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat *)pc->data; 99dfbe8321SBarry Smith PetscErrorCode ierr; 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith PetscFunctionBegin; 1024b9ad928SBarry Smith if (eis->b) {ierr = VecDestroy(eis->b);CHKERRQ(ierr);} 1034b9ad928SBarry Smith if (eis->shell) {ierr = MatDestroy(eis->shell);CHKERRQ(ierr);} 1044b9ad928SBarry Smith if (eis->diag) {ierr = VecDestroy(eis->diag);CHKERRQ(ierr);} 1054b9ad928SBarry Smith ierr = PetscFree(eis);CHKERRQ(ierr); 1064b9ad928SBarry Smith PetscFunctionReturn(0); 1074b9ad928SBarry Smith } 1084b9ad928SBarry Smith 1094b9ad928SBarry Smith #undef __FUNCT__ 1104b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Eisenstat" 1116849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc) 1124b9ad928SBarry Smith { 1134b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 114dfbe8321SBarry Smith PetscErrorCode ierr; 1154b9ad928SBarry Smith PetscTruth flg; 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith PetscFunctionBegin; 1184b9ad928SBarry Smith ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr); 1194b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,0);CHKERRQ(ierr); 1204b9ad928SBarry Smith ierr = PetscOptionsName("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",&flg);CHKERRQ(ierr); 1214b9ad928SBarry Smith if (flg) { 1224b9ad928SBarry Smith ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr); 1234b9ad928SBarry Smith } 1244b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1254b9ad928SBarry Smith PetscFunctionReturn(0); 1264b9ad928SBarry Smith } 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith #undef __FUNCT__ 1294b9ad928SBarry Smith #define __FUNCT__ "PCView_Eisenstat" 1306849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 1314b9ad928SBarry Smith { 1324b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 133dfbe8321SBarry Smith PetscErrorCode ierr; 13432077d6dSBarry Smith PetscTruth iascii; 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith PetscFunctionBegin; 13732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 13832077d6dSBarry Smith if (iascii) { 1394b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",eis->omega);CHKERRQ(ierr); 1404b9ad928SBarry Smith if (eis->usediag) { 1414b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 1424b9ad928SBarry Smith } else { 1434b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 1444b9ad928SBarry Smith } 1454b9ad928SBarry Smith } else { 14679a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name); 1474b9ad928SBarry Smith } 1484b9ad928SBarry Smith PetscFunctionReturn(0); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith #undef __FUNCT__ 1524b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Eisenstat" 1536849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc) 1544b9ad928SBarry Smith { 155dfbe8321SBarry Smith PetscErrorCode ierr; 15613f74950SBarry Smith PetscInt M,N,m,n; 1574b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith PetscFunctionBegin; 1604b9ad928SBarry Smith if (!pc->setupcalled) { 1614b9ad928SBarry Smith ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 1624b9ad928SBarry Smith ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 163f204ca49SKris Buschelman ierr = MatCreate(pc->comm,m,N,M,N,&eis->shell);CHKERRQ(ierr); 164f204ca49SKris Buschelman ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 165f204ca49SKris Buschelman ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 166*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr); 1674b9ad928SBarry Smith ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 1684b9ad928SBarry Smith } 1694b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 1704b9ad928SBarry Smith if (!pc->setupcalled) { 17123ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 172*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr); 1734b9ad928SBarry Smith } 1744b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 1754b9ad928SBarry Smith PetscFunctionReturn(0); 1764b9ad928SBarry Smith } 1774b9ad928SBarry Smith 1784b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith EXTERN_C_BEGIN 1814b9ad928SBarry Smith #undef __FUNCT__ 1824b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat" 183dfbe8321SBarry Smith PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 1844b9ad928SBarry Smith { 1854b9ad928SBarry Smith PC_Eisenstat *eis; 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith PetscFunctionBegin; 1884b9ad928SBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1894b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 1904b9ad928SBarry Smith eis->omega = omega; 1914b9ad928SBarry Smith PetscFunctionReturn(0); 1924b9ad928SBarry Smith } 1934b9ad928SBarry Smith EXTERN_C_END 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith EXTERN_C_BEGIN 1964b9ad928SBarry Smith #undef __FUNCT__ 1974b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat" 198dfbe8321SBarry Smith PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc) 1994b9ad928SBarry Smith { 2004b9ad928SBarry Smith PC_Eisenstat *eis; 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith PetscFunctionBegin; 2034b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2044b9ad928SBarry Smith eis->usediag = PETSC_FALSE; 2054b9ad928SBarry Smith PetscFunctionReturn(0); 2064b9ad928SBarry Smith } 2074b9ad928SBarry Smith EXTERN_C_END 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith #undef __FUNCT__ 2104b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega" 2114b9ad928SBarry Smith /*@ 2124b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 2134b9ad928SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith Collective on PC 2164b9ad928SBarry Smith 2174b9ad928SBarry Smith Input Parameters: 2184b9ad928SBarry Smith + pc - the preconditioner context 2194b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2204b9ad928SBarry Smith 2214b9ad928SBarry Smith Options Database Key: 2224b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith Notes: 2254b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2264b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2274b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2284b9ad928SBarry Smith and right preconditioning. 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 2314b9ad928SBarry Smith which can be chosen with the database options 2324b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith Level: intermediate 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 2374b9ad928SBarry Smith 2384b9ad928SBarry Smith .seealso: PCSORSetOmega() 2394b9ad928SBarry Smith @*/ 240dfbe8321SBarry Smith PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 2414b9ad928SBarry Smith { 242dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith PetscFunctionBegin; 2454482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2464b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr); 2474b9ad928SBarry Smith if (f) { 2484b9ad928SBarry Smith ierr = (*f)(pc,omega);CHKERRQ(ierr); 2494b9ad928SBarry Smith } 2504b9ad928SBarry Smith PetscFunctionReturn(0); 2514b9ad928SBarry Smith } 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith #undef __FUNCT__ 2544b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling" 2554b9ad928SBarry Smith /*@ 2564b9ad928SBarry Smith PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner 2574b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2584b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2594b9ad928SBarry Smith 2604b9ad928SBarry Smith Collective on PC 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith Input Parameter: 2634b9ad928SBarry Smith . pc - the preconditioner context 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith Options Database Key: 2664b9ad928SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith Level: intermediate 2694b9ad928SBarry Smith 2704b9ad928SBarry Smith Note: 2714b9ad928SBarry Smith If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 2724b9ad928SBarry Smith likley want to use this routine since it will save you some unneeded flops. 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith .seealso: PCEisenstatSetOmega() 2774b9ad928SBarry Smith @*/ 278dfbe8321SBarry Smith PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc) 2794b9ad928SBarry Smith { 280dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 2814b9ad928SBarry Smith 2824b9ad928SBarry Smith PetscFunctionBegin; 2834482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2844b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);CHKERRQ(ierr); 2854b9ad928SBarry Smith if (f) { 2864b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 2874b9ad928SBarry Smith } 2884b9ad928SBarry Smith PetscFunctionReturn(0); 2894b9ad928SBarry Smith } 2904b9ad928SBarry Smith 2914b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 2924b9ad928SBarry Smith 2934b9ad928SBarry Smith /*MC 2944b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 2954b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith Options Database Keys: 2984b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 2994b9ad928SBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling() 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith Level: beginner 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith Notes: Only implemented for the SeqAIJ matrix format. 3064b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3074b9ad928SBarry Smith Jacobi with SOR on each block. 3084b9ad928SBarry Smith 3094b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3104b9ad928SBarry Smith PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 3114b9ad928SBarry Smith M*/ 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith EXTERN_C_BEGIN 3144b9ad928SBarry Smith #undef __FUNCT__ 3154b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat" 316dfbe8321SBarry Smith PetscErrorCode PCCreate_Eisenstat(PC pc) 3174b9ad928SBarry Smith { 318dfbe8321SBarry Smith PetscErrorCode ierr; 3194b9ad928SBarry Smith PC_Eisenstat *eis; 3204b9ad928SBarry Smith 3214b9ad928SBarry Smith PetscFunctionBegin; 3224b9ad928SBarry Smith ierr = PetscNew(PC_Eisenstat,&eis);CHKERRQ(ierr); 323*52e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_Eisenstat));CHKERRQ(ierr); 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 326dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 327dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 3284b9ad928SBarry Smith pc->ops->applyrichardson = 0; 3294b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 3304b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 3314b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 3324b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith pc->data = (void*)eis; 3354b9ad928SBarry Smith eis->omega = 1.0; 3364b9ad928SBarry Smith eis->b = 0; 3374b9ad928SBarry Smith eis->diag = 0; 3384b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 3394b9ad928SBarry Smith 3404b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat", 3414b9ad928SBarry Smith PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 3424b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C", 3434b9ad928SBarry Smith "PCEisenstatNoDiagonalScaling_Eisenstat", 3444b9ad928SBarry Smith PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 3454b9ad928SBarry Smith PetscFunctionReturn(0); 3464b9ad928SBarry Smith } 3474b9ad928SBarry Smith EXTERN_C_END 348