xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
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 */
8*b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>           /*I "petscpc.h" I*/
94b9ad928SBarry Smith 
104b9ad928SBarry Smith typedef struct {
114b9ad928SBarry Smith   Mat        shell,A;
124b9ad928SBarry Smith   Vec        b,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;
62e7e72b3dSBarry Smith   if (pc->mat != pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot have different mat and pmat");
634b9ad928SBarry Smith 
644b9ad928SBarry Smith   /* swap shell matrix and true matrix */
654b9ad928SBarry Smith   eis->A    = pc->mat;
664b9ad928SBarry Smith   pc->mat   = eis->shell;
674b9ad928SBarry Smith 
684b9ad928SBarry Smith   if (!eis->b) {
694b9ad928SBarry Smith     ierr = VecDuplicate(b,&eis->b);CHKERRQ(ierr);
7052e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->b);CHKERRQ(ierr);
714b9ad928SBarry Smith   }
724b9ad928SBarry Smith 
734b9ad928SBarry Smith 
744b9ad928SBarry Smith   /* if nonzero initial guess, modify x */
754b9ad928SBarry Smith   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
764b9ad928SBarry Smith   if (nonzero) {
77121471adSBarry Smith     ierr = VecCopy(x,eis->b);CHKERRQ(ierr);
78121471adSBarry Smith     ierr = MatSOR(eis->A,eis->b,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 */
82121471adSBarry Smith   ierr = VecCopy(b,eis->b);CHKERRQ(ierr);
83121471adSBarry Smith 
845c99c7daSBarry Smith   /* modify b by (L + D/omega)^{-1} */
85121471adSBarry Smith   ierr =   MatSOR(eis->A,eis->b,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 */
984b9ad928SBarry Smith   ierr = VecCopy(eis->b,b);CHKERRQ(ierr);
99121471adSBarry Smith 
100121471adSBarry Smith   /* modify x by (U + D/omega)^{-1} */
101121471adSBarry Smith   ierr = VecCopy(x,eis->b);CHKERRQ(ierr);
102121471adSBarry Smith   ierr = MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr);
103121471adSBarry Smith   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;
1156bf464f9SBarry Smith   ierr = VecDestroy(&eis->b);CHKERRQ(ierr);
1166bf464f9SBarry Smith   ierr = MatDestroy(&eis->shell);CHKERRQ(ierr);
1176bf464f9SBarry Smith   ierr = VecDestroy(&eis->diag);CHKERRQ(ierr);
11869d2c0f9SBarry Smith   PetscFunctionReturn(0);
11969d2c0f9SBarry Smith }
12069d2c0f9SBarry Smith 
12169d2c0f9SBarry Smith #undef __FUNCT__
12269d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_Eisenstat"
12369d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc)
12469d2c0f9SBarry Smith {
12569d2c0f9SBarry Smith   PetscErrorCode ierr;
12669d2c0f9SBarry Smith 
12769d2c0f9SBarry Smith   PetscFunctionBegin;
1283f34a243SBarry Smith   ierr = PCReset_Eisenstat(pc);CHKERRQ(ierr);
129c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1304b9ad928SBarry Smith   PetscFunctionReturn(0);
1314b9ad928SBarry Smith }
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith #undef __FUNCT__
1344b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Eisenstat"
1356849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
1364b9ad928SBarry Smith {
1374b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
138dfbe8321SBarry Smith   PetscErrorCode ierr;
139ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1404b9ad928SBarry Smith 
1414b9ad928SBarry Smith   PetscFunctionBegin;
1424b9ad928SBarry Smith   ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr);
14390d69ab7SBarry Smith     ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);CHKERRQ(ierr);
144acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1454b9ad928SBarry Smith     if (flg) {
1464b9ad928SBarry Smith       ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr);
1474b9ad928SBarry Smith     }
1484b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1494b9ad928SBarry Smith   PetscFunctionReturn(0);
1504b9ad928SBarry Smith }
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith #undef __FUNCT__
1534b9ad928SBarry Smith #define __FUNCT__ "PCView_Eisenstat"
1546849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
1554b9ad928SBarry Smith {
1564b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
157dfbe8321SBarry Smith   PetscErrorCode ierr;
158ace3abfcSBarry Smith   PetscBool      iascii;
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith   PetscFunctionBegin;
1612692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
16232077d6dSBarry Smith   if (iascii) {
163a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr);
1644b9ad928SBarry Smith     if (eis->usediag) {
1654b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr);
1664b9ad928SBarry Smith     } else {
1674b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr);
1684b9ad928SBarry Smith     }
1694b9ad928SBarry Smith   } else {
17065e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
1714b9ad928SBarry Smith   }
1724b9ad928SBarry Smith   PetscFunctionReturn(0);
1734b9ad928SBarry Smith }
1744b9ad928SBarry Smith 
1754b9ad928SBarry Smith #undef __FUNCT__
1764b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Eisenstat"
1776849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc)
1784b9ad928SBarry Smith {
179dfbe8321SBarry Smith   PetscErrorCode ierr;
18013f74950SBarry Smith   PetscInt       M,N,m,n;
1814b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith   PetscFunctionBegin;
1844b9ad928SBarry Smith   if (!pc->setupcalled) {
1854b9ad928SBarry Smith     ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr);
1864b9ad928SBarry Smith     ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr);
1877adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)pc)->comm,&eis->shell);CHKERRQ(ierr);
18864aae45aSBarry Smith     ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr);
189f204ca49SKris Buschelman     ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr);
1907ae8954aSSatish Balay     ierr = MatSetUp(eis->shell);CHKERRQ(ierr);
191f204ca49SKris Buschelman     ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr);
19252e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr);
1934b9ad928SBarry Smith     ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr);
1944b9ad928SBarry Smith   }
1954b9ad928SBarry Smith   if (!eis->usediag) PetscFunctionReturn(0);
1964b9ad928SBarry Smith   if (!pc->setupcalled) {
19723ce1328SBarry Smith     ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr);
19852e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr);
1994b9ad928SBarry Smith   }
2004b9ad928SBarry Smith   ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr);
2014b9ad928SBarry Smith   PetscFunctionReturn(0);
2024b9ad928SBarry Smith }
2034b9ad928SBarry Smith 
2044b9ad928SBarry Smith /* --------------------------------------------------------------------*/
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith EXTERN_C_BEGIN
2074b9ad928SBarry Smith #undef __FUNCT__
2084b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat"
2097087cfbeSBarry Smith PetscErrorCode  PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
2104b9ad928SBarry Smith {
2114b9ad928SBarry Smith   PC_Eisenstat  *eis;
2124b9ad928SBarry Smith 
2134b9ad928SBarry Smith   PetscFunctionBegin;
214e7e72b3dSBarry Smith   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
2154b9ad928SBarry Smith   eis = (PC_Eisenstat*)pc->data;
2164b9ad928SBarry Smith   eis->omega = omega;
2174b9ad928SBarry Smith   PetscFunctionReturn(0);
2184b9ad928SBarry Smith }
2194b9ad928SBarry Smith EXTERN_C_END
2204b9ad928SBarry Smith 
2214b9ad928SBarry Smith EXTERN_C_BEGIN
2224b9ad928SBarry Smith #undef __FUNCT__
2234b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat"
2247087cfbeSBarry Smith PetscErrorCode  PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
2254b9ad928SBarry Smith {
2264b9ad928SBarry Smith   PC_Eisenstat *eis;
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith   PetscFunctionBegin;
2294b9ad928SBarry Smith   eis = (PC_Eisenstat*)pc->data;
2304b9ad928SBarry Smith   eis->usediag = PETSC_FALSE;
2314b9ad928SBarry Smith   PetscFunctionReturn(0);
2324b9ad928SBarry Smith }
2334b9ad928SBarry Smith EXTERN_C_END
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith #undef __FUNCT__
2364b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega"
2374b9ad928SBarry Smith /*@
2384b9ad928SBarry Smith    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
2394b9ad928SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
2404b9ad928SBarry Smith 
241ad4df100SBarry Smith    Logically Collective on PC
2424b9ad928SBarry Smith 
2434b9ad928SBarry Smith    Input Parameters:
2444b9ad928SBarry Smith +  pc - the preconditioner context
2454b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2)
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith    Options Database Key:
2484b9ad928SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith    Notes:
2514b9ad928SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
2524b9ad928SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
2534b9ad928SBarry Smith    however, the preconditioned problem must be solved with both left
2544b9ad928SBarry Smith    and right preconditioning.
2554b9ad928SBarry Smith 
2564b9ad928SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
2574b9ad928SBarry Smith    which can be chosen with the database options
2584b9ad928SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
2594b9ad928SBarry Smith 
2604b9ad928SBarry Smith    Level: intermediate
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
2634b9ad928SBarry Smith 
2644b9ad928SBarry Smith .seealso: PCSORSetOmega()
2654b9ad928SBarry Smith @*/
2667087cfbeSBarry Smith PetscErrorCode  PCEisenstatSetOmega(PC pc,PetscReal omega)
2674b9ad928SBarry Smith {
2684ac538c5SBarry Smith   PetscErrorCode ierr;
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith   PetscFunctionBegin;
2710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
272c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
2734ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
2744b9ad928SBarry Smith   PetscFunctionReturn(0);
2754b9ad928SBarry Smith }
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith #undef __FUNCT__
2784b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatNoDiagonalScaling"
2794b9ad928SBarry Smith /*@
2804b9ad928SBarry Smith    PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
2814b9ad928SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
2824b9ad928SBarry Smith    along the diagonal, this may save a small amount of work.
2834b9ad928SBarry Smith 
284ad4df100SBarry Smith    Logically Collective on PC
2854b9ad928SBarry Smith 
2864b9ad928SBarry Smith    Input Parameter:
2874b9ad928SBarry Smith .  pc - the preconditioner context
2884b9ad928SBarry Smith 
2894b9ad928SBarry Smith    Options Database Key:
2904b9ad928SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith    Level: intermediate
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith    Note:
2954b9ad928SBarry Smith      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
2964b9ad928SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith .seealso: PCEisenstatSetOmega()
3014b9ad928SBarry Smith @*/
3027087cfbeSBarry Smith PetscErrorCode  PCEisenstatNoDiagonalScaling(PC pc)
3034b9ad928SBarry Smith {
3044ac538c5SBarry Smith   PetscErrorCode ierr;
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith   PetscFunctionBegin;
3070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3084ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));CHKERRQ(ierr);
3094b9ad928SBarry Smith   PetscFunctionReturn(0);
3104b9ad928SBarry Smith }
3114b9ad928SBarry Smith 
3124b9ad928SBarry Smith /* --------------------------------------------------------------------*/
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith /*MC
3154b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
3164b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
3174b9ad928SBarry Smith 
3184b9ad928SBarry Smith    Options Database Keys:
3194b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
3204b9ad928SBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Level: beginner
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith    Notes: Only implemented for the SeqAIJ matrix format.
3274b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
3284b9ad928SBarry Smith           Jacobi with SOR on each block.
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
3314b9ad928SBarry Smith            PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
3324b9ad928SBarry Smith M*/
3334b9ad928SBarry Smith 
3344b9ad928SBarry Smith EXTERN_C_BEGIN
3354b9ad928SBarry Smith #undef __FUNCT__
3364b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat"
3377087cfbeSBarry Smith PetscErrorCode  PCCreate_Eisenstat(PC pc)
3384b9ad928SBarry Smith {
339dfbe8321SBarry Smith   PetscErrorCode ierr;
3404b9ad928SBarry Smith   PC_Eisenstat   *eis;
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith   PetscFunctionBegin;
34338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr);
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
346dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
347dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
3484b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
3494b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
3504b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
35169d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Eisenstat;
3524b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
3534b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
3544b9ad928SBarry Smith 
3554b9ad928SBarry Smith   pc->data           = (void*)eis;
3564b9ad928SBarry Smith   eis->omega         = 1.0;
3574b9ad928SBarry Smith   eis->b             = 0;
3584b9ad928SBarry Smith   eis->diag          = 0;
3594b9ad928SBarry Smith   eis->usediag       = PETSC_TRUE;
3604b9ad928SBarry Smith 
3614b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
3624b9ad928SBarry Smith                     PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr);
3634b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
3644b9ad928SBarry Smith                     "PCEisenstatNoDiagonalScaling_Eisenstat",
3654b9ad928SBarry Smith                     PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
3664b9ad928SBarry Smith  PetscFunctionReturn(0);
3674b9ad928SBarry Smith }
3684b9ad928SBarry Smith EXTERN_C_END
369