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-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/st/stimpl.h
17: /*@
18: STApply - Applies the spectral transformation operator to a vector, for
19: instance (A - sB)^-1 B in the case of the shift-and-invert tranformation
20: and generalized eigenproblem.
22: Collective on ST and Vec
24: Input Parameters:
25: + st - the spectral transformation context
26: - x - input vector
28: Output Parameter:
29: . y - output vector
31: Level: developer
33: .seealso: STApplyTranspose()
34: @*/
35: PetscErrorCode STApply(ST st,Vec x,Vec y)
36: {
43: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
45: if (!st->setupcalled) { STSetUp(st); }
47: PetscLogEventBegin(ST_Apply,st,x,y,0);
48: st->applys++;
49: (*st->ops->apply)(st,x,y);
50: PetscLogEventEnd(ST_Apply,st,x,y,0);
51: return(0);
52: }
56: /*@
57: STGetBilinearForm - Returns the matrix used in the bilinear form with a semi-definite generalised problem.
59: Collective on ST and Mat
61: Input Parameters:
62: . st - the spectral transformation context
64: Output Parameter:
65: . B - output matrix
67: Level: developer
68: @*/
69: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
70: {
76: (*st->ops->getbilinearform)(st,B);
77: return(0);
78: }
82: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
83: {
85: *B = st->B;
86: return(0);
87: }
91: /*@
92: STApplyTranspose - Applies the transpose of the operator to a vector, for
93: instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
94: and generalized eigenproblem.
96: Collective on ST and Vec
98: Input Parameters:
99: + st - the spectral transformation context
100: - x - input vector
102: Output Parameter:
103: . y - output vector
105: Level: developer
107: .seealso: STApply()
108: @*/
109: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
110: {
117: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
119: if (!st->setupcalled) { STSetUp(st); }
121: PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
122: (*st->ops->applytrans)(st,x,y);
123: PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
124: return(0);
125: }
129: /*@
130: STComputeExplicitOperator - Computes the explicit operator associated
131: to the eigenvalue problem with the specified spectral transformation.
133: Collective on ST
135: Input Parameter:
136: . st - the spectral transform context
138: Output Parameter:
139: . mat - the explicit operator
141: Notes:
142: This routine builds a matrix containing the explicit operator. For
143: example, in generalized problems with shift-and-invert spectral
144: transformation the result would be matrix (A - s B)^-1 B.
145:
146: This computation is done by applying the operator to columns of the
147: identity matrix. Note that the result is a dense matrix.
149: Level: advanced
151: .seealso: STApply()
152: @*/
153: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
154: {
156: Vec in,out;
157: PetscInt i,M,m,*rows,start,end;
158: PetscScalar *array,one = 1.0;
164: MatGetVecs(st->A,&in,&out);
165: VecGetSize(out,&M);
166: VecGetLocalSize(out,&m);
167: VecGetOwnershipRange(out,&start,&end);
168: PetscMalloc(m*sizeof(int),&rows);
169: for (i=0; i<m; i++) rows[i] = start + i;
171: MatCreateMPIDense(st->comm,m,m,M,M,PETSC_NULL,mat);
173: for (i=0; i<M; i++) {
174: VecSet(in,0.0);
175: VecSetValues(in,1,&i,&one,INSERT_VALUES);
176: VecAssemblyBegin(in);
177: VecAssemblyEnd(in);
179: STApply(st,in,out);
180:
181: VecGetArray(out,&array);
182: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
183: VecRestoreArray(out,&array);
184: }
185: PetscFree(rows);
186: VecDestroy(in);
187: VecDestroy(out);
188: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
189: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
190: return(0);
191: }
195: /*@
196: STSetUp - Prepares for the use of a spectral transformation.
198: Collective on ST
200: Input Parameter:
201: . st - the spectral transformation context
203: Level: advanced
205: .seealso: STCreate(), STApply(), STDestroy()
206: @*/
207: PetscErrorCode STSetUp(ST st)
208: {
214: PetscInfo(st,"Setting up new ST\n");
215: if (st->setupcalled) return(0);
216: PetscLogEventBegin(ST_SetUp,st,0,0,0);
217: if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
218: if (!st->type_name) {
219: STSetType(st,STSHIFT);
220: }
221: if (st->w) { VecDestroy(st->w); }
222: MatGetVecs(st->A,&st->w,PETSC_NULL);
223: if (st->ops->setup) {
224: (*st->ops->setup)(st);
225: }
226: st->setupcalled = 1;
227: PetscLogEventEnd(ST_SetUp,st,0,0,0);
228: return(0);
229: }
233: /*@
234: STPostSolve - Optional post-solve phase, intended for any actions that must
235: be performed on the ST object after the eigensolver has finished.
237: Collective on ST
239: Input Parameters:
240: . st - the spectral transformation context
242: Level: developer
244: .seealso: EPSSolve()
245: @*/
246: PetscErrorCode STPostSolve(ST st)
247: {
252: if (st->ops->postsolve) {
253: (*st->ops->postsolve)(st);
254: }
256: return(0);
257: }
261: /*@
262: STBackTransform - Back-transformation phase, intended for
263: spectral transformations which require to transform the computed
264: eigenvalues back to the original eigenvalue problem.
266: Collective on ST
268: Input Parameters:
269: st - the spectral transformation context
270: eigr - real part of a computed eigenvalue
271: eigi - imaginary part of a computed eigenvalue
273: Level: developer
275: .seealso: EPSBackTransform()
276: @*/
277: PetscErrorCode STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
278: {
283: if (st->ops->backtr) {
284: (*st->ops->backtr)(st,eigr,eigi);
285: }
287: return(0);
288: }