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