Actual source code: cayley.c
1: /*
2: Implements the Cayley spectral transform.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <private/stimpl.h> /*I "slepcst.h" I*/
26: typedef struct {
27: PetscScalar nu;
28: PetscBool nu_set;
29: Vec w2;
30: } ST_CAYLEY;
34: PetscErrorCode STApply_Cayley(ST st,Vec x,Vec y)
35: {
37: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
38: PetscScalar nu = ctx->nu;
39:
41: if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
43: if (st->B) {
44: /* generalized eigenproblem: y = (A - sB)^-1 (A + tB)x */
45: MatMult(st->A,x,st->w);
46: MatMult(st->B,x,ctx->w2);
47: VecAXPY(st->w,nu,ctx->w2);
48: STAssociatedKSPSolve(st,st->w,y);
49: }
50: else {
51: /* standard eigenproblem: y = (A - sI)^-1 (A + tI)x */
52: MatMult(st->A,x,st->w);
53: VecAXPY(st->w,nu,x);
54: STAssociatedKSPSolve(st,st->w,y);
55: }
56: return(0);
57: }
61: PetscErrorCode STApplyTranspose_Cayley(ST st,Vec x,Vec y)
62: {
64: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
65: PetscScalar nu = ctx->nu;
66:
68: if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
70: if (st->B) {
71: /* generalized eigenproblem: y = (A + tB)^T (A - sB)^-T x */
72: STAssociatedKSPSolveTranspose(st,x,st->w);
73: MatMultTranspose(st->A,st->w,y);
74: MatMultTranspose(st->B,st->w,ctx->w2);
75: VecAXPY(y,nu,ctx->w2);
76: }
77: else {
78: /* standard eigenproblem: y = (A + tI)^T (A - sI)^-T x */
79: STAssociatedKSPSolveTranspose(st,x,st->w);
80: MatMultTranspose(st->A,st->w,y);
81: VecAXPY(y,nu,st->w);
82: }
83: return(0);
84: }
88: PetscErrorCode STBilinearMatMult_Cayley(Mat B,Vec x,Vec y)
89: {
91: ST st;
92: ST_CAYLEY *ctx;
93: PetscScalar nu;
94:
96: MatShellGetContext(B,(void**)&st);
97: ctx = (ST_CAYLEY*)st->data;
98: nu = ctx->nu;
99:
100: if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
102: if (st->B) {
103: /* generalized eigenproblem: y = (A + tB)x */
104: MatMult(st->A,x,y);
105: MatMult(st->B,x,ctx->w2);
106: VecAXPY(y,nu,ctx->w2);
107: }
108: else {
109: /* standard eigenproblem: y = (A + tI)x */
110: MatMult(st->A,x,y);
111: VecAXPY(y,nu,x);
112: }
113: return(0);
114: }
118: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
119: {
121: PetscInt n,m;
124: MatGetLocalSize(st->B,&n,&m);
125: MatCreateShell(((PetscObject)st)->comm,n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,B);
126: MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))STBilinearMatMult_Cayley);
127: return(0);
128: }
132: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
133: {
134: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
135: PetscInt j;
136: #if !defined(PETSC_USE_COMPLEX)
137: PetscScalar t,i,r;
138: #endif
141: #if !defined(PETSC_USE_COMPLEX)
142: for (j=0;j<n;j++) {
143: if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
144: else {
145: r = eigr[j];
146: i = eigi[j];
147: r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
148: i = - st->sigma * i - ctx->nu * i;
149: t = i * i + r * (r - 2.0) + 1.0;
150: eigr[j] = r / t;
151: eigi[j] = i / t;
152: }
153: }
154: #else
155: for (j=0;j<n;j++) {
156: eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
157: }
158: #endif
159: return(0);
160: }
164: PetscErrorCode STPostSolve_Cayley(ST st)
165: {
169: if (st->shift_matrix == ST_MATMODE_INPLACE) {
170: if (st->B) {
171: MatAXPY(st->A,st->sigma,st->B,st->str);
172: } else {
173: MatShift(st->A,st->sigma);
174: }
175: st->setupcalled = 0;
176: }
177: return(0);
178: }
182: PetscErrorCode STSetUp_Cayley(ST st)
183: {
185: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
188: MatDestroy(&st->mat);
190: /* if the user did not set the shift, use the target value */
191: if (!st->sigma_set) st->sigma = st->defsigma;
193: if (!ctx->nu_set) { ctx->nu = st->sigma; }
194: if (ctx->nu == 0.0 && st->sigma == 0.0) {
195: SETERRQ(((PetscObject)st)->comm,1,"Values of shift and antishift cannot be zero simultaneously");
196: }
198: if (!st->ksp) { STGetKSP(st,&st->ksp); }
199: switch (st->shift_matrix) {
200: case ST_MATMODE_INPLACE:
201: st->mat = PETSC_NULL;
202: if (st->sigma != 0.0) {
203: if (st->B) {
204: MatAXPY(st->A,-st->sigma,st->B,st->str);
205: } else {
206: MatShift(st->A,-st->sigma);
207: }
208: }
209: KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
210: break;
211: case ST_MATMODE_SHELL:
212: STMatShellCreate(st,&st->mat);
213: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
214: break;
215: default:
216: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
217: if (st->sigma != 0.0) {
218: if (st->B) {
219: MatAXPY(st->mat,-st->sigma,st->B,st->str);
220: } else {
221: MatShift(st->mat,-st->sigma);
222: }
223: }
224: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
225: }
226: if (st->B) {
227: VecDestroy(&ctx->w2);
228: MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
229: }
230: KSPSetUp(st->ksp);
231: return(0);
232: }
236: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
237: {
239: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
240: MatStructure flg;
243: if (!ctx->nu_set) { ctx->nu = newshift; }
244: if (ctx->nu == 0.0 && newshift == 0.0) {
245: SETERRQ(((PetscObject)st)->comm,1,"Values of shift and antishift cannot be zero simultaneously");
246: }
248: /* Nothing to be done if STSetUp has not been called yet */
249: if (!st->setupcalled) return(0);
251: /* Check if the new KSP matrix has the same zero structure */
252: if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
253: flg = DIFFERENT_NONZERO_PATTERN;
254: } else {
255: flg = SAME_NONZERO_PATTERN;
256: }
258: switch (st->shift_matrix) {
259: case ST_MATMODE_INPLACE:
260: /* Undo previous operations */
261: if (st->sigma != 0.0) {
262: if (st->B) {
263: MatAXPY(st->A,st->sigma,st->B,st->str);
264: } else {
265: MatShift(st->A,st->sigma);
266: }
267: }
268: /* Apply new shift */
269: if (newshift != 0.0) {
270: if (st->B) {
271: MatAXPY(st->A,-newshift,st->B,st->str);
272: } else {
273: MatShift(st->A,-newshift);
274: }
275: }
276: KSPSetOperators(st->ksp,st->A,st->A,flg);
277: break;
278: case ST_MATMODE_SHELL:
279: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
280: break;
281: default:
282: MatCopy(st->A,st->mat,DIFFERENT_NONZERO_PATTERN);
283: if (newshift != 0.0) {
284: if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
285: else { MatShift(st->mat,-newshift); }
286: }
287: KSPSetOperators(st->ksp,st->mat,st->mat,flg);
288: }
289: st->sigma = newshift;
290: KSPSetUp(st->ksp);
291: return(0);
292: }
296: PetscErrorCode STSetFromOptions_Cayley(ST st)
297: {
299: PetscScalar nu;
300: PetscBool flg;
301: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
302: PC pc;
303: const PCType pctype;
304: const KSPType ksptype;
307: if (!st->ksp) { STGetKSP(st,&st->ksp); }
308: KSPGetPC(st->ksp,&pc);
309: KSPGetType(st->ksp,&ksptype);
310: PCGetType(pc,&pctype);
311: if (!pctype && !ksptype) {
312: if (st->shift_matrix == ST_MATMODE_SHELL) {
313: /* in shell mode use GMRES with Jacobi as the default */
314: KSPSetType(st->ksp,KSPGMRES);
315: PCSetType(pc,PCJACOBI);
316: } else {
317: /* use direct solver as default */
318: KSPSetType(st->ksp,KSPPREONLY);
319: PCSetType(pc,PCREDUNDANT);
320: }
321: }
323: PetscOptionsHead("ST Cayley Options");
324: PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
325: if (flg) {
326: STCayleySetAntishift(st,nu);
327: }
328: PetscOptionsTail();
329: return(0);
330: }
332: EXTERN_C_BEGIN
335: PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
336: {
337: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
340: ctx->nu = newshift;
341: ctx->nu_set = PETSC_TRUE;
342: return(0);
343: }
344: EXTERN_C_END
348: /*@
349: STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
350: spectral transformation.
352: Logically Collective on ST
354: Input Parameters:
355: + st - the spectral transformation context
356: - nu - the anti-shift
358: Options Database Key:
359: . -st_cayley_antishift - Sets the value of the anti-shift
361: Level: intermediate
363: Note:
364: In the generalized Cayley transform, the operator can be expressed as
365: OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
366: Use STSetShift() for setting sigma.
368: .seealso: STSetShift(), STCayleyGetAntishift()
369: @*/
370: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
371: {
377: PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
378: return(0);
379: }
380: EXTERN_C_BEGIN
383: PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
384: {
385: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
388: *nu = ctx->nu;
389: return(0);
390: }
391: EXTERN_C_END
395: /*@
396: STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
397: spectral transformation.
399: Not Collective
401: Input Parameter:
402: . st - the spectral transformation context
404: Output Parameter:
405: . nu - the anti-shift
407: Level: intermediate
409: .seealso: STGetShift(), STCayleySetAntishift()
410: @*/
411: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
412: {
418: PetscTryMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
419: return(0);
420: }
424: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
425: {
427: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
430: #if !defined(PETSC_USE_COMPLEX)
431: PetscViewerASCIIPrintf(viewer," Cayley: antishift: %G\n",ctx->nu);
432: #else
433: PetscViewerASCIIPrintf(viewer," Cayley: antishift: %G+%G i\n",PetscRealPart(ctx->nu),PetscImaginaryPart(ctx->nu));
434: #endif
435: return(0);
436: }
440: PetscErrorCode STReset_Cayley(ST st)
441: {
443: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
446: VecDestroy(&ctx->w2);
447: return(0);
448: }
452: PetscErrorCode STDestroy_Cayley(ST st)
453: {
457: PetscFree(st->data);
458: PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","",PETSC_NULL);
459: return(0);
460: }
462: EXTERN_C_BEGIN
465: PetscErrorCode STCreate_Cayley(ST st)
466: {
470: PetscNewLog(st,ST_CAYLEY,&st->data);
471: st->ops->apply = STApply_Cayley;
472: st->ops->getbilinearform = STGetBilinearForm_Cayley;
473: st->ops->applytrans = STApplyTranspose_Cayley;
474: st->ops->postsolve = STPostSolve_Cayley;
475: st->ops->backtr = STBackTransform_Cayley;
476: st->ops->setfromoptions = STSetFromOptions_Cayley;
477: st->ops->setup = STSetUp_Cayley;
478: st->ops->setshift = STSetShift_Cayley;
479: st->ops->destroy = STDestroy_Cayley;
480: st->ops->reset = STReset_Cayley;
481: st->ops->view = STView_Cayley;
482: st->ops->checknullspace = STCheckNullSpace_Default;
483: PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","STCayleySetAntishift_Cayley",STCayleySetAntishift_Cayley);
484: return(0);
485: }
486: EXTERN_C_END