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