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-2009, Universidad 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
26: typedef struct {
27: PetscScalar tau;
28: PetscTruth tau_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 tau = ctx->tau;
39:
41: if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + 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,tau,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,tau,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 tau = ctx->tau;
66:
68: if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + 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,tau,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,tau,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 tau;
94:
96: MatShellGetContext(B,(void**)&st);
97: ctx = (ST_CAYLEY *) st->data;
98: tau = ctx->tau;
99:
100: if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + 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,tau,ctx->w2);
107: }
108: else {
109: /* standard eigenproblem: y = (A + tI)x */
110: MatMult(st->A,x,y);
111: VecAXPY(y,tau,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,PetscScalar *eigr,PetscScalar *eigi)
133: {
134: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
135: #ifndef PETSC_USE_COMPLEX
136: PetscScalar t,i,r;
140: if (*eigi == 0.0) *eigr = (ctx->tau + *eigr * st->sigma) / (*eigr - 1.0);
141: else {
142: r = *eigr;
143: i = *eigi;
144: r = st->sigma * (r * r + i * i - r) + ctx->tau * (r - 1);
145: i = - st->sigma * i - ctx->tau * i;
146: t = i * i + r * (r - 2.0) + 1.0;
147: *eigr = r / t;
148: *eigi = i / t;
149: }
150: #else
153: *eigr = (ctx->tau + *eigr * st->sigma) / (*eigr - 1.0);
154: #endif
155: return(0);
156: }
160: PetscErrorCode STPostSolve_Cayley(ST st)
161: {
165: if (st->shift_matrix == STMATMODE_INPLACE) {
166: if (st->B) {
167: MatAXPY(st->A,st->sigma,st->B,st->str);
168: } else {
169: MatShift(st->A,st->sigma);
170: }
171: st->setupcalled = 0;
172: }
173: return(0);
174: }
178: PetscErrorCode STSetUp_Cayley(ST st)
179: {
181: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
185: if (st->mat) { MatDestroy(st->mat); }
186: if (!ctx->tau_set) { ctx->tau = st->sigma; }
187: if (ctx->tau == 0.0 && st->sigma == 0.0) {
188: SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
189: }
191: switch (st->shift_matrix) {
192: case STMATMODE_INPLACE:
193: st->mat = PETSC_NULL;
194: if (st->sigma != 0.0) {
195: if (st->B) {
196: MatAXPY(st->A,-st->sigma,st->B,st->str);
197: } else {
198: MatShift(st->A,-st->sigma);
199: }
200: }
201: KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
202: break;
203: case STMATMODE_SHELL:
204: STMatShellCreate(st,&st->mat);
205: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
206: break;
207: default:
208: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
209: if (st->sigma != 0.0) {
210: if (st->B) {
211: MatAXPY(st->mat,-st->sigma,st->B,st->str);
212: } else {
213: MatShift(st->mat,-st->sigma);
214: }
215: }
216: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
217: }
218: if (st->B) {
219: if (ctx->w2) { VecDestroy(ctx->w2); }
220: MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
221: }
222: KSPSetUp(st->ksp);
223: return(0);
224: }
228: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
229: {
231: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
232: MatStructure flg;
235: if (!ctx->tau_set) { ctx->tau = newshift; }
236: if (ctx->tau == 0.0 && newshift == 0.0) {
237: SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
238: }
240: /* Nothing to be done if STSetUp has not been called yet */
241: if (!st->setupcalled) return(0);
243: /* Check if the new KSP matrix has the same zero structure */
244: if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
245: flg = DIFFERENT_NONZERO_PATTERN;
246: } else {
247: flg = SAME_NONZERO_PATTERN;
248: }
250: switch (st->shift_matrix) {
251: case STMATMODE_INPLACE:
252: /* Undo previous operations */
253: if (st->sigma != 0.0) {
254: if (st->B) {
255: MatAXPY(st->A,st->sigma,st->B,st->str);
256: } else {
257: MatShift(st->A,st->sigma);
258: }
259: }
260: /* Apply new shift */
261: if (newshift != 0.0) {
262: if (st->B) {
263: MatAXPY(st->A,-newshift,st->B,st->str);
264: } else {
265: MatShift(st->A,-newshift);
266: }
267: }
268: KSPSetOperators(st->ksp,st->A,st->A,flg);
269: break;
270: case STMATMODE_SHELL:
271: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
272: break;
273: default:
274: MatCopy(st->A, st->mat,SUBSET_NONZERO_PATTERN);
275: if (newshift != 0.0) {
276: if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
277: else { MatShift(st->mat,-newshift); }
278: }
279: KSPSetOperators(st->ksp,st->mat,st->mat,flg);
280: }
281: st->sigma = newshift;
282: KSPSetUp(st->ksp);
283: return(0);
284: }
288: PetscErrorCode STSetFromOptions_Cayley(ST st)
289: {
291: PetscScalar tau;
292: PetscTruth flg;
293: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
296: PetscOptionsHead("ST Cayley Options");
297: PetscOptionsScalar("-st_antishift","Value of the antishift","STSetAntishift",ctx->tau,&tau,&flg);
298: if (flg) {
299: STCayleySetAntishift(st,tau);
300: }
301: PetscOptionsTail();
302: return(0);
303: }
308: PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
309: {
310: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
313: ctx->tau = newshift;
314: ctx->tau_set = PETSC_TRUE;
315: return(0);
316: }
321: /*@
322: STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
323: spectral transformation.
325: Collective on ST
327: Input Parameters:
328: + st - the spectral transformation context
329: - tau - the anti-shift
331: Options Database Key:
332: . -st_antishift - Sets the value of the anti-shift
334: Level: intermediate
336: Note:
337: In the generalized Cayley transform, the operator can be expressed as
338: OP = inv(A - sigma B)*(A + tau B). This function sets the value of tau.
339: Use STSetShift() for setting sigma.
341: .seealso: STSetShift()
342: @*/
343: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar tau)
344: {
345: PetscErrorCode ierr, (*f)(ST,PetscScalar);
349: PetscObjectQueryFunction((PetscObject)st,"STCayleySetAntishift_C",(void (**)(void))&f);
350: if (f) {
351: (*f)(st,tau);
352: }
353: return(0);
354: }
358: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
359: {
361: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
364: #if !defined(PETSC_USE_COMPLEX)
365: PetscViewerASCIIPrintf(viewer," antishift: %g\n",ctx->tau);
366: #else
367: PetscViewerASCIIPrintf(viewer," antishift: %g+%g i\n",PetscRealPart(ctx->tau),PetscImaginaryPart(ctx->tau));
368: #endif
369: STView_Default(st,viewer);
370: return(0);
371: }
375: PetscErrorCode STDestroy_Cayley(ST st)
376: {
378: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
381: if (ctx->w2) { VecDestroy(ctx->w2); }
382: PetscFree(ctx);
383: return(0);
384: }
389: PetscErrorCode STCreate_Cayley(ST st)
390: {
392: ST_CAYLEY *ctx;
395: PetscNew(ST_CAYLEY,&ctx);
396: PetscLogObjectMemory(st,sizeof(ST_CAYLEY));
397: st->data = (void *) ctx;
399: st->ops->apply = STApply_Cayley;
400: st->ops->getbilinearform = STGetBilinearForm_Cayley;
401: st->ops->applytrans = STApplyTranspose_Cayley;
402: st->ops->postsolve = STPostSolve_Cayley;
403: st->ops->backtr = STBackTransform_Cayley;
404: st->ops->setfromoptions = STSetFromOptions_Cayley;
405: st->ops->setup = STSetUp_Cayley;
406: st->ops->setshift = STSetShift_Cayley;
407: st->ops->destroy = STDestroy_Cayley;
408: st->ops->view = STView_Cayley;
409:
410: st->checknullspace = STCheckNullSpace_Default;
412: ctx->tau = 0.0;
413: ctx->tau_set = PETSC_FALSE;
415: PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","STCayleySetAntishift_Cayley",
416: STCayleySetAntishift_Cayley);
418: return(0);
419: }