Actual source code: cayley.c
slepc-3.5.2 2014-10-10
1: /*
2: Implements the Cayley spectral transform.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepc-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: {
39: /* standard eigenproblem: y = (A - sI)^-1 (A + tI)x */
40: /* generalized eigenproblem: y = (A - sB)^-1 (A + tB)x */
41: MatMult(st->T[0],x,st->w);
42: STMatSolve(st,st->w,y);
43: return(0);
44: }
48: PetscErrorCode STApplyTranspose_Cayley(ST st,Vec x,Vec y)
49: {
53: /* standard eigenproblem: y = (A + tI)^T (A - sI)^-T x */
54: /* generalized eigenproblem: y = (A + tB)^T (A - sB)^-T x */
55: STMatSolveTranspose(st,x,st->w);
56: MatMultTranspose(st->T[0],st->w,y);
57: return(0);
58: }
62: static PetscErrorCode MatMult_Cayley(Mat B,Vec x,Vec y)
63: {
65: ST st;
66: ST_CAYLEY *ctx;
67: PetscScalar nu;
70: MatShellGetContext(B,(void**)&st);
71: ctx = (ST_CAYLEY*)st->data;
72: nu = ctx->nu;
74: if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
76: if (st->nmat>1) {
77: /* generalized eigenproblem: y = (A + tB)x */
78: MatMult(st->A[0],x,y);
79: MatMult(st->A[1],x,ctx->w2);
80: VecAXPY(y,nu,ctx->w2);
81: } else {
82: /* standard eigenproblem: y = (A + tI)x */
83: MatMult(st->A[0],x,y);
84: VecAXPY(y,nu,x);
85: }
86: return(0);
87: }
91: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
92: {
96: STSetUp(st);
97: *B = st->T[0];
98: PetscObjectReference((PetscObject)*B);
99: return(0);
100: }
104: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
105: {
106: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
107: PetscInt j;
108: #if !defined(PETSC_USE_COMPLEX)
109: PetscScalar t,i,r;
110: #endif
113: #if !defined(PETSC_USE_COMPLEX)
114: for (j=0;j<n;j++) {
115: if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
116: else {
117: r = eigr[j];
118: i = eigi[j];
119: r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
120: i = - st->sigma * i - ctx->nu * i;
121: t = i * i + r * (r - 2.0) + 1.0;
122: eigr[j] = r / t;
123: eigi[j] = i / t;
124: }
125: }
126: #else
127: for (j=0;j<n;j++) {
128: eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
129: }
130: #endif
131: return(0);
132: }
136: PetscErrorCode STPostSolve_Cayley(ST st)
137: {
141: if (st->shift_matrix == ST_MATMODE_INPLACE) {
142: if (st->nmat>1) {
143: MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
144: } else {
145: MatShift(st->A[0],st->sigma);
146: }
147: st->Astate[0] = ((PetscObject)st->A[0])->state;
148: st->setupcalled = 0;
149: }
150: return(0);
151: }
155: PetscErrorCode STSetUp_Cayley(ST st)
156: {
158: PetscInt n,m;
159: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
162: /* if the user did not set the shift, use the target value */
163: if (!st->sigma_set) st->sigma = st->defsigma;
165: if (!ctx->nu_set) ctx->nu = st->sigma;
166: if (ctx->nu == 0.0 && st->sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)st),1,"Values of shift and antishift cannot be zero simultaneously");
168: /* T[0] = A+nu*B */
169: if (st->shift_matrix==ST_MATMODE_INPLACE) {
170: MatGetLocalSize(st->A[0],&n,&m);
171: MatCreateShell(PetscObjectComm((PetscObject)st),n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,&st->T[0]);
172: MatShellSetOperation(st->T[0],MATOP_MULT,(void(*)(void))MatMult_Cayley);
173: PetscLogObjectParent((PetscObject)st,(PetscObject)st->T[0]);
174: } else {
175: STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PETSC_TRUE,&st->T[0]);
176: }
178: /* T[1] = A-sigma*B */
179: STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,&st->T[1]);
180: st->P = st->T[1];
181: PetscObjectReference((PetscObject)st->P);
182: if (st->nmat>1) {
183: VecDestroy(&ctx->w2);
184: MatGetVecs(st->A[1],&ctx->w2,NULL);
185: PetscLogObjectParent((PetscObject)st,(PetscObject)ctx->w2);
186: }
187: if (!st->ksp) { STGetKSP(st,&st->ksp); }
188: KSPSetOperators(st->ksp,st->P,st->P);
189: KSPSetUp(st->ksp);
190: return(0);
191: }
195: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
196: {
198: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
201: if (newshift==0.0 && (!ctx->nu_set || (ctx->nu_set && ctx->nu==0.0))) SETERRQ(PetscObjectComm((PetscObject)st),1,"Values of shift and antishift cannot be zero simultaneously");
203: /* Nothing to be done if STSetUp has not been called yet */
204: if (!st->setupcalled) return(0);
206: if (!ctx->nu_set) {
207: if (st->shift_matrix!=ST_MATMODE_INPLACE) {
208: STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
209: }
210: ctx->nu = newshift;
211: }
212: STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,&st->T[1]);
213: if (st->P!=st->T[1]) {
214: MatDestroy(&st->P);
215: st->P = st->T[1];
216: PetscObjectReference((PetscObject)st->P);
217: }
218: KSPSetOperators(st->ksp,st->P,st->P);
219: KSPSetUp(st->ksp);
220: return(0);
221: }
225: PetscErrorCode STSetFromOptions_Cayley(ST st)
226: {
228: PetscScalar nu;
229: PetscBool flg;
230: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
231: PC pc;
232: PCType pctype;
233: KSPType ksptype;
236: if (!st->ksp) { STGetKSP(st,&st->ksp); }
237: KSPGetPC(st->ksp,&pc);
238: KSPGetType(st->ksp,&ksptype);
239: PCGetType(pc,&pctype);
240: if (!pctype && !ksptype) {
241: if (st->shift_matrix == ST_MATMODE_SHELL) {
242: /* in shell mode use GMRES with Jacobi as the default */
243: KSPSetType(st->ksp,KSPGMRES);
244: PCSetType(pc,PCJACOBI);
245: } else {
246: /* use direct solver as default */
247: KSPSetType(st->ksp,KSPPREONLY);
248: PCSetType(pc,PCREDUNDANT);
249: }
250: }
252: PetscOptionsHead("ST Cayley Options");
253: PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
254: if (flg) {
255: STCayleySetAntishift(st,nu);
256: }
257: PetscOptionsTail();
258: return(0);
259: }
263: static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
264: {
266: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
269: if (st->setupcalled && st->shift_matrix!=ST_MATMODE_INPLACE) {
270: STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
271: }
272: ctx->nu = newshift;
273: ctx->nu_set = PETSC_TRUE;
274: return(0);
275: }
279: /*@
280: STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
281: spectral transformation.
283: Logically Collective on ST
285: Input Parameters:
286: + st - the spectral transformation context
287: - nu - the anti-shift
289: Options Database Key:
290: . -st_cayley_antishift - Sets the value of the anti-shift
292: Level: intermediate
294: Note:
295: In the generalized Cayley transform, the operator can be expressed as
296: OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
297: Use STSetShift() for setting sigma.
299: .seealso: STSetShift(), STCayleyGetAntishift()
300: @*/
301: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
302: {
308: PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
309: return(0);
310: }
313: static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
314: {
315: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
318: *nu = ctx->nu;
319: return(0);
320: }
324: /*@
325: STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
326: spectral transformation.
328: Not Collective
330: Input Parameter:
331: . st - the spectral transformation context
333: Output Parameter:
334: . nu - the anti-shift
336: Level: intermediate
338: .seealso: STGetShift(), STCayleySetAntishift()
339: @*/
340: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
341: {
347: PetscTryMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
348: return(0);
349: }
353: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
354: {
356: char str[50];
357: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
360: SlepcSNPrintfScalar(str,50,ctx->nu,PETSC_FALSE);
361: PetscViewerASCIIPrintf(viewer," Cayley: antishift: %s\n",str);
362: return(0);
363: }
367: PetscErrorCode STReset_Cayley(ST st)
368: {
370: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
373: VecDestroy(&ctx->w2);
374: return(0);
375: }
379: PetscErrorCode STDestroy_Cayley(ST st)
380: {
384: PetscFree(st->data);
385: PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL);
386: PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL);
387: return(0);
388: }
392: PETSC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
393: {
395: ST_CAYLEY *ctx;
398: PetscNewLog(st,&ctx);
399: st->data = (void*)ctx;
401: st->ops->apply = STApply_Cayley;
402: st->ops->getbilinearform = STGetBilinearForm_Cayley;
403: st->ops->applytrans = STApplyTranspose_Cayley;
404: st->ops->postsolve = STPostSolve_Cayley;
405: st->ops->backtransform = STBackTransform_Cayley;
406: st->ops->setfromoptions = STSetFromOptions_Cayley;
407: st->ops->setup = STSetUp_Cayley;
408: st->ops->setshift = STSetShift_Cayley;
409: st->ops->destroy = STDestroy_Cayley;
410: st->ops->reset = STReset_Cayley;
411: st->ops->view = STView_Cayley;
412: st->ops->checknullspace = STCheckNullSpace_Default;
413: PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley);
414: PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley);
415: return(0);
416: }