Actual source code: stsolve.c
1: /*
2: The ST (spectral transformation) interface routines, callable by users.
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
28: /*@
29: STApply - Applies the spectral transformation operator to a vector, for
30: instance (A - sB)^-1 B in the case of the shift-and-invert tranformation
31: and generalized eigenproblem.
33: Collective on ST and Vec
35: Input Parameters:
36: + st - the spectral transformation context
37: - x - input vector
39: Output Parameter:
40: . y - output vector
42: Level: developer
44: .seealso: STApplyTranspose()
45: @*/
46: PetscErrorCode STApply(ST st,Vec x,Vec y)
47: {
54: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
56: if (!st->setupcalled) { STSetUp(st); }
58: PetscLogEventBegin(ST_Apply,st,x,y,0);
59: st->applys++;
60: if (st->D) { /* with balancing */
61: VecPointwiseDivide(st->wb,x,st->D);
62: (*st->ops->apply)(st,st->wb,y);
63: VecPointwiseMult(y,y,st->D);
64: }
65: else {
66: (*st->ops->apply)(st,x,y);
67: }
68: PetscLogEventEnd(ST_Apply,st,x,y,0);
69: return(0);
70: }
74: /*@
75: STGetBilinearForm - Returns the matrix used in the bilinear form with a
76: generalized problem with semi-definite B.
78: Collective on ST and Mat
80: Input Parameters:
81: . st - the spectral transformation context
83: Output Parameter:
84: . B - output matrix
86: Note:
87: The output matrix B must be destroyed after use.
88:
89: Level: developer
90: @*/
91: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
92: {
98: (*st->ops->getbilinearform)(st,B);
99: return(0);
100: }
104: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
105: {
109: *B = st->B;
110: if (*B) {
111: PetscObjectReference((PetscObject)*B);
112: }
113: return(0);
114: }
118: /*@
119: STApplyTranspose - Applies the transpose of the operator to a vector, for
120: instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
121: and generalized eigenproblem.
123: Collective on ST and Vec
125: Input Parameters:
126: + st - the spectral transformation context
127: - x - input vector
129: Output Parameter:
130: . y - output vector
132: Level: developer
134: .seealso: STApply()
135: @*/
136: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
137: {
144: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
146: if (!st->setupcalled) { STSetUp(st); }
148: PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
149: st->applys++;
150: if (st->D) { /* with balancing */
151: VecPointwiseMult(st->wb,x,st->D);
152: (*st->ops->applytrans)(st,st->wb,y);
153: VecPointwiseDivide(y,y,st->D);
154: }
155: else {
156: (*st->ops->applytrans)(st,x,y);
157: }
158: PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
159: return(0);
160: }
164: /*@
165: STComputeExplicitOperator - Computes the explicit operator associated
166: to the eigenvalue problem with the specified spectral transformation.
168: Collective on ST
170: Input Parameter:
171: . st - the spectral transform context
173: Output Parameter:
174: . mat - the explicit operator
176: Notes:
177: This routine builds a matrix containing the explicit operator. For
178: example, in generalized problems with shift-and-invert spectral
179: transformation the result would be matrix (A - s B)^-1 B.
180:
181: This computation is done by applying the operator to columns of the
182: identity matrix. Note that the result is a dense matrix.
184: Level: advanced
186: .seealso: STApply()
187: @*/
188: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
189: {
191: Vec in,out;
192: PetscInt i,M,m,*rows,start,end;
193: PetscScalar *array,one = 1.0;
199: MatGetVecs(st->A,&in,&out);
200: VecGetSize(out,&M);
201: VecGetLocalSize(out,&m);
202: VecGetOwnershipRange(out,&start,&end);
203: PetscMalloc(m*sizeof(PetscInt),&rows);
204: for (i=0; i<m; i++) rows[i] = start + i;
206: MatCreateMPIDense(((PetscObject)st)->comm,m,m,M,M,PETSC_NULL,mat);
208: for (i=0; i<M; i++) {
209: VecSet(in,0.0);
210: VecSetValues(in,1,&i,&one,INSERT_VALUES);
211: VecAssemblyBegin(in);
212: VecAssemblyEnd(in);
214: STApply(st,in,out);
215:
216: VecGetArray(out,&array);
217: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
218: VecRestoreArray(out,&array);
219: }
220: PetscFree(rows);
221: VecDestroy(in);
222: VecDestroy(out);
223: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
224: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
225: return(0);
226: }
230: /*@
231: STSetUp - Prepares for the use of a spectral transformation.
233: Collective on ST
235: Input Parameter:
236: . st - the spectral transformation context
238: Level: advanced
240: .seealso: STCreate(), STApply(), STDestroy()
241: @*/
242: PetscErrorCode STSetUp(ST st)
243: {
249: PetscInfo(st,"Setting up new ST\n");
250: if (st->setupcalled) return(0);
251: PetscLogEventBegin(ST_SetUp,st,0,0,0);
252: if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
253: if (!((PetscObject)st)->type_name) {
254: STSetType(st,STSHIFT);
255: }
256: if (st->w) { VecDestroy(st->w); }
257: MatGetVecs(st->A,&st->w,PETSC_NULL);
258: if (st->ops->setup) {
259: (*st->ops->setup)(st);
260: }
261: st->setupcalled = 1;
262: PetscLogEventEnd(ST_SetUp,st,0,0,0);
263: return(0);
264: }
268: /*@
269: STPostSolve - Optional post-solve phase, intended for any actions that must
270: be performed on the ST object after the eigensolver has finished.
272: Collective on ST
274: Input Parameters:
275: . st - the spectral transformation context
277: Level: developer
279: .seealso: EPSSolve()
280: @*/
281: PetscErrorCode STPostSolve(ST st)
282: {
287: if (st->ops->postsolve) {
288: (*st->ops->postsolve)(st);
289: }
291: return(0);
292: }
296: /*@
297: STBackTransform - Back-transformation phase, intended for
298: spectral transformations which require to transform the computed
299: eigenvalues back to the original eigenvalue problem.
301: Collective on ST
303: Input Parameters:
304: st - the spectral transformation context
305: eigr - real part of a computed eigenvalue
306: eigi - imaginary part of a computed eigenvalue
308: Level: developer
310: .seealso: EPSBackTransform()
311: @*/
312: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
313: {
318: if (st->ops->backtr) {
319: (*st->ops->backtr)(st,n,eigr,eigi);
320: }
322: return(0);
323: }