Actual source code: stsolve.c
2: /*
3: The ST (spectral transformation) interface routines, callable by users.
4: */
6: #include src/st/stimpl.h
10: /*@
11: STApply - Applies the spectral transformation operator to a vector, for
12: instance (A - sB)^-1 B in the case of the shift-and-invert tranformation
13: and generalized eigenproblem.
15: Collective on ST and Vec
17: Input Parameters:
18: + st - the spectral transformation context
19: - x - input vector
21: Output Parameter:
22: . y - output vector
24: Level: developer
26: .seealso: STApplyB(), STApplyNoB()
27: @*/
28: PetscErrorCode STApply(ST st,Vec x,Vec y)
29: {
30: PetscErrorCode ierr;
36: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
38: if (!st->setupcalled) { STSetUp(st); }
40: PetscLogEventBegin(ST_Apply,st,x,y,0);
41: (*st->ops->apply)(st,x,y);
42: PetscLogEventEnd(ST_Apply,st,x,y,0);
43: return(0);
44: }
48: /*@
49: STApplyB - Applies the B matrix to a vector.
51: Collective on ST and Vec
53: Input Parameters:
54: + st - the spectral transformation context
55: - x - input vector
57: Output Parameter:
58: . y - output vector
60: Level: developer
62: .seealso: STApply(), STApplyNoB()
63: @*/
64: PetscErrorCode STApplyB(ST st,Vec x,Vec y)
65: {
66: PetscErrorCode ierr;
72: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
74: if (!st->setupcalled) { STSetUp(st); }
76: if (x == st->x && x->state == st->xstate) {
77: VecCopy(st->Bx, y);
78: return(0);
79: }
81: PetscLogEventBegin(ST_ApplyB,st,x,y,0);
82: (*st->ops->applyB)(st,x,y);
83: PetscLogEventEnd(ST_ApplyB,st,x,y,0);
84:
85: st->x = x;
86: st->xstate = x->state;
87: VecCopy(y,st->Bx);
88: return(0);
89: }
93: /*@
94: STApplyNoB - Applies the spectral transformation operator to a vector
95: which has already been multiplied by matrix B. For instance, this routine
96: would perform the operation y =(A - sB)^-1 x in the case of the
97: shift-and-invert tranformation and generalized eigenproblem.
99: Collective on ST and Vec
101: Input Parameters:
102: + st - the spectral transformation context
103: - x - input vector, where it is assumed that x=Bw for some vector w
105: Output Parameter:
106: . y - output vector
108: Level: developer
110: .seealso: STApply(), STApplyB()
111: @*/
112: PetscErrorCode STApplyNoB(ST st,Vec x,Vec y)
113: {
114: PetscErrorCode ierr;
120: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
122: if (!st->setupcalled) { STSetUp(st); }
124: PetscLogEventBegin(ST_ApplyNoB,st,x,y,0);
125: (*st->ops->applynoB)(st,x,y);
126: PetscLogEventEnd(ST_ApplyNoB,st,x,y,0);
127: return(0);
128: }
132: /*@
133: STNorm - Computes de norm of a vector as the square root of the inner
134: product (x,x) as defined by STInnerProduct().
136: Collective on ST and Vec
138: Input Parameters:
139: + st - the spectral transformation context
140: - x - input vector
142: Output Parameter:
143: + norm - the computed norm
145: Notes:
146: This function will usually compute the 2-norm of a vector, ||x||_2. But
147: this behaviour may be different if using a non-standard inner product changed
148: via STSetBilinearForm(). For example, if using the B-inner product for
149: positive definite B, (x,y)_B=y^H Bx, then the computed norm is ||x||_B =
150: sqrt( x^H Bx ).
152: Level: developer
154: .seealso: STInnerProduct()
155: @*/
156: PetscErrorCode STNorm(ST st,Vec x,PetscReal *norm)
157: {
158: PetscErrorCode ierr;
159: PetscScalar p;
165:
166: STInnerProduct(st,x,x,&p);
168: if (PetscAbsScalar(p)<PETSC_MACHINE_EPSILON)
169: PetscLogInfo(st,"STNorm: Zero norm, either the vector is zero or a semi-inner product is being used\n" );
171: #if defined(PETSC_USE_COMPLEX)
172: if (PetscRealPart(p)<0.0 || PetscAbsReal(PetscImaginaryPart(p))>PETSC_MACHINE_EPSILON)
173: SETERRQ(1,"STNorm: The inner product is not well defined");
174: *norm = PetscSqrtScalar(PetscRealPart(p));
175: #else
176: if (p<0.0) SETERRQ(1,"STNorm: The inner product is not well defined");
177: *norm = PetscSqrtScalar(p);
178: #endif
180: return(0);
181: }
185: /*@
186: STInnerProduct - Computes de inner product of two vectors.
188: Collective on ST and Vec
190: Input Parameters:
191: + st - the spectral transformation context
192: . x - input vector
193: - y - input vector
195: Output Parameter:
196: + p - result of the inner product
198: Notes:
199: This function will usually compute the standard dot product of vectors
200: x and y, (x,y)=y^H x. However this behaviour may be different if changed
201: via STSetBilinearForm(). This allows use of other inner products such as
202: the indefinite product y^T x for complex symmetric problems or the
203: B-inner product for positive definite B, (x,y)_B=y^H Bx.
205: Level: developer
207: .seealso: STSetBilinearForm(), STApplyB(), VecDot()
208: @*/
209: PetscErrorCode STInnerProduct(ST st,Vec x,Vec y,PetscScalar *p)
210: {
211: PetscErrorCode ierr;
218:
219: PetscLogEventBegin(ST_InnerProduct,st,x,0,0);
220: switch (st->bilinear_form) {
221: case STINNER_HERMITIAN:
222: case STINNER_SYMMETRIC:
223: VecCopy(x,st->w);
224: break;
225: case STINNER_B_HERMITIAN:
226: case STINNER_B_SYMMETRIC:
227: STApplyB(st,x,st->w);
228: break;
229: }
230: switch (st->bilinear_form) {
231: case STINNER_HERMITIAN:
232: case STINNER_B_HERMITIAN:
233: VecDot(st->w,y,p);
234: break;
235: case STINNER_SYMMETRIC:
236: case STINNER_B_SYMMETRIC:
237: VecTDot(st->w,y,p);
238: break;
239: }
240: PetscLogEventEnd(ST_InnerProduct,st,x,0,0);
241: return(0);
242: }
246: /*@
247: STSetUp - Prepares for the use of a spectral transformation.
249: Collective on ST
251: Input Parameter:
252: . st - the spectral transformation context
254: Level: advanced
256: .seealso: STCreate(), STApply(), STDestroy()
257: @*/
258: PetscErrorCode STSetUp(ST st)
259: {
260: PetscErrorCode ierr;
265: PetscLogInfo(st,"STSetUp:Setting up new ST\n");
266: if (st->setupcalled) return(0);
267: PetscLogEventBegin(ST_SetUp,st,0,0,0);
268: if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
269: if (!st->type_name) {
270: STSetType(st,STSHIFT);
271: }
272: if (st->w) { VecDestroy(st->w); }
273: if (st->Bx) { VecDestroy(st->Bx); }
274: MatGetVecs(st->A,&st->w,&st->Bx);
275: st->x = PETSC_NULL;
276: if (st->ops->setup) {
277: (*st->ops->setup)(st);
278: }
279: st->setupcalled = 1;
280: PetscLogEventEnd(ST_SetUp,st,0,0,0);
281: return(0);
282: }
286: /*
287: STPreSolve - Optional pre-solve phase, intended for any actions that
288: must be performed on the ST object before the eigensolver starts iterating.
290: Collective on ST
292: Input Parameters:
293: st - the spectral transformation context
294: eps - the eigenproblem solver context
296: Level: developer
298: Sample of Usage:
300: STPreSolve(st,eps);
301: EPSSolve(eps,its);
302: STPostSolve(st,eps);
303: */
304: PetscErrorCode STPreSolve(ST st,EPS eps)
305: {
306: PetscErrorCode ierr;
311: if (st->ops->presolve) {
312: (*st->ops->presolve)(st);
313: }
314: return(0);
315: }
319: /*
320: STPostSolve - Optional post-solve phase, intended for any actions that must
321: be performed on the ST object after the eigensolver has finished.
323: Collective on ST
325: Input Parameters:
326: st - the spectral transformation context
327: eps - the eigenproblem solver context
329: Sample of Usage:
331: STPreSolve(st,eps);
332: EPSSolve(eps,its);
333: STPostSolve(st,eps);
334: */
335: PetscErrorCode STPostSolve(ST st,EPS eps)
336: {
337: PetscErrorCode ierr;
341: if (st->ops->postsolve) {
342: (*st->ops->postsolve)(st);
343: }
345: return(0);
346: }
350: /*@
351: STBackTransform - Back-transformation phase, intended for
352: spectral transformations which require to transform the computed
353: eigenvalues back to the original eigenvalue problem.
355: Collective on ST
357: Input Parameters:
358: st - the spectral transformation context
359: eigr - real part of a computed eigenvalue
360: eigi - imaginary part of a computed eigenvalue
362: Level: developer
364: .seealso: EPSBackTransform()
365: @*/
366: PetscErrorCode STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
367: {
368: PetscErrorCode ierr;
372: if (st->ops->backtr) {
373: (*st->ops->backtr)(st,eigr,eigi);
374: }
376: return(0);
377: }