xref: /petsc/src/ksp/pc/impls/eisens/eisen.c (revision f3fe499b4cc4d64bf04aa4f5e4963dcc4eb56541)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a  Eisenstat trick SSOR  preconditioner. This uses about
5  %50 of the usual amount of floating point ops used for SSOR + Krylov
6  method. But it requires actually solving the preconditioned problem
7  with both left and right preconditioning.
8 */
9 #include "private/pcimpl.h"           /*I "petscpc.h" I*/
10 
11 typedef struct {
12   Mat        shell,A;
13   Vec        b,diag;     /* temporary storage for true right hand side */
14   PetscReal  omega;
15   PetscTruth usediag;    /* indicates preconditioner should include diagonal scaling*/
16 } PC_Eisenstat;
17 
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "PCMult_Eisenstat"
21 static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
22 {
23   PetscErrorCode ierr;
24   PC             pc;
25   PC_Eisenstat   *eis;
26 
27   PetscFunctionBegin;
28   ierr = MatShellGetContext(mat,(void **)&pc);CHKERRQ(ierr);
29   eis = (PC_Eisenstat*)pc->data;
30   ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 #undef __FUNCT__
35 #define __FUNCT__ "PCApply_Eisenstat"
36 static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
37 {
38   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
39   PetscErrorCode ierr;
40   PetscTruth     hasop;
41 
42   PetscFunctionBegin;
43   if (eis->usediag)  {
44     ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
45     if (hasop) {
46       ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr);
47     } else {
48       ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr);
49     }
50   } else {ierr = VecCopy(x,y);CHKERRQ(ierr);}
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "PCPreSolve_Eisenstat"
56 static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
57 {
58   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
59   PetscTruth     nonzero;
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   if (pc->mat != pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot have different mat and pmat");
64 
65   /* swap shell matrix and true matrix */
66   eis->A    = pc->mat;
67   pc->mat   = eis->shell;
68 
69   if (!eis->b) {
70     ierr = VecDuplicate(b,&eis->b);CHKERRQ(ierr);
71     ierr = PetscLogObjectParent(pc,eis->b);CHKERRQ(ierr);
72   }
73 
74 
75   /* if nonzero initial guess, modify x */
76   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
77   if (nonzero) {
78     ierr = VecCopy(x,eis->b);CHKERRQ(ierr);
79     ierr = MatSOR(eis->A,eis->b,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr);
80   }
81 
82   /* save true b, other option is to swap pointers */
83   ierr = VecCopy(b,eis->b);CHKERRQ(ierr);
84 
85   /* modify b by (L + D/omega)^{-1} */
86   ierr =   MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "PCPostSolve_Eisenstat"
92 static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
93 {
94   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
95   PetscErrorCode ierr;
96 
97   PetscFunctionBegin;
98   /* get back true b */
99   ierr = VecCopy(eis->b,b);CHKERRQ(ierr);
100 
101   /* modify x by (U + D/omega)^{-1} */
102   ierr = VecCopy(x,eis->b);CHKERRQ(ierr);
103   ierr = MatSOR(eis->A,eis->b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr);
104   pc->mat = eis->A;
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "PCDestroy_Eisenstat"
110 static PetscErrorCode PCDestroy_Eisenstat(PC pc)
111 {
112   PC_Eisenstat   *eis = (PC_Eisenstat *)pc->data;
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin;
116   if (eis->b)     {ierr = VecDestroy(eis->b);CHKERRQ(ierr);}
117   if (eis->shell) {ierr = MatDestroy(eis->shell);CHKERRQ(ierr);}
118   if (eis->diag)  {ierr = VecDestroy(eis->diag);CHKERRQ(ierr);}
119   ierr = PetscFree(eis);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "PCSetFromOptions_Eisenstat"
125 static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
126 {
127   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
128   PetscErrorCode ierr;
129   PetscTruth     flg = PETSC_FALSE;
130 
131   PetscFunctionBegin;
132   ierr = PetscOptionsHead("Eisenstat SSOR options");CHKERRQ(ierr);
133     ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,PETSC_NULL);CHKERRQ(ierr);
134     ierr = PetscOptionsTruth("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
135     if (flg) {
136       ierr = PCEisenstatNoDiagonalScaling(pc);CHKERRQ(ierr);
137     }
138   ierr = PetscOptionsTail();CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "PCView_Eisenstat"
144 static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
145 {
146   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
147   PetscErrorCode ierr;
148   PetscTruth     iascii;
149 
150   PetscFunctionBegin;
151   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
152   if (iascii) {
153     ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %G\n",eis->omega);CHKERRQ(ierr);
154     if (eis->usediag) {
155       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr);
156     } else {
157       ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr);
158     }
159   } else {
160     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
161   }
162   PetscFunctionReturn(0);
163 }
164 
165 #undef __FUNCT__
166 #define __FUNCT__ "PCSetUp_Eisenstat"
167 static PetscErrorCode PCSetUp_Eisenstat(PC pc)
168 {
169   PetscErrorCode ierr;
170   PetscInt       M,N,m,n;
171   PC_Eisenstat   *eis = (PC_Eisenstat*)pc->data;
172 
173   PetscFunctionBegin;
174   if (!pc->setupcalled) {
175     ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr);
176     ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr);
177     ierr = MatCreate(((PetscObject)pc)->comm,&eis->shell);CHKERRQ(ierr);
178     ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr);
179     ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr);
180     ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr);
181     ierr = PetscLogObjectParent(pc,eis->shell);CHKERRQ(ierr);
182     ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);CHKERRQ(ierr);
183   }
184   if (!eis->usediag) PetscFunctionReturn(0);
185   if (!pc->setupcalled) {
186     ierr = MatGetVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr);
187     ierr = PetscLogObjectParent(pc,eis->diag);CHKERRQ(ierr);
188   }
189   ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr);
190   PetscFunctionReturn(0);
191 }
192 
193 /* --------------------------------------------------------------------*/
194 
195 EXTERN_C_BEGIN
196 #undef __FUNCT__
197 #define __FUNCT__ "PCEisenstatSetOmega_Eisenstat"
198 PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
199 {
200   PC_Eisenstat  *eis;
201 
202   PetscFunctionBegin;
203   if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
204   eis = (PC_Eisenstat*)pc->data;
205   eis->omega = omega;
206   PetscFunctionReturn(0);
207 }
208 EXTERN_C_END
209 
210 EXTERN_C_BEGIN
211 #undef __FUNCT__
212 #define __FUNCT__ "PCEisenstatNoDiagonalScaling_Eisenstat"
213 PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
214 {
215   PC_Eisenstat *eis;
216 
217   PetscFunctionBegin;
218   eis = (PC_Eisenstat*)pc->data;
219   eis->usediag = PETSC_FALSE;
220   PetscFunctionReturn(0);
221 }
222 EXTERN_C_END
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "PCEisenstatSetOmega"
226 /*@
227    PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
228    to use with Eisenstat's trick (where omega = 1.0 by default).
229 
230    Logically Collective on PC
231 
232    Input Parameters:
233 +  pc - the preconditioner context
234 -  omega - relaxation coefficient (0 < omega < 2)
235 
236    Options Database Key:
237 .  -pc_eisenstat_omega <omega> - Sets omega
238 
239    Notes:
240    The Eisenstat trick implementation of SSOR requires about 50% of the
241    usual amount of floating point operations used for SSOR + Krylov method;
242    however, the preconditioned problem must be solved with both left
243    and right preconditioning.
244 
245    To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
246    which can be chosen with the database options
247 $    -pc_type  sor  -pc_sor_symmetric
248 
249    Level: intermediate
250 
251 .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
252 
253 .seealso: PCSORSetOmega()
254 @*/
255 PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC pc,PetscReal omega)
256 {
257   PetscErrorCode ierr,(*f)(PC,PetscReal);
258 
259   PetscFunctionBegin;
260   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
261   PetscValidLogicalCollectiveReal(pc,omega,2);
262   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr);
263   if (f) {
264     ierr = (*f)(pc,omega);CHKERRQ(ierr);
265   }
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "PCEisenstatNoDiagonalScaling"
271 /*@
272    PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
273    not to do additional diagonal preconditioning. For matrices with a constant
274    along the diagonal, this may save a small amount of work.
275 
276    Logically Collective on PC
277 
278    Input Parameter:
279 .  pc - the preconditioner context
280 
281    Options Database Key:
282 .  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
283 
284    Level: intermediate
285 
286    Note:
287      If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
288    likley want to use this routine since it will save you some unneeded flops.
289 
290 .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
291 
292 .seealso: PCEisenstatSetOmega()
293 @*/
294 PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC pc)
295 {
296   PetscErrorCode ierr,(*f)(PC);
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
300   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);CHKERRQ(ierr);
301   if (f) {
302     ierr = (*f)(pc);CHKERRQ(ierr);
303   }
304   PetscFunctionReturn(0);
305 }
306 
307 /* --------------------------------------------------------------------*/
308 
309 /*MC
310      PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
311            preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
312 
313    Options Database Keys:
314 +  -pc_eisenstat_omega <omega> - Sets omega
315 -  -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
316 
317    Level: beginner
318 
319   Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
320 
321    Notes: Only implemented for the SeqAIJ matrix format.
322           Not a true parallel SOR, in parallel this implementation corresponds to block
323           Jacobi with SOR on each block.
324 
325 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
326            PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
327 M*/
328 
329 EXTERN_C_BEGIN
330 #undef __FUNCT__
331 #define __FUNCT__ "PCCreate_Eisenstat"
332 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Eisenstat(PC pc)
333 {
334   PetscErrorCode ierr;
335   PC_Eisenstat   *eis;
336 
337   PetscFunctionBegin;
338   ierr = PetscNewLog(pc,PC_Eisenstat,&eis);CHKERRQ(ierr);
339 
340   pc->ops->apply           = PCApply_Eisenstat;
341   pc->ops->presolve        = PCPreSolve_Eisenstat;
342   pc->ops->postsolve       = PCPostSolve_Eisenstat;
343   pc->ops->applyrichardson = 0;
344   pc->ops->setfromoptions  = PCSetFromOptions_Eisenstat;
345   pc->ops->destroy         = PCDestroy_Eisenstat;
346   pc->ops->view            = PCView_Eisenstat;
347   pc->ops->setup           = PCSetUp_Eisenstat;
348 
349   pc->data           = (void*)eis;
350   eis->omega         = 1.0;
351   eis->b             = 0;
352   eis->diag          = 0;
353   eis->usediag       = PETSC_TRUE;
354 
355   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
356                     PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr);
357   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
358                     "PCEisenstatNoDiagonalScaling_Eisenstat",
359                     PCEisenstatNoDiagonalScaling_Eisenstat);CHKERRQ(ierr);
360  PetscFunctionReturn(0);
361 }
362 EXTERN_C_END
363