Actual source code: eisen.c
1: /*
2: Defines a Eisenstat trick SSOR preconditioner. This uses about
3: %50 of the usual amount of floating point ops used for SSOR + Krylov
4: method. But it requires actually solving the preconditioned problem
5: with both left and right preconditioning.
6: */
7: #include src/ksp/pc/pcimpl.h
9: typedef struct {
10: Mat shell,A;
11: Vec b,diag; /* temporary storage for true right hand side */
12: PetscReal omega;
13: PetscTruth usediag; /* indicates preconditioner should include diagonal scaling*/
14: } PC_Eisenstat;
19: static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
20: {
22: PC pc;
23: PC_Eisenstat *eis;
26: MatShellGetContext(mat,(void **)&pc);
27: eis = (PC_Eisenstat*)pc->data;
28: MatRelax(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);
29: return(0);
30: }
34: static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
35: {
36: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
40: if (eis->usediag) {VecPointwiseMult(x,eis->diag,y);}
41: else {VecCopy(x,y);}
42: return(0);
43: }
47: static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec x,Vec b)
48: {
49: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
50: PetscTruth nonzero;
54: if (pc->mat != pc->pmat) SETERRQ(PETSC_ERR_SUP,"Cannot have different mat and pmat");
55:
56: /* swap shell matrix and true matrix */
57: eis->A = pc->mat;
58: pc->mat = eis->shell;
60: if (!eis->b) {
61: VecDuplicate(b,&eis->b);
62: PetscLogObjectParent(pc,eis->b);
63: }
64:
65: /* save true b, other option is to swap pointers */
66: VecCopy(b,eis->b);
68: /* if nonzero initial guess, modify x */
69: KSPGetInitialGuessNonzero(ksp,&nonzero);
70: if (nonzero) {
71: MatRelax(eis->A,x,eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);
72: }
74: /* modify b by (L + D)^{-1} */
75: MatRelax(eis->A,b,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_FORWARD_SWEEP),0.0,1,1,b);
76: return(0);
77: }
81: static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec x,Vec b)
82: {
83: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
87: MatRelax(eis->A,x,eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_BACKWARD_SWEEP),0.0,1,1,x);
88: pc->mat = eis->A;
89: /* get back true b */
90: VecCopy(eis->b,b);
91: return(0);
92: }
96: static PetscErrorCode PCDestroy_Eisenstat(PC pc)
97: {
98: PC_Eisenstat *eis = (PC_Eisenstat *)pc->data;
102: if (eis->b) {VecDestroy(eis->b);}
103: if (eis->shell) {MatDestroy(eis->shell);}
104: if (eis->diag) {VecDestroy(eis->diag);}
105: PetscFree(eis);
106: return(0);
107: }
111: static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
112: {
113: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
115: PetscTruth flg;
118: PetscOptionsHead("Eisenstat SSOR options");
119: PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,0);
120: PetscOptionsName("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",&flg);
121: if (flg) {
122: PCEisenstatNoDiagonalScaling(pc);
123: }
124: PetscOptionsTail();
125: return(0);
126: }
130: static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
131: {
132: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
134: PetscTruth iascii;
137: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
138: if (iascii) {
139: PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",eis->omega);
140: if (eis->usediag) {
141: PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");
142: } else {
143: PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");
144: }
145: } else {
146: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for Eisenstat PC",((PetscObject)viewer)->type_name);
147: }
148: return(0);
149: }
153: static PetscErrorCode PCSetUp_Eisenstat(PC pc)
154: {
156: PetscInt M,N,m,n;
157: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
160: if (!pc->setupcalled) {
161: MatGetSize(pc->mat,&M,&N);
162: MatGetLocalSize(pc->mat,&m,&n);
163: MatCreate(pc->comm,m,N,M,N,&eis->shell);
164: MatSetType(eis->shell,MATSHELL);
165: MatShellSetContext(eis->shell,(void*)pc);
166: PetscLogObjectParent(pc,eis->shell);
167: MatShellSetOperation(eis->shell,MATOP_MULT,(void(*)(void))PCMult_Eisenstat);
168: }
169: if (!eis->usediag) return(0);
170: if (!pc->setupcalled) {
171: MatGetVecs(pc->pmat,&eis->diag,0);
172: PetscLogObjectParent(pc,eis->diag);
173: }
174: MatGetDiagonal(pc->pmat,eis->diag);
175: return(0);
176: }
178: /* --------------------------------------------------------------------*/
183: PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
184: {
185: PC_Eisenstat *eis;
188: if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
189: eis = (PC_Eisenstat*)pc->data;
190: eis->omega = omega;
191: return(0);
192: }
198: PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
199: {
200: PC_Eisenstat *eis;
203: eis = (PC_Eisenstat*)pc->data;
204: eis->usediag = PETSC_FALSE;
205: return(0);
206: }
211: /*@
212: PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
213: to use with Eisenstat's trick (where omega = 1.0 by default).
215: Collective on PC
217: Input Parameters:
218: + pc - the preconditioner context
219: - omega - relaxation coefficient (0 < omega < 2)
221: Options Database Key:
222: . -pc_eisenstat_omega <omega> - Sets omega
224: Notes:
225: The Eisenstat trick implementation of SSOR requires about 50% of the
226: usual amount of floating point operations used for SSOR + Krylov method;
227: however, the preconditioned problem must be solved with both left
228: and right preconditioning.
230: To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
231: which can be chosen with the database options
232: $ -pc_type sor -pc_sor_symmetric
234: Level: intermediate
236: .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
238: .seealso: PCSORSetOmega()
239: @*/
240: PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega)
241: {
242: PetscErrorCode ierr,(*f)(PC,PetscReal);
246: PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatSetOmega_C",(void (**)(void))&f);
247: if (f) {
248: (*f)(pc,omega);
249: }
250: return(0);
251: }
255: /*@
256: PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
257: not to do additional diagonal preconditioning. For matrices with a constant
258: along the diagonal, this may save a small amount of work.
260: Collective on PC
262: Input Parameter:
263: . pc - the preconditioner context
265: Options Database Key:
266: . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
268: Level: intermediate
270: Note:
271: If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
272: likley want to use this routine since it will save you some unneeded flops.
274: .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
276: .seealso: PCEisenstatSetOmega()
277: @*/
278: PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc)
279: {
280: PetscErrorCode ierr,(*f)(PC);
284: PetscObjectQueryFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",(void (**)(void))&f);
285: if (f) {
286: (*f)(pc);
287: }
288: return(0);
289: }
291: /* --------------------------------------------------------------------*/
293: /*MC
294: PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
295: preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
297: Options Database Keys:
298: + -pc_eisenstat_omega <omega> - Sets omega
299: - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
301: Level: beginner
303: Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
305: Notes: Only implemented for the SeqAIJ matrix format.
306: Not a true parallel SOR, in parallel this implementation corresponds to block
307: Jacobi with SOR on each block.
309: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
310: PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
311: M*/
316: PetscErrorCode PCCreate_Eisenstat(PC pc)
317: {
319: PC_Eisenstat *eis;
322: PetscNew(PC_Eisenstat,&eis);
323: PetscLogObjectMemory(pc,sizeof(PC_Eisenstat));
325: pc->ops->apply = PCApply_Eisenstat;
326: pc->ops->presolve = PCPreSolve_Eisenstat;
327: pc->ops->postsolve = PCPostSolve_Eisenstat;
328: pc->ops->applyrichardson = 0;
329: pc->ops->setfromoptions = PCSetFromOptions_Eisenstat;
330: pc->ops->destroy = PCDestroy_Eisenstat;
331: pc->ops->view = PCView_Eisenstat;
332: pc->ops->setup = PCSetUp_Eisenstat;
334: pc->data = (void*)eis;
335: eis->omega = 1.0;
336: eis->b = 0;
337: eis->diag = 0;
338: eis->usediag = PETSC_TRUE;
340: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatSetOmega_C","PCEisenstatSetOmega_Eisenstat",
341: PCEisenstatSetOmega_Eisenstat);
342: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",
343: "PCEisenstatNoDiagonalScaling_Eisenstat",
344: PCEisenstatNoDiagonalScaling_Eisenstat);
345: return(0);
346: }