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-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

 58:   if (!st->ops->apply) SETERRQ(PetscObjectComm((PetscObject)st),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:   } 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:    Not collective, though a parallel Mat may be returned

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

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

 86:    Notes:
 87:    The output matrix B must be destroyed after use. It will be NULL in
 88:    case of standard eigenproblems.

 90:    Level: developer
 91: @*/
 92: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
 93: {

 99:   if (!st->A) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Matrices must be set first");
100:   (*st->ops->getbilinearform)(st,B);
101:   return(0);
102: }

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

111:   if (st->nmat==1) *B = NULL;
112:   else {
113:     *B = st->A[1];
114:     PetscObjectReference((PetscObject)*B);
115:   }
116:   return(0);
117: }

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

126:    Collective on ST and Vec

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

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

135:    Level: developer

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

147:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");

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

151:   if (!st->ops->applytrans) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST does not have applytrans");
152:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
153:   st->applys++;
154:   if (st->D) { /* with balancing */
155:     VecPointwiseMult(st->wb,x,st->D);
156:     (*st->ops->applytrans)(st,st->wb,y);
157:     VecPointwiseDivide(y,y,st->D);
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. This is analogous to MatComputeExplicitOperator().

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;
198:   PetscMPIInt       size;

203:   if (!st->A) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Matrices must be set first");
204:   MPI_Comm_size(PetscObjectComm((PetscObject)st),&size);

206:   MatGetVecs(st->A[0],&in,&out);
207:   VecGetSize(out,&M);
208:   VecGetLocalSize(out,&m);
209:   VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
210:   VecGetOwnershipRange(out,&start,&end);
211:   PetscMalloc(m*sizeof(PetscInt),&rows);
212:   for (i=0;i<m;i++) rows[i] = start + i;

214:   MatCreate(PetscObjectComm((PetscObject)st),mat);
215:   MatSetSizes(*mat,m,m,M,M);
216:   if (size == 1) {
217:     MatSetType(*mat,MATSEQDENSE);
218:     MatSeqDenseSetPreallocation(*mat,NULL);
219:   } else {
220:     MatSetType(*mat,MATMPIAIJ);
221:     MatMPIAIJSetPreallocation(*mat,m,NULL,M-m,NULL);
222:   }

224:   for (i=0;i<M;i++) {
225:     VecSet(in,0.0);
226:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
227:     VecAssemblyBegin(in);
228:     VecAssemblyEnd(in);

230:     STApply(st,in,out);

232:     VecGetArrayRead(out,&array);
233:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
234:     VecRestoreArrayRead(out,&array);
235:   }
236:   PetscFree(rows);
237:   VecDestroy(&in);
238:   VecDestroy(&out);
239:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
240:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
241:   return(0);
242: }

246: /*@
247:    STSetUp - Prepares for the use of a spectral transformation.

249:    Collective on ST

251:    Input Parameter:
252: .  st - the spectral transformation context

254:    Level: advanced

256: .seealso: STCreate(), STApply(), STDestroy()
257: @*/
258: PetscErrorCode STSetUp(ST st)
259: {
260:   PetscInt       i,n,k;

265:   if (!st->A) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Matrices must be set first");
266:   if (st->setupcalled) return(0);
267:   PetscInfo(st,"Setting up new ST\n");
268:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
269:   if (!((PetscObject)st)->type_name) {
270:     STSetType(st,STSHIFT);
271:   }
272:   if (!st->T) {
273:     PetscMalloc(PetscMax(2,st->nmat)*sizeof(Mat),&st->T);
274:     PetscLogObjectMemory(st,PetscMax(2,st->nmat)*sizeof(Mat));
275:   }
276:   for (i=0;i<PetscMax(2,st->nmat);i++) st->T[i] = NULL;
277:   if (!st->w) {
278:     MatGetVecs(st->A[0],&st->w,NULL);
279:     PetscLogObjectParent(st,st->w);
280:   }
281:   if (st->D) {
282:     MatGetLocalSize(st->A[0],NULL,&n);
283:     VecGetLocalSize(st->D,&k);
284:     if (n != k) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Balance matrix has wrong dimension %D (should be %D)",k,n);
285:     if (!st->wb) {
286:       VecDuplicate(st->D,&st->wb);
287:       PetscLogObjectParent(st,st->wb);
288:     }
289:   }
290:   if (st->ops->setup) { (*st->ops->setup)(st); }
291:   st->setupcalled = 1;
292:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
293:   return(0);
294: }

