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