xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith    Defines a  Eisenstat trick SSOR  preconditioner. This uses about
54b9ad928SBarry Smith  %50 of the usual amount of floating point ops used for SSOR + Krylov
64b9ad928SBarry Smith  method. But it requires actually solving the preconditioned problem
74b9ad928SBarry Smith  with both left and right preconditioning.
84b9ad928SBarry Smith */
96356e834SBarry Smith #include "private/pcimpl.h"           /*I "petscpc.h" I*/
104b9ad928SBarry Smith 
114b9ad928SBarry Smith typedef struct {
124b9ad928SBarry Smith   Mat        shell,A;
134b9ad928SBarry Smith   Vec        b,diag;     /* temporary storage for true right hand side */
144b9ad928SBarry Smith   PetscReal  omega;
154b9ad928SBarry Smith   PetscTruth usediag;    /* indicates preconditioner should include diagonal scaling*/
164b9ad928SBarry Smith } PC_Eisenstat;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith 
194b9ad928SBarry Smith #undef __FUNCT__
204b9ad928SBarry Smith #define __FUNCT__ "PCMult_Eisenstat"
216849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
224b9ad928SBarry Smith {
23dfbe8321SBarry Smith   PetscErrorCode ierr;
244b9ad928SBarry Smith   PC             pc;
254b9ad928SBarry Smith   PC_Eisenstat   *eis;
264b9ad928SBarry Smith 
274b9ad928SBarry Smith   PetscFunctionBegin;
284b9ad928SBarry Smith   ierr = MatShellGetContext(mat,(void **)&pc);CHKERRQ(ierr);
294b9ad928SBarry Smith   eis = (PC_Eisenstat*)pc->data;
3041f059aeSBarry Smith   ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr);
314b9ad928SBarry Smith   PetscFunctionReturn(0);
324b9ad928SBarry Smith }
334b9ad928SBarry Smith 
344b9ad928SBarry Smith #undef __FUNCT__
354b9ad928SBarry Smith #define __FUNCT__ "PCApply_Eisenstat"
366849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
374b9ad928SBarry Smith {
384b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
39dfbe8321SBarry Smith   PetscErrorCode ierr;
4089c6957cSBarry Smith   PetscTruth     hasop;
414b9ad928SBarry Smith 
424b9ad928SBarry Smith   PetscFunctionBegin;
4389c6957cSBarry Smith   if (eis->usediag)  {
4489c6957cSBarry Smith     ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
4589c6957cSBarry Smith     if (hasop) {
4689c6957cSBarry Smith       ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr);
4789c6957cSBarry Smith     } else {
4889c6957cSBarry Smith       ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr);
4989c6957cSBarry Smith     }
5089c6957cSBarry Smith   } else {ierr = VecCopy(x,y);CHKERRQ(ierr);}
514b9ad928SBarry Smith   PetscFunctionReturn(0);
524b9ad928SBarry Smith }
534b9ad928SBarry Smith 
544b9ad928SBarry Smith #undef __FUNCT__
55dc231df0SBarry Smith #define __FUNCT__ "PCPreSolve_Eisenstat"
56946967b8SBarry Smith static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
574b9ad928SBarry Smith {
584b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
594b9ad928SBarry Smith   PetscTruth     nonzero;
60dfbe8321SBarry Smith   PetscErrorCode ierr;
614b9ad928SBarry Smith 
624b9ad928SBarry Smith   PetscFunctionBegin;
634b9ad928SBarry Smith   if (pc->mat != pc->pmat) SETERRQ(PETSC_ERR_SUP,"Cannot have different mat and pmat");
644b9ad928SBarry Smith 
654b9ad928SBarry Smith   /* swap shell matrix and true matrix */
664b9ad928SBarry Smith   eis->A    = pc->mat;
674b9ad928SBarry Smith   pc->mat   = eis->shell;
684b9ad928SBarry Smith 
694b9ad928SBarry Smith   if (!eis->b) {
704b9ad928SBarry Smith     ierr = VecDuplicate(b,&eis->b);CHKERRQ(ierr);
7152e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->b);CHKERRQ(ierr);
724b9ad928SBarry Smith   }
734b9ad928SBarry Smith 
744b9ad928SBarry Smith   /* save true b, other option is to swap pointers */
754b9ad928SBarry Smith   ierr = VecCopy(b,eis->b);CHKERRQ(ierr);
764b9ad928SBarry Smith 
774b9ad928SBarry Smith   /* if nonzero initial guess, modify x */
784b9ad928SBarry Smith   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
794b9ad928SBarry Smith   if (nonzero) {
8041f059aeSBarry Smith     ierr = MatSOR(eis->A,x,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr);
814b9ad928SBarry Smith   }
824b9ad928SBarry Smith 
835c99c7daSBarry Smith   /* modify b by (L + D/omega)^{-1} */
8441f059aeSBarry Smith   ierr =   MatSOR(eis->A,b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr);
854b9ad928SBarry Smith   PetscFunctionReturn(0);
864b9ad928SBarry Smith }
874b9ad928SBarry Smith 
884b9ad928SBarry Smith #undef __FUNCT__
89dc231df0SBarry Smith #define __FUNCT__ "PCPostSolve_Eisenstat"
90946967b8SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
914b9ad928SBarry Smith {
924b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
93dfbe8321SBarry Smith   PetscErrorCode ierr;
944b9ad928SBarry Smith 
954b9ad928SBarry Smith   PetscFunctionBegin;
965c99c7daSBarry Smith   /* modify x by (U + D/omega)^{-1} */
9741f059aeSBarry Smith   ierr =   MatSOR(eis->A,x,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr);
984b9ad928SBarry Smith   pc->mat = eis->A;
994b9ad928SBarry Smith   /* get back true b */
1004b9ad928SBarry Smith   ierr = VecCopy(eis->b,b);CHKERRQ(ierr);
1014b9ad928SBarry Smith   PetscFunctionReturn(0);
1024b9ad928SBarry Smith }
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith #undef __FUNCT__
1054b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Eisenstat"
1066849ba73SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc)
1074b9ad928SBarry Smith {
1084b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat *)pc->data;
109dfbe8321SBarry Smith   PetscErrorCode ierr;
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith   PetscFunctionBegin;
1124b9ad928SBarry Smith   if (eis->b)     {ierr = VecDestroy(eis->b);CHKERRQ(ierr);}
1134b9ad928SBarry Smith   if (eis->shell) {ierr = MatDestroy(eis->shell);CHKERRQ(ierr);}
1144b9ad928SBarry Smith   if (eis->diag)  {ierr = VecDestroy(eis->diag);CHKERRQ(ierr);}
1154b9ad928SBarry Smith   ierr = PetscFree(eis);CHKERRQ(ierr);
1164b9ad928SBarry Smith   PetscFunctionReturn(0);
1174b9ad928SBarry Smith }
1184b9ad928SBarry Smith 
1194b9ad928SBarry Smith #undef __FUNCT__
1204b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Eisenstat"
1216849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
1224b9ad928SBarry Smith {
1234b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
124dfbe8321SBarry Smith   PetscErrorCode ierr;
12590d69ab7SBarry Smith   PetscTruth     flg = PETSC_FALSE;
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith   PetscFunctionBegin;
1284b9ad928SBarry Smith   ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr);
12990d69ab7SBarry Smith     ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);CHKERRQ(ierr);
13090d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1314b9ad928SBarry Smith     if (flg) {
1324b9ad928SBarry Smith       ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr);
1334b9ad928SBarry Smith     }
1344b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1354b9ad928SBarry Smith   PetscFunctionReturn(0);
1364b9ad928SBarry Smith }
1374b9ad928SBarry Smith 
1384b9ad928SBarry Smith #undef __FUNCT__
1394b9ad928SBarry Smith #define __FUNCT__ "PCView_Eisenstat"
1406849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
1414b9ad928SBarry Smith {
1424b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
143dfbe8321SBarry Smith   PetscErrorCode ierr;
14432077d6dSBarry Smith   PetscTruth     iascii;
1454b9ad928SBarry Smith 
1464b9ad928SBarry Smith   PetscFunctionBegin;
14732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
14832077d6dSBarry Smith   if (iascii) {
149a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr);
1504b9ad928SBarry Smith     if (eis->usediag) {
1514b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr);
1524b9ad928SBarry Smith     } else {
1534b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr);
1544b9ad928SBarry Smith     }
1554b9ad928SBarry Smith   } else {
15679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
1574b9ad928SBarry Smith   }
1584b9ad928SBarry Smith   PetscFunctionReturn(0);
1594b9ad928SBarry Smith }
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith #undef __FUNCT__
1624b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Eisenstat"
1636849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc)
1644b9ad928SBarry Smith {
165dfbe8321SBarry Smith   PetscErrorCode ierr;
16613f74950SBarry Smith   PetscInt       M,N,m,n;
1674b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
1684b9ad928SBarry Smith 
1694b9ad928SBarry Smith   PetscFunctionBegin;
1704b9ad928SBarry Smith   if (!pc->setupcalled) {
1714b9ad928SBarry Smith     ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr);
1724b9ad928SBarry Smith     ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr);
1737adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)pc)->comm,&eis->shell);CHKERRQ(ierr);
17464aae45aSBarry Smith     ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr);
175f204ca49SKris Buschelman     ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr);
176f204ca49SKris Buschelman     ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr);
17752e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr);
1784b9ad928SBarry Smith     ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr);
1794b9ad928SBarry Smith   }
1804b9ad928SBarry Smith   if (!eis->usediag) PetscFunctionReturn(0);
1814b9ad928SBarry Smith   if (!pc->setupcalled) {
18223ce1328SBarry Smith     ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr);
18352e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr);
1844b9ad928SBarry Smith   }
1854b9ad928SBarry Smith   ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr);
1864b9ad928SBarry Smith   PetscFunctionReturn(0);
1874b9ad928SBarry Smith }
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith /* --------------------------------------------------------------------*/
1904b9ad928SBarry Smith 
1914b9ad928SBarry Smith EXTERN_C_BEGIN
1924b9ad928SBarry Smith #undef __FUNCT__
1934b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat"
194dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
1954b9ad928SBarry Smith {
1964b9ad928SBarry Smith   PC_Eisenstat  *eis;
1974b9ad928SBarry Smith 
1984b9ad928SBarry Smith   PetscFunctionBegin;
1994b9ad928SBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
2004b9ad928SBarry Smith   eis = (PC_Eisenstat*)pc->data;
2014b9ad928SBarry Smith   eis->omega = omega;
2024b9ad928SBarry Smith   PetscFunctionReturn(0);
2034b9ad928SBarry Smith }
2044b9ad928SBarry Smith EXTERN_C_END
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith EXTERN_C_BEGIN
2074b9ad928SBarry Smith #undef __FUNCT__
2084b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat"
209dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
2104b9ad928SBarry Smith {
2114b9ad928SBarry Smith   PC_Eisenstat *eis;
2124b9ad928SBarry Smith 
2134b9ad928SBarry Smith   PetscFunctionBegin;
2144b9ad928SBarry Smith   eis = (PC_Eisenstat*)pc->data;
2154b9ad928SBarry Smith   eis->usediag = PETSC_FALSE;
2164b9ad928SBarry Smith   PetscFunctionReturn(0);
2174b9ad928SBarry Smith }
2184b9ad928SBarry Smith EXTERN_C_END
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith #undef __FUNCT__
2214b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega"
2224b9ad928SBarry Smith /*@
2234b9ad928SBarry Smith    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
2244b9ad928SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
2254b9ad928SBarry Smith 
2264b9ad928SBarry Smith    Collective on PC
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith    Input Parameters:
2294b9ad928SBarry Smith +  pc - the preconditioner context
2304b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2)
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith    Options Database Key:
2334b9ad928SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith    Notes:
2364b9ad928SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
2374b9ad928SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
2384b9ad928SBarry Smith    however, the preconditioned problem must be solved with both left
2394b9ad928SBarry Smith    and right preconditioning.
2404b9ad928SBarry Smith 
2414b9ad928SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
2424b9ad928SBarry Smith    which can be chosen with the database options
2434b9ad928SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith    Level: intermediate
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith .seealso: PCSORSetOmega()
2504b9ad928SBarry Smith @*/
251dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC pc,PetscReal omega)
2524b9ad928SBarry Smith {
253dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscReal);
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith   PetscFunctionBegin;
256*0700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2574b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr);
2584b9ad928SBarry Smith   if (f) {
2594b9ad928SBarry Smith     ierr = (*f)(pc,omega);CHKERRQ(ierr);
2604b9ad928SBarry Smith   }
2614b9ad928SBarry Smith   PetscFunctionReturn(0);
2624b9ad928SBarry Smith }
2634b9ad928SBarry Smith 
2644b9ad928SBarry Smith #undef __FUNCT__
2654b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling"
2664b9ad928SBarry Smith /*@
2674b9ad928SBarry Smith    PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
2684b9ad928SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
2694b9ad928SBarry Smith    along the diagonal, this may save a small amount of work.
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Collective on PC
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith    Input Parameter:
2744b9ad928SBarry Smith .  pc - the preconditioner context
2754b9ad928SBarry Smith 
2764b9ad928SBarry Smith    Options Database Key:
2774b9ad928SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith    Level: intermediate
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith    Note:
2824b9ad928SBarry Smith      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
2834b9ad928SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith .seealso: PCEisenstatSetOmega()
2884b9ad928SBarry Smith @*/
289dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC pc)
2904b9ad928SBarry Smith {
291dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC);
2924b9ad928SBarry Smith 
2934b9ad928SBarry Smith   PetscFunctionBegin;
294*0700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2954b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);CHKERRQ(ierr);
2964b9ad928SBarry Smith   if (f) {
2974b9ad928SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
2984b9ad928SBarry Smith   }
2994b9ad928SBarry Smith   PetscFunctionReturn(0);
3004b9ad928SBarry Smith }
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith /* --------------------------------------------------------------------*/
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith /*MC
3054b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
3064b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith    Options Database Keys:
3094b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
3104b9ad928SBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
3114b9ad928SBarry Smith 
3124b9ad928SBarry Smith    Level: beginner
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith    Notes: Only implemented for the SeqAIJ matrix format.
3174b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
3184b9ad928SBarry Smith           Jacobi with SOR on each block.
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
3214b9ad928SBarry Smith            PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
3224b9ad928SBarry Smith M*/
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith EXTERN_C_BEGIN
3254b9ad928SBarry Smith #undef __FUNCT__
3264b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat"
327dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Eisenstat(PC pc)
3284b9ad928SBarry Smith {
329dfbe8321SBarry Smith   PetscErrorCode ierr;
3304b9ad928SBarry Smith   PC_Eisenstat   *eis;
3314b9ad928SBarry Smith 
3324b9ad928SBarry Smith   PetscFunctionBegin;
33338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr);
3344b9ad928SBarry Smith 
3354b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
336dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
337dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
3384b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
3394b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
3404b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
3414b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
3424b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith   pc->data           = (void*)eis;
3454b9ad928SBarry Smith   eis->omega         = 1.0;
3464b9ad928SBarry Smith   eis->b             = 0;
3474b9ad928SBarry Smith   eis->diag          = 0;
3484b9ad928SBarry Smith   eis->usediag       = PETSC_TRUE;
3494b9ad928SBarry Smith 
3504b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
3514b9ad928SBarry Smith                     PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr);
3524b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
3534b9ad928SBarry Smith                     "PCEisenstatNoDiagonalScaling_Eisenstat",
3544b9ad928SBarry Smith                     PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
3554b9ad928SBarry Smith  PetscFunctionReturn(0);
3564b9ad928SBarry Smith }
3574b9ad928SBarry Smith EXTERN_C_END
358