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