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-2010, 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 nu;
28: PetscTruth 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: #ifndef PETSC_USE_COMPLEX
137: PetscScalar t,i,r;
141: for (j=0;j<n;j++) {
142: if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
143: else {
144: r = eigr[j];
145: i = eigi[j];
146: r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
147: i = - st->sigma * i - ctx->nu * i;
148: t = i * i + r * (r - 2.0) + 1.0;
149: eigr[j] = r / t;
150: eigi[j] = i / t;
151: }
152: }
153: #else
156: for (j=0;j<n;j++) {
157: eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
158: }
159: #endif
160: return(0);
161: }
165: PetscErrorCode STPostSolve_Cayley(ST st)
166: {
170: if (st->shift_matrix == ST_MATMODE_INPLACE) {
171: if (st->B) {
172: MatAXPY(st->A,st->sigma,st->B,st->str);
173: } else {
174: MatShift(st->A,st->sigma);
175: }
176: st->setupcalled = 0;
177: }
178: return(0);
179: }
183: PetscErrorCode STSetUp_Cayley(ST st)
184: {
186: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
190: if (st->mat) { MatDestroy(st->mat); }
192: /* if the user did not set the shift, use the target value */
193: if (!st->sigma_set) st->sigma = st->defsigma;
195: if (!ctx->nu_set) { ctx->nu = st->sigma; }
196: if (ctx->nu == 0.0 && st->sigma == 0.0) {
197: SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
198: }
200: switch (st->shift_matrix) {
201: case ST_MATMODE_INPLACE:
202: st->mat = PETSC_NULL;
203: if (st->sigma != 0.0) {
204: if (st->B) {
205: MatAXPY(st->A,-st->sigma,st->B,st->str);
206: } else {
207: MatShift(st->A,-st->sigma);
208: }
209: }
210: KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
211: break;
212: case ST_MATMODE_SHELL:
213: STMatShellCreate(st,&st->mat);
214: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
215: break;
216: default:
217: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
218: if (st->sigma != 0.0) {
219: if (st->B) {
220: MatAXPY(st->mat,-st->sigma,st->B,st->str);
221: } else {
222: MatShift(st->mat,-st->sigma);
223: }
224: }
225: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
226: }
227: if (st->B) {
228: if (ctx->w2) { VecDestroy(ctx->w2); }
229: MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
230: }
231: KSPSetUp(st->ksp);
232: return(0);
233: }
237: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
238: {
240: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
241: MatStructure flg;
244: if (!ctx->nu_set) { ctx->nu = newshift; }
245: if (ctx->nu == 0.0 && newshift == 0.0) {
246: SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
247: }
249: /* Nothing to be done if STSetUp has not been called yet */
250: if (!st->setupcalled) return(0);
252: /* Check if the new KSP matrix has the same zero structure */
253: if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
254: flg = DIFFERENT_NONZERO_PATTERN;
255: } else {
256: flg = SAME_NONZERO_PATTERN;
257: }
259: switch (st->shift_matrix) {
260: case ST_MATMODE_INPLACE:
261: /* Undo previous operations */
262: if (st->sigma != 0.0) {
263: if (st->B) {
264: MatAXPY(st->A,st->sigma,st->B,st->str);
265: } else {
266: MatShift(st->A,st->sigma);
267: }
268: }
269: /* Apply new shift */
270: if (newshift != 0.0) {
271: if (st->B) {
272: MatAXPY(st->A,-newshift,st->B,st->str);
273: } else {
274: MatShift(st->A,-newshift);
275: }
276: }
277: KSPSetOperators(st->ksp,st->A,st->A,flg);
278: break;
279: case ST_MATMODE_SHELL:
280: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
281: break;
282: default:
283: MatCopy(st->A, st->mat,SUBSET_NONZERO_PATTERN);
284: if (newshift != 0.0) {
285: if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
286: else { MatShift(st->mat,-newshift); }
287: }
288: KSPSetOperators(st->ksp,st->mat,st->mat,flg);
289: }
290: st->sigma = newshift;
291: KSPSetUp(st->ksp);
292: return(0);
293: }
297: PetscErrorCode STSetFromOptions_Cayley(ST st)
298: {
300: PetscScalar nu;
301: PetscTruth flg;
302: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
303: PC pc;
304: const PCType pctype;
305: const KSPType ksptype;
309: KSPGetPC(st->ksp,&pc);
310: KSPGetType(st->ksp,&ksptype);
311: PCGetType(pc,&pctype);
312: if (!pctype && !ksptype) {
313: if (st->shift_matrix == ST_MATMODE_SHELL) {
314: /* in shell mode use GMRES with Jacobi as the default */
315: KSPSetType(st->ksp,KSPGMRES);
316: PCSetType(pc,PCJACOBI);
317: } else {
318: /* use direct solver as default */
319: KSPSetType(st->ksp,KSPPREONLY);
320: PCSetType(pc,PCREDUNDANT);
321: }
322: }
324: PetscOptionsBegin(((PetscObject)st)->comm,((PetscObject)st)->prefix,"ST Cayley Options","ST");
325: PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
326: if (flg) {
327: STCayleySetAntishift(st,nu);
328: }
329: PetscOptionsEnd();
330: return(0);
331: }
336: PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
337: {
338: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
341: ctx->nu = newshift;
342: ctx->nu_set = PETSC_TRUE;
343: return(0);
344: }
349: /*@
350: STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
351: spectral transformation.
353: Collective on ST
355: Input Parameters:
356: + st - the spectral transformation context
357: - nu - the anti-shift
359: Options Database Key:
360: . -st_cayley_antishift - Sets the value of the anti-shift
362: Level: intermediate
364: Note:
365: In the generalized Cayley transform, the operator can be expressed as
366: OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
367: Use STSetShift() for setting sigma.
369: .seealso: STSetShift()
370: @*/
371: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
372: {
373: PetscErrorCode ierr, (*f)(ST,PetscScalar);
377: PetscObjectQueryFunction((PetscObject)st,"STCayleySetAntishift_C",(void (**)(void))&f);
378: if (f) {
379: (*f)(st,nu);
380: }
381: return(0);
382: }
386: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
387: {
389: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
392: #if !defined(PETSC_USE_COMPLEX)
393: PetscViewerASCIIPrintf(viewer," antishift: %g\n",ctx->nu);
394: #else
395: PetscViewerASCIIPrintf(viewer," antishift: %g+%g i\n",PetscRealPart(ctx->nu),PetscImaginaryPart(ctx->nu));
396: #endif
397: STView_Default(st,viewer);
398: return(0);
399: }
403: PetscErrorCode STDestroy_Cayley(ST st)
404: {
406: ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;
409: if (ctx->w2) { VecDestroy(ctx->w2); }
410: PetscFree(ctx);
411: PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","",PETSC_NULL);
412: return(0);
413: }
418: PetscErrorCode STCreate_Cayley(ST st)
419: {
421: ST_CAYLEY *ctx;
424: PetscNew(ST_CAYLEY,&ctx);
425: PetscLogObjectMemory(st,sizeof(ST_CAYLEY));
426: st->data = (void *) ctx;
428: st->ops->apply = STApply_Cayley;
429: st->ops->getbilinearform = STGetBilinearForm_Cayley;
430: st->ops->applytrans = STApplyTranspose_Cayley;
431: st->ops->postsolve = STPostSolve_Cayley;
432: st->ops->backtr = STBackTransform_Cayley;
433: st->ops->setfromoptions = STSetFromOptions_Cayley;
434: st->ops->setup = STSetUp_Cayley;
435: st->ops->setshift = STSetShift_Cayley;
436: st->ops->destroy = STDestroy_Cayley;
437: st->ops->view = STView_Cayley;
438:
439: st->checknullspace = STCheckNullSpace_Default;
441: ctx->nu = 0.0;
442: ctx->nu_set = PETSC_FALSE;
444: PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","STCayleySetAntishift_Cayley",STCayleySetAntishift_Cayley);
446: return(0);
447: }