xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 */
8af0996ceSBarry Smith #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 
176849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
184b9ad928SBarry Smith {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
204b9ad928SBarry Smith   PC             pc;
214b9ad928SBarry Smith   PC_Eisenstat   *eis;
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   PetscFunctionBegin;
243ec1f749SStefano Zampini   ierr = MatShellGetContext(mat,&pc);CHKERRQ(ierr);
254b9ad928SBarry Smith   eis  = (PC_Eisenstat*)pc->data;
2641f059aeSBarry Smith   ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr);
274b9ad928SBarry Smith   PetscFunctionReturn(0);
284b9ad928SBarry Smith }
294b9ad928SBarry Smith 
306849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
314b9ad928SBarry Smith {
324b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
33dfbe8321SBarry Smith   PetscErrorCode ierr;
34ace3abfcSBarry Smith   PetscBool      hasop;
354b9ad928SBarry Smith 
364b9ad928SBarry Smith   PetscFunctionBegin;
3789c6957cSBarry Smith   if (eis->usediag) {
3889c6957cSBarry Smith     ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
3989c6957cSBarry Smith     if (hasop) {
4089c6957cSBarry Smith       ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr);
4189c6957cSBarry Smith     } else {
4289c6957cSBarry Smith       ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr);
4389c6957cSBarry Smith     }
4489c6957cSBarry Smith   } else {ierr = VecCopy(x,y);CHKERRQ(ierr);}
454b9ad928SBarry Smith   PetscFunctionReturn(0);
464b9ad928SBarry Smith }
474b9ad928SBarry Smith 
48946967b8SBarry Smith static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
494b9ad928SBarry Smith {
504b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
51ace3abfcSBarry Smith   PetscBool      nonzero;
52dfbe8321SBarry Smith   PetscErrorCode ierr;
534b9ad928SBarry Smith 
544b9ad928SBarry Smith   PetscFunctionBegin;
5578c391d7SBarry Smith   if (pc->presolvedone < 2) {
56*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(pc->mat != pc->pmat,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot have different mat and pmat");
574b9ad928SBarry Smith     /* swap shell matrix and true matrix */
584b9ad928SBarry Smith     eis->A  = pc->mat;
594b9ad928SBarry Smith     pc->mat = eis->shell;
604b9ad928SBarry Smith   }
614b9ad928SBarry Smith 
6278c391d7SBarry Smith   if (!eis->b[pc->presolvedone-1]) {
6378c391d7SBarry Smith     ierr = VecDuplicate(b,&eis->b[pc->presolvedone-1]);CHKERRQ(ierr);
643bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->b[pc->presolvedone-1]);CHKERRQ(ierr);
6578c391d7SBarry Smith   }
664b9ad928SBarry Smith 
674b9ad928SBarry Smith   /* if nonzero initial guess, modify x */
684b9ad928SBarry Smith   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
694b9ad928SBarry Smith   if (nonzero) {
7078c391d7SBarry Smith     ierr = VecCopy(x,eis->b[pc->presolvedone-1]);CHKERRQ(ierr);
7178c391d7SBarry Smith     ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr);
724b9ad928SBarry Smith   }
734b9ad928SBarry Smith 
74121471adSBarry Smith   /* save true b, other option is to swap pointers */
7578c391d7SBarry Smith   ierr = VecCopy(b,eis->b[pc->presolvedone-1]);CHKERRQ(ierr);
76121471adSBarry Smith 
775c99c7daSBarry Smith   /* modify b by (L + D/omega)^{-1} */
7878c391d7SBarry 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);
794b9ad928SBarry Smith   PetscFunctionReturn(0);
804b9ad928SBarry Smith }
814b9ad928SBarry Smith 
82946967b8SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
834b9ad928SBarry Smith {
844b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
85dfbe8321SBarry Smith   PetscErrorCode ierr;
864b9ad928SBarry Smith 
874b9ad928SBarry Smith   PetscFunctionBegin;
884b9ad928SBarry Smith   /* get back true b */
8978c391d7SBarry Smith   ierr = VecCopy(eis->b[pc->presolvedone],b);CHKERRQ(ierr);
90121471adSBarry Smith 
91121471adSBarry Smith   /* modify x by (U + D/omega)^{-1} */
9278c391d7SBarry Smith   ierr = VecCopy(x,eis->b[pc->presolvedone]);CHKERRQ(ierr);
9378c391d7SBarry 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);
942fa5cd67SKarl Rupp   if (!pc->presolvedone) pc->mat = eis->A;
954b9ad928SBarry Smith   PetscFunctionReturn(0);
964b9ad928SBarry Smith }
974b9ad928SBarry Smith 
9869d2c0f9SBarry Smith static PetscErrorCode PCReset_Eisenstat(PC pc)
994b9ad928SBarry Smith {
1004b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
101dfbe8321SBarry Smith   PetscErrorCode ierr;
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith   PetscFunctionBegin;
10478c391d7SBarry Smith   ierr = VecDestroy(&eis->b[0]);CHKERRQ(ierr);
10578c391d7SBarry Smith   ierr = VecDestroy(&eis->b[1]);CHKERRQ(ierr);
1066bf464f9SBarry Smith   ierr = MatDestroy(&eis->shell);CHKERRQ(ierr);
1076bf464f9SBarry Smith   ierr = VecDestroy(&eis->diag);CHKERRQ(ierr);
10869d2c0f9SBarry Smith   PetscFunctionReturn(0);
10969d2c0f9SBarry Smith }
11069d2c0f9SBarry Smith 
11169d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc)
11269d2c0f9SBarry Smith {
11369d2c0f9SBarry Smith   PetscErrorCode ierr;
11469d2c0f9SBarry Smith 
11569d2c0f9SBarry Smith   PetscFunctionBegin;
1163f34a243SBarry Smith   ierr = PCReset_Eisenstat(pc);CHKERRQ(ierr);
117c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1184b9ad928SBarry Smith   PetscFunctionReturn(0);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith 
1214416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PetscOptionItems *PetscOptionsObject,PC pc)
1224b9ad928SBarry Smith {
1234b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
124dfbe8321SBarry Smith   PetscErrorCode ierr;
1258afaa268SBarry Smith   PetscBool      set,flg;
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith   PetscFunctionBegin;
128e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Eisenstat SSOR options");CHKERRQ(ierr);
1290298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,NULL);CHKERRQ(ierr);
130c60c7ad4SBarry 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);
131c60c7ad4SBarry Smith   if (set) {
132c60c7ad4SBarry Smith     ierr = PCEisenstatSetNoDiagonalScaling(pc,flg);CHKERRQ(ierr);
1334b9ad928SBarry Smith   }
1344b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1354b9ad928SBarry Smith   PetscFunctionReturn(0);
1364b9ad928SBarry Smith }
1374b9ad928SBarry Smith 
1386849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
1394b9ad928SBarry Smith {
1404b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
141dfbe8321SBarry Smith   PetscErrorCode ierr;
142ace3abfcSBarry Smith   PetscBool      iascii;
1434b9ad928SBarry Smith 
1444b9ad928SBarry Smith   PetscFunctionBegin;
145251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
14632077d6dSBarry Smith   if (iascii) {
147efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  omega = %g\n",(double)eis->omega);CHKERRQ(ierr);
1484b9ad928SBarry Smith     if (eis->usediag) {
149efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Using diagonal scaling (default)\n");CHKERRQ(ierr);
1504b9ad928SBarry Smith     } else {
151efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Not using diagonal scaling\n");CHKERRQ(ierr);
1524b9ad928SBarry Smith     }
1534b9ad928SBarry Smith   }
1544b9ad928SBarry Smith   PetscFunctionReturn(0);
1554b9ad928SBarry Smith }
1564b9ad928SBarry Smith 
1576849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc)
1584b9ad928SBarry Smith {
159dfbe8321SBarry Smith   PetscErrorCode ierr;
16013f74950SBarry Smith   PetscInt       M,N,m,n;
1614b9ad928SBarry Smith   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith   PetscFunctionBegin;
1644b9ad928SBarry Smith   if (!pc->setupcalled) {
1654b9ad928SBarry Smith     ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr);
1664b9ad928SBarry Smith     ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr);
167ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)pc),&eis->shell);CHKERRQ(ierr);
16864aae45aSBarry Smith     ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr);
169f204ca49SKris Buschelman     ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr);
1707ae8954aSSatish Balay     ierr = MatSetUp(eis->shell);CHKERRQ(ierr);
1713ec1f749SStefano Zampini     ierr = MatShellSetContext(eis->shell,pc);CHKERRQ(ierr);
1723bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->shell);CHKERRQ(ierr);
1734b9ad928SBarry Smith     ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void (*)(void))PCMult_Eisenstat);CHKERRQ(ierr);
1744b9ad928SBarry Smith   }
1754b9ad928SBarry Smith   if (!eis->usediag) PetscFunctionReturn(0);
1764b9ad928SBarry Smith   if (!pc->setupcalled) {
1770a545947SLisandro Dalcin     ierr = MatCreateVecs(pc->pmat,&eis->diag,NULL);CHKERRQ(ierr);
1783bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->diag);CHKERRQ(ierr);
1794b9ad928SBarry Smith   }
1804b9ad928SBarry Smith   ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr);
1814b9ad928SBarry Smith   PetscFunctionReturn(0);
1824b9ad928SBarry Smith }
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith /* --------------------------------------------------------------------*/
1854b9ad928SBarry Smith 
1861e6b0712SBarry Smith static PetscErrorCode  PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
1874b9ad928SBarry Smith {
188c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   PetscFunctionBegin;
191*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(omega >= 2.0 || omega <= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
1924b9ad928SBarry Smith   eis->omega = omega;
1934b9ad928SBarry Smith   PetscFunctionReturn(0);
1944b9ad928SBarry Smith }
1954b9ad928SBarry Smith 
196c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc,PetscBool flg)
1974b9ad928SBarry Smith {
198c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith   PetscFunctionBegin;
201c60c7ad4SBarry Smith   eis->usediag = flg;
202c60c7ad4SBarry Smith   PetscFunctionReturn(0);
203c60c7ad4SBarry Smith }
204c60c7ad4SBarry Smith 
205c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatGetOmega_Eisenstat(PC pc,PetscReal *omega)
206c60c7ad4SBarry Smith {
207c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
208c60c7ad4SBarry Smith 
209c60c7ad4SBarry Smith   PetscFunctionBegin;
210c60c7ad4SBarry Smith   *omega = eis->omega;
211c60c7ad4SBarry Smith   PetscFunctionReturn(0);
212c60c7ad4SBarry Smith }
213c60c7ad4SBarry Smith 
214c60c7ad4SBarry Smith static PetscErrorCode  PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc,PetscBool *flg)
215c60c7ad4SBarry Smith {
216c60c7ad4SBarry Smith   PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
217c60c7ad4SBarry Smith 
218c60c7ad4SBarry Smith   PetscFunctionBegin;
219c60c7ad4SBarry Smith   *flg = eis->usediag;
2204b9ad928SBarry Smith   PetscFunctionReturn(0);
2214b9ad928SBarry Smith }
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith /*@
2244b9ad928SBarry Smith    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
2254b9ad928SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
2264b9ad928SBarry Smith 
227ad4df100SBarry Smith    Logically Collective on PC
2284b9ad928SBarry Smith 
2294b9ad928SBarry Smith    Input Parameters:
2304b9ad928SBarry Smith +  pc - the preconditioner context
2314b9ad928SBarry Smith -  omega - relaxation coefficient (0 < omega < 2)
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith    Options Database Key:
2344b9ad928SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
2354b9ad928SBarry Smith 
2364b9ad928SBarry Smith    Notes:
2374b9ad928SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
2384b9ad928SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
2394b9ad928SBarry Smith    however, the preconditioned problem must be solved with both left
2404b9ad928SBarry Smith    and right preconditioning.
2414b9ad928SBarry Smith 
2424b9ad928SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
2434b9ad928SBarry Smith    which can be chosen with the database options
2444b9ad928SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
2454b9ad928SBarry Smith 
2464b9ad928SBarry Smith    Level: intermediate
2474b9ad928SBarry Smith 
2484b9ad928SBarry Smith .seealso: PCSORSetOmega()
2494b9ad928SBarry Smith @*/
2507087cfbeSBarry Smith PetscErrorCode  PCEisenstatSetOmega(PC pc,PetscReal omega)
2514b9ad928SBarry Smith {
2524ac538c5SBarry Smith   PetscErrorCode ierr;
2534b9ad928SBarry Smith 
2544b9ad928SBarry Smith   PetscFunctionBegin;
2550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
256c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(pc,omega,2);
2574ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr);
2584b9ad928SBarry Smith   PetscFunctionReturn(0);
2594b9ad928SBarry Smith }
2604b9ad928SBarry Smith 
2614b9ad928SBarry Smith /*@
262c60c7ad4SBarry Smith    PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner
2634b9ad928SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
2644b9ad928SBarry Smith    along the diagonal, this may save a small amount of work.
2654b9ad928SBarry Smith 
266ad4df100SBarry Smith    Logically Collective on PC
2674b9ad928SBarry Smith 
268c60c7ad4SBarry Smith    Input Parameters:
269c60c7ad4SBarry Smith +  pc - the preconditioner context
270c60c7ad4SBarry Smith -  flg - PETSC_TRUE turns off diagonal scaling inside the algorithm
2714b9ad928SBarry Smith 
2724b9ad928SBarry Smith    Options Database Key:
273c60c7ad4SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith    Level: intermediate
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith    Note:
278f9ff08acSPierre Jolivet      If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will
2794b9ad928SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith .seealso: PCEisenstatSetOmega()
2824b9ad928SBarry Smith @*/
283c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatSetNoDiagonalScaling(PC pc,PetscBool flg)
2844b9ad928SBarry Smith {
2854ac538c5SBarry Smith   PetscErrorCode ierr;
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith   PetscFunctionBegin;
2880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
289c60c7ad4SBarry Smith   ierr = PetscTryMethod(pc,"PCEisenstatSetNoDiagonalScaling_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr);
290c60c7ad4SBarry Smith   PetscFunctionReturn(0);
291c60c7ad4SBarry Smith }
292c60c7ad4SBarry Smith 
293c60c7ad4SBarry Smith /*@
294c60c7ad4SBarry Smith    PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega,
295c60c7ad4SBarry Smith    to use with Eisenstat's trick (where omega = 1.0 by default).
296c60c7ad4SBarry Smith 
297c60c7ad4SBarry Smith    Logically Collective on PC
298c60c7ad4SBarry Smith 
299c60c7ad4SBarry Smith    Input Parameter:
300c60c7ad4SBarry Smith .  pc - the preconditioner context
301c60c7ad4SBarry Smith 
302c60c7ad4SBarry Smith    Output Parameter:
303c60c7ad4SBarry Smith .  omega - relaxation coefficient (0 < omega < 2)
304c60c7ad4SBarry Smith 
305c60c7ad4SBarry Smith    Options Database Key:
306c60c7ad4SBarry Smith .  -pc_eisenstat_omega <omega> - Sets omega
307c60c7ad4SBarry Smith 
308c60c7ad4SBarry Smith    Notes:
309c60c7ad4SBarry Smith    The Eisenstat trick implementation of SSOR requires about 50% of the
310c60c7ad4SBarry Smith    usual amount of floating point operations used for SSOR + Krylov method;
311c60c7ad4SBarry Smith    however, the preconditioned problem must be solved with both left
312c60c7ad4SBarry Smith    and right preconditioning.
313c60c7ad4SBarry Smith 
314c60c7ad4SBarry Smith    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
315c60c7ad4SBarry Smith    which can be chosen with the database options
316c60c7ad4SBarry Smith $    -pc_type  sor  -pc_sor_symmetric
317c60c7ad4SBarry Smith 
318c60c7ad4SBarry Smith    Level: intermediate
319c60c7ad4SBarry Smith 
320c60c7ad4SBarry Smith .seealso: PCSORGetOmega(), PCEisenstatSetOmega()
321c60c7ad4SBarry Smith @*/
322c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatGetOmega(PC pc,PetscReal *omega)
323c60c7ad4SBarry Smith {
324c60c7ad4SBarry Smith   PetscErrorCode ierr;
325c60c7ad4SBarry Smith 
326c60c7ad4SBarry Smith   PetscFunctionBegin;
327c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
328c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCEisenstatGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr);
329c60c7ad4SBarry Smith   PetscFunctionReturn(0);
330c60c7ad4SBarry Smith }
331c60c7ad4SBarry Smith 
332c60c7ad4SBarry Smith /*@
333163d334eSBarry Smith    PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner
334c60c7ad4SBarry Smith    not to do additional diagonal preconditioning. For matrices with a constant
335c60c7ad4SBarry Smith    along the diagonal, this may save a small amount of work.
336c60c7ad4SBarry Smith 
337c60c7ad4SBarry Smith    Logically Collective on PC
338c60c7ad4SBarry Smith 
339c60c7ad4SBarry Smith    Input Parameter:
340c60c7ad4SBarry Smith .  pc - the preconditioner context
341c60c7ad4SBarry Smith 
342c60c7ad4SBarry Smith    Output Parameter:
343c60c7ad4SBarry Smith .  flg - PETSC_TRUE means there is no diagonal scaling applied
344c60c7ad4SBarry Smith 
345c60c7ad4SBarry Smith    Options Database Key:
346c60c7ad4SBarry Smith .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
347c60c7ad4SBarry Smith 
348c60c7ad4SBarry Smith    Level: intermediate
349c60c7ad4SBarry Smith 
350c60c7ad4SBarry Smith    Note:
351f9ff08acSPierre Jolivet      If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will
352c60c7ad4SBarry Smith    likley want to use this routine since it will save you some unneeded flops.
353c60c7ad4SBarry Smith 
354c60c7ad4SBarry Smith .seealso: PCEisenstatGetOmega()
355c60c7ad4SBarry Smith @*/
356c60c7ad4SBarry Smith PetscErrorCode  PCEisenstatGetNoDiagonalScaling(PC pc,PetscBool *flg)
357c60c7ad4SBarry Smith {
358c60c7ad4SBarry Smith   PetscErrorCode ierr;
359c60c7ad4SBarry Smith 
360c60c7ad4SBarry Smith   PetscFunctionBegin;
361c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
362163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCEisenstatGetNoDiagonalScaling_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr);
3634b9ad928SBarry Smith   PetscFunctionReturn(0);
3644b9ad928SBarry Smith }
3654b9ad928SBarry Smith 
3668066bbecSBarry Smith static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool* change)
3678066bbecSBarry Smith {
3688066bbecSBarry Smith   PetscFunctionBegin;
3698066bbecSBarry Smith   *change = PETSC_TRUE;
3708066bbecSBarry Smith   PetscFunctionReturn(0);
3718066bbecSBarry Smith }
3728066bbecSBarry Smith 
3734b9ad928SBarry Smith /* --------------------------------------------------------------------*/
3744b9ad928SBarry Smith 
3754b9ad928SBarry Smith /*MC
3764b9ad928SBarry Smith      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
3774b9ad928SBarry Smith            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
3784b9ad928SBarry Smith 
3794b9ad928SBarry Smith    Options Database Keys:
3804b9ad928SBarry Smith +  -pc_eisenstat_omega <omega> - Sets omega
381c60c7ad4SBarry Smith -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling()
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith    Level: beginner
3844b9ad928SBarry Smith 
38595452b02SPatrick Sanan    Notes:
38695452b02SPatrick Sanan     Only implemented for the SeqAIJ matrix format.
3874b9ad928SBarry Smith           Not a true parallel SOR, in parallel this implementation corresponds to block
3884b9ad928SBarry Smith           Jacobi with SOR on each block.
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
391c60c7ad4SBarry Smith            PCEisenstatSetNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
3924b9ad928SBarry Smith M*/
3934b9ad928SBarry Smith 
3948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc)
3954b9ad928SBarry Smith {
396dfbe8321SBarry Smith   PetscErrorCode ierr;
3974b9ad928SBarry Smith   PC_Eisenstat   *eis;
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith   PetscFunctionBegin;
400b00a9115SJed Brown   ierr = PetscNewLog(pc,&eis);CHKERRQ(ierr);
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith   pc->ops->apply           = PCApply_Eisenstat;
403dc231df0SBarry Smith   pc->ops->presolve        = PCPreSolve_Eisenstat;
404dc231df0SBarry Smith   pc->ops->postsolve       = PCPostSolve_Eisenstat;
4050a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
4064b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
4074b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Eisenstat;
40869d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Eisenstat;
4094b9ad928SBarry Smith   pc->ops->view            = PCView_Eisenstat;
4104b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Eisenstat;
4114b9ad928SBarry Smith 
4123ec1f749SStefano Zampini   pc->data     = eis;
4134b9ad928SBarry Smith   eis->omega   = 1.0;
4140a545947SLisandro Dalcin   eis->b[0]    = NULL;
4150a545947SLisandro Dalcin   eis->b[1]    = NULL;
4160a545947SLisandro Dalcin   eis->diag    = NULL;
4174b9ad928SBarry Smith   eis->usediag = PETSC_TRUE;
4184b9ad928SBarry Smith 
419bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr);
420c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetNoDiagonalScaling_C",PCEisenstatSetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
421c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetOmega_C",PCEisenstatGetOmega_Eisenstat);CHKERRQ(ierr);
422c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetNoDiagonalScaling_C",PCEisenstatGetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
4238066bbecSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Eisenstat);CHKERRQ(ierr);
4244b9ad928SBarry Smith   PetscFunctionReturn(0);
4254b9ad928SBarry Smith }
426