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-2011, Universitat 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>            /*I "slepcst.h" I*/

 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(((PetscObject)st)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");

 56:   if (!st->setupcalled) { STSetUp(st); }

 58:   if (!st->ops->apply) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST does not have apply");
 59:   PetscLogEventBegin(ST_Apply,st,x,y,0);
 60:   st->applys++;
 61:   if (st->D) { /* with balancing */
 62:     VecPointwiseDivide(st->wb,x,st->D);
 63:     (*st->ops->apply)(st,st->wb,y);
 64:     VecPointwiseMult(y,y,st->D);
 65:   }
 66:   else {
 67:     (*st->ops->apply)(st,x,y);
 68:   }
 69:   PetscLogEventEnd(ST_Apply,st,x,y,0);
 70:   return(0);
 71: }

 75: /*@
 76:    STGetBilinearForm - Returns the matrix used in the bilinear form with a 
 77:    generalized problem with semi-definite B.

 79:    Not collective, though a parallel Mat may be returned

 81:    Input Parameters:
 82: .  st - the spectral transformation context

 84:    Output Parameter:
 85: .  B - output matrix

 87:    Notes:
 88:    The output matrix B must be destroyed after use. It will be PETSC_NULL in
 89:    case of standard eigenproblems.
 90:    
 91:    Level: developer
 92: @*/
 93: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
 94: {

100:   (*st->ops->getbilinearform)(st,B);
101:   return(0);
102: }

106: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
107: {

111:   *B = st->B;
112:   if (*B) {
113:      PetscObjectReference((PetscObject)*B);
114:   }
115:   return(0);
116: }

120: /*@
121:    STApplyTranspose - Applies the transpose of the operator to a vector, for
122:    instance B^T(A - sB)^-T in the case of the shift-and-invert tranformation
123:    and generalized eigenproblem.

125:    Collective on ST and Vec

127:    Input Parameters:
128: +  st - the spectral transformation context
129: -  x  - input vector

131:    Output Parameter:
132: .  y - output vector

134:    Level: developer

136: .seealso: STApply()
137: @*/
138: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
139: {

146:   if (x == y) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");

148:   if (!st->setupcalled) { STSetUp(st); }

150:   if (!st->ops->applytrans) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST does not have applytrans");
151:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
152:   st->applys++;
153:   if (st->D) { /* with balancing */
154:     VecPointwiseMult(st->wb,x,st->D);
155:     (*st->ops->applytrans)(st,st->wb,y);
156:     VecPointwiseDivide(y,y,st->D);
157:   }
158:   else {
159:     (*st->ops->applytrans)(st,x,y);
160:   }
161:   PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
162:   return(0);
163: }

167: /*@
168:    STComputeExplicitOperator - Computes the explicit operator associated
169:    to the eigenvalue problem with the specified spectral transformation.  

171:    Collective on ST

173:    Input Parameter:
174: .  st - the spectral transform context

176:    Output Parameter:
177: .  mat - the explicit operator

179:    Notes:
180:    This routine builds a matrix containing the explicit operator. For 
181:    example, in generalized problems with shift-and-invert spectral
182:    transformation the result would be matrix (A - s B)^-1 B.

184:    This computation is done by applying the operator to columns of the 
185:    identity matrix. Note that the result is a dense matrix.

187:    Level: advanced

189: .seealso: STApply()   
190: @*/
191: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
192: {
193:   PetscErrorCode    ierr;
194:   Vec               in,out;
195:   PetscInt          i,M,m,*rows,start,end;
196:   const PetscScalar *array;
197:   PetscScalar       one = 1.0;


203:   MatGetVecs(st->A,&in,&out);
204:   VecGetSize(out,&M);
205:   VecGetLocalSize(out,&m);
206:   VecGetOwnershipRange(out,&start,&end);
207:   PetscMalloc(m*sizeof(PetscInt),&rows);
208:   for (i=0; i<m; i++) rows[i] = start + i;

210:   MatCreateMPIDense(((PetscObject)st)->comm,m,m,M,M,PETSC_NULL,mat);

212:   for (i=0; i<M; i++) {
213:     VecSet(in,0.0);
214:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
215:     VecAssemblyBegin(in);
216:     VecAssemblyEnd(in);

218:     STApply(st,in,out);
219: 
220:     VecGetArrayRead(out,&array);
221:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
222:     VecRestoreArrayRead(out,&array);
223:   }
224:   PetscFree(rows);
225:   VecDestroy(&in);
226:   VecDestroy(&out);
227:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
228:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
229:   return(0);
230: }

234: /*@
235:    STSetUp - Prepares for the use of a spectral transformation.

237:    Collective on ST

239:    Input Parameter:
240: .  st - the spectral transformation context

242:    Level: advanced

244: .seealso: STCreate(), STApply(), STDestroy()
245: @*/
246: PetscErrorCode STSetUp(ST st)
247: {
248:   PetscInt       n,k;

253:   if (!st->A) {SETERRQ(((PetscObject)st)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
254:   if (st->setupcalled) return(0);
255:   PetscInfo(st,"Setting up new ST\n");
256:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
257:   if (!((PetscObject)st)->type_name) {
258:     STSetType(st,STSHIFT);
259:   }
260:   VecDestroy(&st->w);
261:   MatGetVecs(st->A,&st->w,PETSC_NULL);
262:   if (st->D) {
263:     MatGetLocalSize(st->A,PETSC_NULL,&n);
264:     VecGetLocalSize(st->D,&k);
265:     if (n != k) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Balance matrix has wrong dimension %D (should be %D)",k,n);
266:     VecDestroy(&st->wb);
267:     VecDuplicate(st->D,&st->wb);
268:   }
269:   if (st->ops->setup) { (*st->ops->setup)(st); }
270:   st->setupcalled = 1;
271:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
272:   return(0);
273: }

277: /*@
278:    STPostSolve - Optional post-solve phase, intended for any actions that must 
279:    be performed on the ST object after the eigensolver has finished.

281:    Collective on ST

283:    Input Parameters:
284: .  st  - the spectral transformation context

286:    Level: developer

288: .seealso: EPSSolve()
289: @*/
290: PetscErrorCode STPostSolve(ST st)
291: {

296:   if (st->ops->postsolve) {
297:     (*st->ops->postsolve)(st);
298:   }
299:   return(0);
300: }

304: /*@
305:    STBackTransform - Back-transformation phase, intended for 
306:    spectral transformations which require to transform the computed 
307:    eigenvalues back to the original eigenvalue problem.

309:    Not Collective

311:    Input Parameters:
312:    st   - the spectral transformation context
313:    eigr - real part of a computed eigenvalue
314:    eigi - imaginary part of a computed eigenvalue

316:    Level: developer

318: .seealso: EPSBackTransform()
319: @*/
320: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
321: {

326:   if (st->ops->backtr) {
327:     (*st->ops->backtr)(st,n,eigr,eigi);
328:   }
329:   return(0);
330: }