xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision c60c7ad4eb3baa38c2a8deee7f23b741b9653612)
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);
713bb1ff40SBarry 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;
1408afaa268SBarry Smith   PetscBool      set,flg;
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);
145*c60c7ad4SBarry Smith   ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatSetNoDiagonalScaling",eis->usediag ? PETSC_FALSE : PETSC_TRUE,&flg,&set);CHKERRQ(ierr);
146*c60c7ad4SBarry Smith   if (set) {
147*c60c7ad4SBarry Smith     ierr = PCEisenstatSetNoDiagonalScaling(pc,flg);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) {
16457622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",(double)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);
1913bb1ff40SBarry 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) {
1962a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr);
1973bb1ff40SBarry 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 {
209*c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
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->omega = omega;
2144b9ad928SBarry Smith   PetscFunctionReturn(0);
2154b9ad928SBarry Smith }
2164b9ad928SBarry Smith 
2174b9ad928SBarry Smith #undef __FUNCT__
218*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatSetNoDiagonalScaling_Eisenstat"
219*c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc,PetscBool flg)
2204b9ad928SBarry Smith {
221*c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith   PetscFunctionBegin;
224*c60c7ad4SBarry Smith   eis->usediag = flg;
225*c60c7ad4SBarry Smith   PetscFunctionReturn(0);
226*c60c7ad4SBarry Smith }
227*c60c7ad4SBarry Smith 
228*c60c7ad4SBarry Smith #undef __FUNCT__
229*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatGetOmega_Eisenstat"
230*c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatGetOmega_Eisenstat(PC pc,PetscReal *omega)
231*c60c7ad4SBarry Smith {
232*c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
233*c60c7ad4SBarry Smith 
234*c60c7ad4SBarry Smith   PetscFunctionBegin;
235*c60c7ad4SBarry Smith   *omega = eis->omega;
236*c60c7ad4SBarry Smith   PetscFunctionReturn(0);
237*c60c7ad4SBarry Smith }
238*c60c7ad4SBarry Smith 
239*c60c7ad4SBarry Smith #undef __FUNCT__
240*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatGetNoDiagonalScaling_Eisenstat"
241*c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc,PetscBool *flg)
242*c60c7ad4SBarry Smith {
243*c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
244*c60c7ad4SBarry Smith 
245*c60c7ad4SBarry Smith   PetscFunctionBegin;
246*c60c7ad4SBarry Smith   *flg = eis->usediag;
2474b9ad928SBarry Smith   PetscFunctionReturn(0);
2484b9ad928SBarry Smith }
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith #undef __FUNCT__
2514b9ad928SBarry Smith #define __FUNCT__ "PCEisenstatSetOmega"
2524b9ad928SBarry Smith /*@
2534b9ad928SBarry Smith    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
2544b9ad928SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
2554b9ad928SBarry Smith 
256ad4df100SBarry Smith    Logically Collective on PC
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith    Input Parameters:
2594b9ad928SBarry Smith +  pc - the preconditioner context
2604b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2)
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith    Options Database Key:
2634b9ad928SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Notes:
2664b9ad928SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
2674b9ad928SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
2684b9ad928SBarry Smith    however, the preconditioned problem must be solved with both left
2694b9ad928SBarry Smith    and right preconditioning.
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
2724b9ad928SBarry Smith    which can be chosen with the database options
2734b9ad928SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith    Level: intermediate
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith .seealso: PCSORSetOmega()
2804b9ad928SBarry Smith @*/
2817087cfbeSBarry Smith PetscErrorCode  PCEisenstatSetOmega(PC pc,PetscReal omega)
2824b9ad928SBarry Smith {
2834ac538c5SBarry Smith   PetscErrorCode ierr;
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith   PetscFunctionBegin;
2860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
287c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
2884ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
2894b9ad928SBarry Smith   PetscFunctionReturn(0);
2904b9ad928SBarry Smith }
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith #undef __FUNCT__
293*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatSetNoDiagonalScaling"
2944b9ad928SBarry Smith /*@
295*c60c7ad4SBarry Smith    PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner
2964b9ad928SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
2974b9ad928SBarry Smith    along the diagonal, this may save a small amount of work.
2984b9ad928SBarry Smith 
299ad4df100SBarry Smith    Logically Collective on PC
3004b9ad928SBarry Smith 
301*c60c7ad4SBarry Smith    Input Parameters:
302*c60c7ad4SBarry Smith +  pc - the preconditioner context
303*c60c7ad4SBarry Smith -  flg - PETSC_TRUE turns off diagonal scaling inside the algorithm
3044b9ad928SBarry Smith 
3054b9ad928SBarry Smith    Options Database Key:
306*c60c7ad4SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith    Level: intermediate
3094b9ad928SBarry Smith 
3104b9ad928SBarry Smith    Note:
3114b9ad928SBarry Smith      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
3124b9ad928SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith .seealso: PCEisenstatSetOmega()
3174b9ad928SBarry Smith @*/
318*c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatSetNoDiagonalScaling(PC pc,PetscBool flg)
3194b9ad928SBarry Smith {
3204ac538c5SBarry Smith   PetscErrorCode ierr;
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith   PetscFunctionBegin;
3230700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
324*c60c7ad4SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatSetNoDiagonalScaling_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
325*c60c7ad4SBarry Smith   PetscFunctionReturn(0);
326*c60c7ad4SBarry Smith }
327*c60c7ad4SBarry Smith 
328*c60c7ad4SBarry Smith #undef __FUNCT__
329*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatGetOmega"
330*c60c7ad4SBarry Smith /*@
331*c60c7ad4SBarry Smith    PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega,
332*c60c7ad4SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
333*c60c7ad4SBarry Smith 
334*c60c7ad4SBarry Smith    Logically Collective on PC
335*c60c7ad4SBarry Smith 
336*c60c7ad4SBarry Smith    Input Parameter:
337*c60c7ad4SBarry Smith .  pc - the preconditioner context
338*c60c7ad4SBarry Smith 
339*c60c7ad4SBarry Smith    Output Parameter:
340*c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2)
341*c60c7ad4SBarry Smith 
342*c60c7ad4SBarry Smith    Options Database Key:
343*c60c7ad4SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
344*c60c7ad4SBarry Smith 
345*c60c7ad4SBarry Smith    Notes:
346*c60c7ad4SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
347*c60c7ad4SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
348*c60c7ad4SBarry Smith    however, the preconditioned problem must be solved with both left
349*c60c7ad4SBarry Smith    and right preconditioning.
350*c60c7ad4SBarry Smith 
351*c60c7ad4SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
352*c60c7ad4SBarry Smith    which can be chosen with the database options
353*c60c7ad4SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
354*c60c7ad4SBarry Smith 
355*c60c7ad4SBarry Smith    Level: intermediate
356*c60c7ad4SBarry Smith 
357*c60c7ad4SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
358*c60c7ad4SBarry Smith 
359*c60c7ad4SBarry Smith .seealso: PCSORGetOmega(), PCEisenstatSetOmega()
360*c60c7ad4SBarry Smith @*/
361*c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatGetOmega(PC pc,PetscReal *omega)
362*c60c7ad4SBarry Smith {
363*c60c7ad4SBarry Smith   PetscErrorCode ierr;
364*c60c7ad4SBarry Smith 
365*c60c7ad4SBarry Smith   PetscFunctionBegin;
366*c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
367*c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCEisenstatGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
368*c60c7ad4SBarry Smith   PetscFunctionReturn(0);
369*c60c7ad4SBarry Smith }
370*c60c7ad4SBarry Smith 
371*c60c7ad4SBarry Smith #undef __FUNCT__
372*c60c7ad4SBarry Smith #define __FUNCT__ "PCEisenstatGetNoDiagonalScaling"
373*c60c7ad4SBarry Smith /*@
374*c60c7ad4SBarry Smith    PCEisenstatGetNoDiagonalScaling - Causes the Eisenstat preconditioner
375*c60c7ad4SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
376*c60c7ad4SBarry Smith    along the diagonal, this may save a small amount of work.
377*c60c7ad4SBarry Smith 
378*c60c7ad4SBarry Smith    Logically Collective on PC
379*c60c7ad4SBarry Smith 
380*c60c7ad4SBarry Smith    Input Parameter:
381*c60c7ad4SBarry Smith .  pc - the preconditioner context
382*c60c7ad4SBarry Smith 
383*c60c7ad4SBarry Smith    Output Parameter:
384*c60c7ad4SBarry Smith .  flg - PETSC_TRUE means there is no diagonal scaling applied
385*c60c7ad4SBarry Smith 
386*c60c7ad4SBarry Smith    Options Database Key:
387*c60c7ad4SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
388*c60c7ad4SBarry Smith 
389*c60c7ad4SBarry Smith    Level: intermediate
390*c60c7ad4SBarry Smith 
391*c60c7ad4SBarry Smith    Note:
392*c60c7ad4SBarry Smith      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
393*c60c7ad4SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
394*c60c7ad4SBarry Smith 
395*c60c7ad4SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
396*c60c7ad4SBarry Smith 
397*c60c7ad4SBarry Smith .seealso: PCEisenstatGetOmega()
398*c60c7ad4SBarry Smith @*/
399*c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatGetNoDiagonalScaling(PC pc,PetscBool *flg)
400*c60c7ad4SBarry Smith {
401*c60c7ad4SBarry Smith   PetscErrorCode ierr;
402*c60c7ad4SBarry Smith 
403*c60c7ad4SBarry Smith   PetscFunctionBegin;
404*c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
405*c60c7ad4SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatGetNoDiagonalScaling_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
4064b9ad928SBarry Smith   PetscFunctionReturn(0);
4074b9ad928SBarry Smith }
4084b9ad928SBarry Smith 
4094b9ad928SBarry Smith /* --------------------------------------------------------------------*/
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith /*MC
4124b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
4134b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith    Options Database Keys:
4164b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
417*c60c7ad4SBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith    Level: beginner
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith   Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith    Notes: Only implemented for the SeqAIJ matrix format.
4244b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
4254b9ad928SBarry Smith           Jacobi with SOR on each block.
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
428*c60c7ad4SBarry Smith            PCEisenstatSetNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
4294b9ad928SBarry Smith M*/
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith #undef __FUNCT__
4324b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Eisenstat"
4338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc)
4344b9ad928SBarry Smith {
435dfbe8321SBarry Smith   PetscErrorCode ierr;
4364b9ad928SBarry Smith   PC_Eisenstat   *eis;
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith   PetscFunctionBegin;
439b00a9115SJed Brown   ierr = PetscNewLog(pc,&eis);CHKERRQ(ierr);
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
442dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
443dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
4444b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
4454b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
4464b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
44769d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Eisenstat;
4484b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
4494b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
4504b9ad928SBarry Smith 
4514b9ad928SBarry Smith   pc->data     = (void*)eis;
4524b9ad928SBarry Smith   eis->omega   = 1.0;
45378c391d7SBarry Smith   eis->b[0]    = 0;
45478c391d7SBarry Smith   eis->b[1]    = 0;
4554b9ad928SBarry Smith   eis->diag    = 0;
4564b9ad928SBarry Smith   eis->usediag = PETSC_TRUE;
4574b9ad928SBarry Smith 
458bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr);
459*c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetNoDiagonalScaling_C",PCEisenstatSetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
460*c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetOmega_C",PCEisenstatGetOmega_Eisenstat);CHKERRQ(ierr);
461*c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetNoDiagonalScaling_C",PCEisenstatGetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
4624b9ad928SBarry Smith   PetscFunctionReturn(0);
4634b9ad928SBarry Smith }
464