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: int STApply(ST st,Vec x,Vec y)
 29: {
 30:   int        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: int STApplyB(ST st,Vec x,Vec y)
 65: {
 66:   int        ierr;

 72:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

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

 76:   PetscLogEventBegin(ST_ApplyB,st,x,y,0);
 77:   (*st->ops->applyB)(st,x,y);
 78:   PetscLogEventEnd(ST_ApplyB,st,x,y,0);
 79:   return(0);
 80: }

 84: /*@
 85:    STApplyNoB - Applies the spectral transformation operator to a vector 
 86:    which has already been multiplied by matrix B. For instance, this routine
 87:    would perform the operation y =(A - sB)^-1 x in the case of the 
 88:    shift-and-invert tranformation and generalized eigenproblem.

 90:    Collective on ST and Vec

 92:    Input Parameters:
 93: +  st - the spectral transformation context
 94: -  x  - input vector, where it is assumed that x=Bw for some vector w

 96:    Output Parameter:
 97: .  y - output vector

 99:    Level: developer

101: .seealso: STApply(), STApplyB()
102: @*/
103: int STApplyNoB(ST st,Vec x,Vec y)
104: {
105:   int        ierr;
106:   PetscTruth isSinv;

112:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");

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

116:   PetscTypeCompare((PetscObject)st,STSINV,&isSinv);
117:   if (!isSinv) { SETERRQ(1,"Function only available in Shift-and-Invert"); }

119:   PetscLogEventBegin(ST_ApplyNoB,st,x,y,0);
120:   (*st->ops->applynoB)(st,x,y);
121:   PetscLogEventEnd(ST_ApplyNoB,st,x,y,0);
122:   return(0);
123: }

127: /*@
128:    STSetUp - Prepares for the use of a spectral transformation.

130:    Collective on ST

132:    Input Parameter:
133: .  st - the spectral transformation context

135:    Level: advanced

137: .seealso: STCreate(), STApply(), STDestroy()
138: @*/
139: int STSetUp(ST st)
140: {


146:   PetscLogInfo(st,"STSetUp:Setting up new ST\n");
147:   if (st->setupcalled) return(0);
148:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
149:   if (!st->vec) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Vector must be set first");}
150:   if (!st->A) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
151:   if (!st->type_name) {
152:     STSetType(st,STSHIFT);
153:   }
154:   if (st->ops->setup) {
155:     (*st->ops->setup)(st);
156:   }
157:   st->setupcalled = 1;
158:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
159:   return(0);
160: }

164: /*
165:    STPreSolve - Optional pre-solve phase, intended for any actions that 
166:    must be performed on the ST object before the eigensolver starts iterating.

168:    Collective on ST

170:    Input Parameters:
171:    st  - the spectral transformation context
172:    eps - the eigenproblem solver context

174:    Level: developer

176:    Sample of Usage:

178:     STPreSolve(st,eps);
179:     EPSSolve(eps,its);
180:     STPostSolve(st,eps);
181: */
182: int STPreSolve(ST st,EPS eps)
183: {
184:   int         ierr;


189:   if (st->ops->presolve) {
190:     (*st->ops->presolve)(st);
191:   }
192:   return(0);
193: }

197: /*
198:    STPostSolve - Optional post-solve phase, intended for any actions that must 
199:    be performed on the ST object after the eigensolver has finished.

201:    Collective on ST

203:    Input Parameters:
204:    st  - the spectral transformation context
205:    eps - the eigenproblem solver context

207:    Sample of Usage:

209:     STPreSolve(st,eps);
210:     EPSSolve(eps,its);
211:     STPostSolve(st,eps);
212: */
213: int STPostSolve(ST st,EPS eps)
214: {
215:   int         ierr;

219:   if (st->ops->postsolve) {
220:     (*st->ops->postsolve)(st);
221:   }

223:   return(0);
224: }

228: /*
229:    STBackTransform - Optional back-transformation phase, intended for 
230:    spectral transformation which require to transform the computed 
231:    eigenvalues back to the original eigenvalue problem.

233:    Collective on ST

235:    Input Parameters:
236:    st   - the spectral transformation context
237:    eigr - real part of a computed eigenvalue
238:    eigi - imaginary part of a computed eigenvalue

240:    Level: developer

242: .seealso: EPSBackTransform()
243: */
244: int STBackTransform(ST st,PetscScalar* eigr,PetscScalar* eigi)
245: {
246:   int         ierr;

250:   if (st->ops->backtr) {
251:     (*st->ops->backtr)(st,eigr,eigi);
252:   }

254:   return(0);
255: }

259: int STDefaultApplyB(ST st,Vec x,Vec y)
260: {
261:   int        ierr;

264:   if( st->B ) {
265:     MatMult( st->B, x, y );
266:   }
267:   else {
268:     VecCopy( x, y );
269:   }
270:   return(0);
271: }