Actual source code: eisen.c
petsc-3.5.4 2015-05-23
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*/
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;
20: static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x)
21: {
23: PC pc;
24: PC_Eisenstat *eis;
27: MatShellGetContext(mat,(void**)&pc);
28: eis = (PC_Eisenstat*)pc->data;
29: MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);
30: return(0);
31: }
35: static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y)
36: {
37: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
39: PetscBool hasop;
42: if (eis->usediag) {
43: MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
44: if (hasop) {
45: MatMultDiagonalBlock(pc->pmat,x,y);
46: } else {
47: VecPointwiseMult(y,x,eis->diag);
48: }
49: } else {VecCopy(x,y);}
50: return(0);
51: }
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;
62: if (pc->presolvedone < 2) {
63: if (pc->mat != pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),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: }
69: if (!eis->b[pc->presolvedone-1]) {
70: VecDuplicate(b,&eis->b[pc->presolvedone-1]);
71: PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->b[pc->presolvedone-1]);
72: }
74: /* if nonzero initial guess, modify x */
75: KSPGetInitialGuessNonzero(ksp,&nonzero);
76: if (nonzero) {
77: VecCopy(x,eis->b[pc->presolvedone-1]);
78: MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);
79: }
81: /* save true b, other option is to swap pointers */
82: VecCopy(b,eis->b[pc->presolvedone-1]);
84: /* modify b by (L + D/omega)^{-1} */
85: MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);
86: return(0);
87: }
91: static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x)
92: {
93: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
97: /* get back true b */
98: VecCopy(eis->b[pc->presolvedone],b);
100: /* modify x by (U + D/omega)^{-1} */
101: VecCopy(x,eis->b[pc->presolvedone]);
102: MatSOR(eis->A,eis->b[pc->presolvedone],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);
103: if (!pc->presolvedone) pc->mat = eis->A;
104: return(0);
105: }
109: static PetscErrorCode PCReset_Eisenstat(PC pc)
110: {
111: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
115: VecDestroy(&eis->b[0]);
116: VecDestroy(&eis->b[1]);
117: MatDestroy(&eis->shell);
118: VecDestroy(&eis->diag);
119: return(0);
120: }
124: static PetscErrorCode PCDestroy_Eisenstat(PC pc)
125: {
129: PCReset_Eisenstat(pc);
130: PetscFree(pc->data);
131: return(0);
132: }
136: static PetscErrorCode PCSetFromOptions_Eisenstat(PC pc)
137: {
138: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
140: PetscBool flg = PETSC_FALSE;
143: PetscOptionsHead("Eisenstat SSOR options");
144: PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,NULL);
145: PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatNoDiagonalScaling",flg,&flg,NULL);
146: if (flg) {
147: PCEisenstatNoDiagonalScaling(pc);
148: }
149: PetscOptionsTail();
150: return(0);
151: }
155: static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer)
156: {
157: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
159: PetscBool iascii;
162: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
163: if (iascii) {
164: PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",(double)eis->omega);
165: if (eis->usediag) {
166: PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");
167: } else {
168: PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");
169: }
170: }
171: return(0);
172: }
176: static PetscErrorCode PCSetUp_Eisenstat(PC pc)
177: {
179: PetscInt M,N,m,n;
180: PC_Eisenstat *eis = (PC_Eisenstat*)pc->data;
183: if (!pc->setupcalled) {
184: MatGetSize(pc->mat,&M,&N);
185: MatGetLocalSize(pc->mat,&m,&n);
186: MatCreate(PetscObjectComm((PetscObject)pc),&eis->shell);
187: MatSetSizes(eis->shell,m,n,M,N);
188: MatSetType(eis->shell,MATSHELL);
189: MatSetUp(eis->shell);
190: MatShellSetContext(eis->shell,(void*)pc);
191: PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->shell);
192: MatShellSetOperation(eis->shell,MATOP_MULT,(void (*)(void))PCMult_Eisenstat);
193: }
194: if (!eis->usediag) return(0);
195: if (!pc->setupcalled) {
196: MatGetVecs(pc->pmat,&eis->diag,0);
197: PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->diag);
198: }
199: MatGetDiagonal(pc->pmat,eis->diag);
200: return(0);
201: }
203: /* --------------------------------------------------------------------*/
207: static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega)
208: {
209: PC_Eisenstat *eis;
212: if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range");
213: eis = (PC_Eisenstat*)pc->data;
214: eis->omega = omega;
215: return(0);
216: }
220: static PetscErrorCode PCEisenstatNoDiagonalScaling_Eisenstat(PC pc)
221: {
222: PC_Eisenstat *eis;
225: eis = (PC_Eisenstat*)pc->data;
226: eis->usediag = PETSC_FALSE;
227: return(0);
228: }
232: /*@
233: PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega,
234: to use with Eisenstat's trick (where omega = 1.0 by default).
236: Logically Collective on PC
238: Input Parameters:
239: + pc - the preconditioner context
240: - omega - relaxation coefficient (0 < omega < 2)
242: Options Database Key:
243: . -pc_eisenstat_omega <omega> - Sets omega
245: Notes:
246: The Eisenstat trick implementation of SSOR requires about 50% of the
247: usual amount of floating point operations used for SSOR + Krylov method;
248: however, the preconditioned problem must be solved with both left
249: and right preconditioning.
251: To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner,
252: which can be chosen with the database options
253: $ -pc_type sor -pc_sor_symmetric
255: Level: intermediate
257: .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega
259: .seealso: PCSORSetOmega()
260: @*/
261: PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega)
262: {
268: PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));
269: return(0);
270: }
274: /*@
275: PCEisenstatNoDiagonalScaling - Causes the Eisenstat preconditioner
276: not to do additional diagonal preconditioning. For matrices with a constant
277: along the diagonal, this may save a small amount of work.
279: Logically Collective on PC
281: Input Parameter:
282: . pc - the preconditioner context
284: Options Database Key:
285: . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
287: Level: intermediate
289: Note:
290: If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will
291: likley want to use this routine since it will save you some unneeded flops.
293: .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR
295: .seealso: PCEisenstatSetOmega()
296: @*/
297: PetscErrorCode PCEisenstatNoDiagonalScaling(PC pc)
298: {
303: PetscTryMethod(pc,"PCEisenstatNoDiagonalScaling_C",(PC),(pc));
304: return(0);
305: }
307: /* --------------------------------------------------------------------*/
309: /*MC
310: PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel)
311: preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed.
313: Options Database Keys:
314: + -pc_eisenstat_omega <omega> - Sets omega
315: - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatNoDiagonalScaling()
317: Level: beginner
319: Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick
321: Notes: Only implemented for the SeqAIJ matrix format.
322: Not a true parallel SOR, in parallel this implementation corresponds to block
323: Jacobi with SOR on each block.
325: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
326: PCEisenstatNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR
327: M*/
331: PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc)
332: {
334: PC_Eisenstat *eis;
337: PetscNewLog(pc,&eis);
339: pc->ops->apply = PCApply_Eisenstat;
340: pc->ops->presolve = PCPreSolve_Eisenstat;
341: pc->ops->postsolve = PCPostSolve_Eisenstat;
342: pc->ops->applyrichardson = 0;
343: pc->ops->setfromoptions = PCSetFromOptions_Eisenstat;
344: pc->ops->destroy = PCDestroy_Eisenstat;
345: pc->ops->reset = PCReset_Eisenstat;
346: pc->ops->view = PCView_Eisenstat;
347: pc->ops->setup = PCSetUp_Eisenstat;
349: pc->data = (void*)eis;
350: eis->omega = 1.0;
351: eis->b[0] = 0;
352: eis->b[1] = 0;
353: eis->diag = 0;
354: eis->usediag = PETSC_TRUE;
356: PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat);
357: PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatNoDiagonalScaling_C",PCEisenstatNoDiagonalScaling_Eisenstat);
358: return(0);
359: }