298: /*
299:    Computes a generalized AXPY operation on the A[:] matrices provided by the user,
300:    and stores the result in one of the T[:] matrices. The computation is different
301:    depending on whether the ST has two or three matrices. Also, in quadratic problems
302:    this function is used to compute two matrices, those corresponding to degree 1 and 2.

304:    Builds matrix in T[k] as follows:
305:    deg=1: T[k] = A+alpha*B (if st->nmat==2) or B+2*alpha*C (if st->nmat==3)
306:    deg=2: T[k] = A+alpha*B+alpha*alpha*C
307: */
308: PetscErrorCode STMatGAXPY_Private(ST st,PetscScalar alpha,PetscScalar beta,PetscInt deg,PetscInt k,PetscBool initial)
309: {
311:   PetscScalar    gamma;
312:   PetscInt       matIdx[3],t,i;

315:   if (st->nmat==3 && deg==1) t = 1;
316:   else t = 0;
317:   switch (st->shift_matrix) {
318:   case ST_MATMODE_INPLACE:
319:     if (initial) {
320:       PetscObjectReference((PetscObject)st->A[t]);
321:       st->T[k] = st->A[t];
322:       gamma = alpha;
323:     } else gamma = alpha-beta;
324:     if (gamma != 0.0) {
325:       if (st->nmat>1) {
326:         MatAXPY(st->T[k],gamma,st->A[t+1],st->str);
327:         if (st->nmat==3 && deg==2) {
328:           MatAXPY(st->T[k],gamma*gamma,st->A[2],st->str);
329:         }
330:       } else {
331:         MatShift(st->T[k],gamma);
332:       }
333:     }
334:     break;
335:   case ST_MATMODE_SHELL:
336:     if (initial) {
337:       if (st->nmat>1) {
338:         for (i=0;i<=deg;i++) {
339:           matIdx[i] = t+i;
340:         }
341:         STMatShellCreate(st,alpha,deg+1,matIdx,&st->T[k]);
342:       } else {
343:         STMatShellCreate(st,alpha,deg,NULL,&st->T[k]);
344:       }
345:       PetscLogObjectParent(st,st->T[k]);
346:     } else {
347:       STMatShellShift(st->T[k],alpha);
348:     }
349:     break;
350:   default:
351:     if (alpha == 0.0) {
352:       if (!initial) {
353:         MatDestroy(&st->T[k]);
354:       }
355:       PetscObjectReference((PetscObject)st->A[t]);
356:       st->T[k] = st->A[t];
357:     } else {
358:       if (initial) {
359:         MatDuplicate(st->A[t],MAT_COPY_VALUES,&st->T[k]);
360:         PetscLogObjectParent(st,st->T[k]);
361:       } else {
362:         if (beta==0.0) {
363:           MatDestroy(&st->T[k]);
364:           MatDuplicate(st->A[t],MAT_COPY_VALUES,&st->T[k]);
365:           PetscLogObjectParent(st,st->T[k]);
366:         } else {
367:           MatCopy(st->A[t],st->T[k],SAME_NONZERO_PATTERN);
368:         }
369:       }
370:       if (st->nmat>1) {
371:         MatAXPY(st->T[k],alpha,st->A[t+1],st->str);
372:         if (st->nmat==3 && deg==2) {
373:           MatAXPY(st->T[k],alpha*alpha,st->A[2],st->str);
374:         }
375:       } else {
376:         MatShift(st->T[k],alpha);
377:       }
378:     }
379:   }
380:   STMatSetHermitian(st,st->T[k]);
381:   return(0);
382: }

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

390:    Collective on ST

392:    Input Parameters:
393: .  st  - the spectral transformation context

395:    Level: developer

397: .seealso: EPSSolve()
398: @*/
399: PetscErrorCode STPostSolve(ST st)
400: {

405:   if (st->ops->postsolve) {
406:     (*st->ops->postsolve)(st);
407:   }
408:   return(0);
409: }

413: /*@
414:    STBackTransform - Back-transformation phase, intended for
415:    spectral transformations which require to transform the computed
416:    eigenvalues back to the original eigenvalue problem.

418:    Not Collective

420:    Input Parameters:
421:    st   - the spectral transformation context
422:    eigr - real part of a computed eigenvalue
423:    eigi - imaginary part of a computed eigenvalue

425:    Level: developer
426: @*/
427: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
428: {

433:   if (st->ops->backtransform) {
434:     (*st->ops->backtransform)(st,n,eigr,eigi);
435:   }
436:   return(0);
437: }