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-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
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: (*st->ops->apply)(st,x,y);
61: PetscLogEventEnd(ST_Apply,st,x,y,0);
62: return(0);
63: }
67: /*@
68: STGetBilinearForm - Returns the matrix used in the bilinear form with a semi-definite generalised problem.
70: Collective on ST and Mat
72: Input Parameters:
73: . st - the spectral transformation context
75: Output Parameter:
76: . B - output matrix
78: Note:
79: The output matrix B must be destroyed after use.
80:
81: Level: developer
82: @*/
83: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
84: {
90: (*st->ops->getbilinearform)(st,B);
91: return(0);
92: }
96: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
97: {
101: *B = st->B;
102: if (*B) {
103: PetscObjectReference((PetscObject)*B);
104: }
105: return(0);
106: }
110: /*@
111: STApplyTranspose - Applies the transpose of the operator to a vector, for
112: instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
113: and generalized eigenproblem.
115: Collective on ST and Vec
117: Input Parameters:
118: + st - the spectral transformation context
119: - x - input vector
121: Output Parameter:
122: . y - output vector
124: Level: developer
126: .seealso: STApply()
127: @*/
128: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
129: {
136: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
138: if (!st->setupcalled) { STSetUp(st); }
140: PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
141: (*st->ops->applytrans)(st,x,y);
142: PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
143: return(0);
144: }
148: /*@
149: STComputeExplicitOperator - Computes the explicit operator associated
150: to the eigenvalue problem with the specified spectral transformation.
152: Collective on ST
154: Input Parameter:
155: . st - the spectral transform context
157: Output Parameter:
158: . mat - the explicit operator
160: Notes:
161: This routine builds a matrix containing the explicit operator. For
162: example, in generalized problems with shift-and-invert spectral
163: transformation the result would be matrix (A - s B)^-1 B.
164:
165: This computation is done by applying the operator to columns of the
166: identity matrix. Note that the result is a dense matrix.
168: Level: advanced
170: .seealso: STApply()
171: @*/
172: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
173: {
175: Vec in,out;
176: PetscInt i,M,m,*rows,start,end;
177: PetscScalar *array,one = 1.0;
183: MatGetVecs(st->A,&in,&out);
184: VecGetSize(out,&M);
185: VecGetLocalSize(out,&m);
186: VecGetOwnershipRange(out,&start,&end);
187: PetscMalloc(m*sizeof(PetscInt),&rows);
188: for (i=0; i<m; i++) rows[i] = start + i;
190: MatCreateMPIDense(((PetscObject)st)->comm,m,m,M,M,PETSC_NULL,mat);
192: for (i=0; i<M; i++) {
193: VecSet(in,0.0);
194: VecSetValues(in,1,&i,&one,INSERT_VALUES);
195: VecAssemblyBegin(in);
196: VecAssemblyEnd(in);
198: STApply(st,in,out);
199:
200: VecGetArray(out,&array);
201: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
202: VecRestoreArray(out,&array);
203: }
204: PetscFree(rows);
205: VecDestroy(in);
206: VecDestroy(out);
207: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
208: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
209: return(0);
210: }
214: /*@
215: STSetUp - Prepares for the use of a spectral transformation.
217: Collective on ST
219: Input Parameter:
220: . st - the spectral transformation context
222: Level: advanced
224: .seealso: STCreate(), STApply(), STDestroy()
225: @*/
226: PetscErrorCode STSetUp(ST st)
227: {
233: PetscInfo(st,"Setting up new ST\n");
234: if (st->setupcalled) return(0);
235: PetscLogEventBegin(ST_SetUp,st,0,0,0);
236: if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
237: if (!((PetscObject)st)->type_name) {
238: STSetType(st,STSHIFT);
239: }
240: if (st->w) { VecDestroy(st->w); }
241: MatGetVecs(st->A,&st->w,PETSC_NULL);
242: if (st->ops->setup) {
243: (*st->ops->setup)(st);
244: }
245: st->setupcalled = 1;
246: PetscLogEventEnd(ST_SetUp,st,0,0,0);
247: return(0);
248: }
252: /*@
253: STPostSolve - Optional post-solve phase, intended for any actions that must
254: be performed on the ST object after the eigensolver has finished.
256: Collective on ST
258: Input Parameters:
259: . st - the spectral transformation context
261: Level: developer
263: .seealso: EPSSolve()
264: @*/
265: PetscErrorCode STPostSolve(ST st)
266: {
271: if (st->ops->postsolve) {
272: (*st->ops->postsolve)(st);
273: }
275: return(0);
276: }
280: /*@
281: STBackTransform - Back-transformation phase, intended for
282: spectral transformations which require to transform the computed
283: eigenvalues back to the original eigenvalue problem.
285: Collective on ST
287: Input Parameters:
288: st - the spectral transformation context
289: eigr - real part of a computed eigenvalue
290: eigi - imaginary part of a computed eigenvalue
292: Level: developer
294: .seealso: EPSBackTransform()
295: @*/
296: PetscErrorCode STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
297: {
302: if (st->ops->backtr) {
303: (*st->ops->backtr)(st,eigr,eigi);
304: }
306: return(0);
307: }