Actual source code: stsolve.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:     The ST (spectral transformation) interface routines, callable by users.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, 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: {

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

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

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

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

 81:    Collective on ST and Vec

 83:    Input Parameters:
 84: +  st - the spectral transformation context
 85: -  x  - input vector

 87:    Output Parameter:
 88: .  y - output vector

 90:    Level: developer

 92: .seealso: STApply()
 93: @*/
 94: PetscErrorCode STApplyTranspose(ST st,Vec x,Vec y)
 95: {

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

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

108:   if (!st->ops->applytrans) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST does not have applytrans");
109:   PetscLogEventBegin(ST_ApplyTranspose,st,x,y,0);
110:   st->applys++;
111:   if (st->D) { /* with balancing */
112:     VecPointwiseMult(st->wb,x,st->D);
113:     (*st->ops->applytrans)(st,st->wb,y);
114:     VecPointwiseDivide(y,y,st->D);
115:   } else {
116:     (*st->ops->applytrans)(st,x,y);
117:   }
118:   PetscLogEventEnd(ST_ApplyTranspose,st,x,y,0);
119:   return(0);
120: }

124: /*@
125:    STGetBilinearForm - Returns the matrix used in the bilinear form with a
126:    generalized problem with semi-definite B.

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

130:    Input Parameters:
131: .  st - the spectral transformation context

133:    Output Parameter:
134: .  B - output matrix

136:    Notes:
137:    The output matrix B must be destroyed after use. It will be NULL in
138:    case of standard eigenproblems.

140:    Level: developer
141: @*/
142: PetscErrorCode STGetBilinearForm(ST st,Mat *B)
143: {

150:   STCheckMatrices(st,1);
151:   (*st->ops->getbilinearform)(st,B);
152:   return(0);
153: }

157: PetscErrorCode STGetBilinearForm_Default(ST st,Mat *B)
158: {

162:   if (st->nmat==1) *B = NULL;
163:   else {
164:     *B = st->A[1];
165:     PetscObjectReference((PetscObject)*B);
166:   }
167:   return(0);
168: }

172: /*@
173:    STComputeExplicitOperator - Computes the explicit operator associated
174:    to the eigenvalue problem with the specified spectral transformation.

176:    Collective on ST

178:    Input Parameter:
179: .  st - the spectral transform context

181:    Output Parameter:
182: .  mat - the explicit operator

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

189:    This computation is done by applying the operator to columns of the
190:    identity matrix. This is analogous to MatComputeExplicitOperator().

192:    Level: advanced

194: .seealso: STApply()
195: @*/
196: PetscErrorCode STComputeExplicitOperator(ST st,Mat *mat)
197: {
198:   PetscErrorCode    ierr;
199:   Vec               in,out;
200:   PetscInt          i,M,m,*rows,start,end;
201:   const PetscScalar *array;
202:   PetscScalar       one = 1.0;
203:   PetscMPIInt       size;

208:   STCheckMatrices(st,1);
209:   if (st->nmat>2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Can only be used with 1 or 2 matrices");
210:   MPI_Comm_size(PetscObjectComm((PetscObject)st),&size);

212:   MatGetVecs(st->A[0],&in,&out);
213:   VecGetSize(out,&M);
214:   VecGetLocalSize(out,&m);
215:   VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
216:   VecGetOwnershipRange(out,&start,&end);
217:   PetscMalloc1(m,&rows);
218:   for (i=0;i<m;i++) rows[i] = start + i;

220:   MatCreate(PetscObjectComm((PetscObject)st),mat);
221:   MatSetSizes(*mat,m,m,M,M);
222:   if (size == 1) {
223:     MatSetType(*mat,MATSEQDENSE);
224:     MatSeqDenseSetPreallocation(*mat,NULL);
225:   } else {
226:     MatSetType(*mat,MATMPIAIJ);
227:     MatMPIAIJSetPreallocation(*mat,m,NULL,M-m,NULL);
228:   }

230:   for (i=0;i<M;i++) {
231:     VecSet(in,0.0);
232:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
233:     VecAssemblyBegin(in);
234:     VecAssemblyEnd(in);

236:     STApply(st,in,out);

238:     VecGetArrayRead(out,&array);
239:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
240:     VecRestoreArrayRead(out,&array);
241:   }
242:   PetscFree(rows);
243:   VecDestroy(&in);
244:   VecDestroy(&out);
245:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
246:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
247:   return(0);
248: }

252: /*@
253:    STSetUp - Prepares for the use of a spectral transformation.

255:    Collective on ST

257:    Input Parameter:
258: .  st - the spectral transformation context

260:    Level: advanced

262: .seealso: STCreate(), STApply(), STDestroy()
263: @*/
264: PetscErrorCode STSetUp(ST st)
265: {
266:   PetscInt       i,n,k;

271:   STCheckMatrices(st,1);
272:   if (st->setupcalled) return(0);
273:   PetscInfo(st,"Setting up new ST\n");
274:   PetscLogEventBegin(ST_SetUp,st,0,0,0);
275:   if (!((PetscObject)st)->type_name) {
276:     STSetType(st,STSHIFT);
277:   }
278:   if (!st->T) {
279:     PetscMalloc(PetscMax(2,st->nmat)*sizeof(Mat),&st->T);
280:     PetscLogObjectMemory((PetscObject)st,PetscMax(2,st->nmat)*sizeof(Mat));
281:     for (i=0;i<PetscMax(2,st->nmat);i++) st->T[i] = NULL;
282:   } else {
283:     for (i=0;i<PetscMax(2,st->nmat);i++) {
284:       MatDestroy(&st->T[i]);
285:     }
286:   }
287:   MatDestroy(&st->P);
288:   if (!st->w) {
289:     MatGetVecs(st->A[0],&st->w,NULL);
290:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->w);
291:   }
292:   if (st->D) {
293:     MatGetLocalSize(st->A[0],NULL,&n);
294:     VecGetLocalSize(st->D,&k);
295:     if (n != k) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Balance matrix has wrong dimension %D (should be %D)",k,n);
296:     if (!st->wb) {
297:       VecDuplicate(st->D,&st->wb);
298:       PetscLogObjectParent((PetscObject)st,(PetscObject)st->wb);
299:     }
300:   }
301:   if (st->ops->setup) { (*st->ops->setup)(st); }
302:   st->setupcalled = 1;
303:   PetscLogEventEnd(ST_SetUp,st,0,0,0);
304:   return(0);
305: }

309: /*
310:    Computes coefficients for the transformed polynomial,
311:    and stores the result in argument S.

313:    alpha - value of the parameter of the transformed polynomial
314:    beta - value of the previous shift (only used in inplace mode)
315:    k - number of A matrices involved in the computation
316:    coeffs - coefficients of the expansion
317:    initial - true if this is the first time (only relevant for shell mode)
318: */
319: PetscErrorCode STMatMAXPY_Private(ST st,PetscScalar alpha,PetscScalar beta,PetscInt k,PetscScalar *coeffs,PetscBool initial,Mat *S)
320: {
322:   PetscInt       *matIdx=NULL,nmat,i,ini=-1;
323:   PetscScalar    t=1.0,ta,gamma;
324:   PetscBool      nz=PETSC_FALSE;

327:   nmat = st->nmat-k;
328:   switch (st->shift_matrix) {
329:   case ST_MATMODE_INPLACE:
330:     if (st->nmat>2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"ST_MATMODE_INPLACE not supported for polynomial eigenproblems");
331:     if (initial) {
332:       PetscObjectReference((PetscObject)st->A[0]);
333:       *S = st->A[0];
334:       gamma = alpha;
335:     } else gamma = alpha-beta;
336:     if (gamma != 0.0) {
337:       if (st->nmat>1) {
338:         MatAXPY(*S,gamma,st->A[1],st->str);
339:       } else {
340:         MatShift(*S,gamma);
341:       }
342:     }
343:     break;
344:   case ST_MATMODE_SHELL:
345:     if (initial) {
346:       if (st->nmat>2) {
347:         PetscMalloc(nmat*sizeof(PetscInt),&matIdx);
348:         for (i=0;i<nmat;i++) matIdx[i] = k+i;
349:       }
350:       STMatShellCreate(st,alpha,nmat,matIdx,coeffs,S);
351:       PetscLogObjectParent((PetscObject)st,(PetscObject)*S);
352:       if (st->nmat>2) { PetscFree(matIdx); }
353:     } else {
354:       STMatShellShift(*S,alpha);
355:     }
356:     break;
357:   case ST_MATMODE_COPY:
358:     if (coeffs) {
359:       for (i=0;i<nmat && ini==-1;i++) {
360:         if (coeffs[i]!=0.0) ini = i;
361:         else t *= alpha;
362:       }
363:       if (coeffs[ini] != 1.0) nz = PETSC_TRUE;
364:       for (i=ini+1;i<nmat&&!nz;i++) if (coeffs[i]!=0.0) nz = PETSC_TRUE;
365:     } else { nz = PETSC_TRUE; ini = 0; }
366:     if ((alpha == 0.0 || !nz) && t==1.0) {
367:       MatDestroy(S);
368:       PetscObjectReference((PetscObject)st->A[k+ini]);
369:       *S = st->A[k+ini];
370:     } else {
371:       if (*S && *S!=st->A[k+ini]) {
372:         MatCopy(st->A[k+ini],*S,DIFFERENT_NONZERO_PATTERN);
373:       } else {
374:         MatDestroy(S);
375:         MatDuplicate(st->A[k+ini],MAT_COPY_VALUES,S);
376:         PetscLogObjectParent((PetscObject)st,(PetscObject)*S);
377:       }
378:       if (coeffs && coeffs[ini]!=1.0) {
379:         MatScale(*S,coeffs[ini]);
380:       }
381:       for (i=ini+k+1;i<PetscMax(2,st->nmat);i++) {
382:         t *= alpha;
383:         ta = t;
384:         if (coeffs) ta *= coeffs[i-k];
385:         if (ta!=0.0) {
386:           if (st->nmat>1) {
387:             MatAXPY(*S,ta,st->A[i],st->str);
388:           } else {
389:             MatShift(*S,ta);
390:           }
391:         }
392:       }
393:     }
394:   }
395:   STMatSetHermitian(st,*S);
396:   return(0);
397: }

401: /*
402:    Computes the values of the coefficients required by STMatMAXPY_Private
403:    for the case of monomial basis.
404: */
405: PetscErrorCode STCoeffs_Monomial(ST st, PetscScalar *coeffs)
406: {
407:   PetscInt  k,i,ini,inip;  

410:   /* Compute binomial coefficients */
411:   ini = (st->nmat*(st->nmat-1))/2;
412:   for (i=0;i<st->nmat;i++) coeffs[ini+i]=1.0;
413:   for (k=st->nmat-1;k>=1;k--) {
414:     inip = ini+1;
415:     ini = (k*(k-1))/2;
416:     coeffs[ini] = 1.0;
417:     for (i=1;i<k;i++) coeffs[ini+i] = coeffs[ini+i-1]+coeffs[inip+i-1];
418:   }
419:   return(0);
420: }

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

428:    Collective on ST

430:    Input Parameters:
431: .  st  - the spectral transformation context

433:    Level: developer

435: .seealso: EPSSolve()
436: @*/
437: PetscErrorCode STPostSolve(ST st)
438: {

444:   if (st->ops->postsolve) {
445:     (*st->ops->postsolve)(st);
446:   }
447:   return(0);
448: }

452: /*@
453:    STBackTransform - Back-transformation phase, intended for
454:    spectral transformations which require to transform the computed
455:    eigenvalues back to the original eigenvalue problem.

457:    Not Collective

459:    Input Parameters:
460:    st   - the spectral transformation context
461:    eigr - real part of a computed eigenvalue
462:    eigi - imaginary part of a computed eigenvalue

464:    Level: developer
465: @*/
466: PetscErrorCode STBackTransform(ST st,PetscInt n,PetscScalar* eigr,PetscScalar* eigi)
467: {

473:   if (st->ops->backtransform) {
474:     (*st->ops->backtransform)(st,n,eigr,eigi);
475:   }
476:   return(0);
477: }

481: /*@
482:    STMatSetUp - Build the preconditioner matrix used in STMatSolve().

484:    Collective on ST

486:    Input Parameters:
487: +  st     - the spectral transformation context
488: .  sigma  - the shift
489: -  coeffs - the coefficients

491:    Note:
492:    This function is not intended to be called by end users, but by SLEPc
493:    solvers that use ST. It builds matrix st->P as follows, then calls KSPSetUp().
494: .vb
495:     If (coeffs):  st->P = Sum_{i=0:nmat-1} coeffs[i]*sigma^i*A_i.
496:     else          st->P = Sum_{i=0:nmat-1} sigma^i*A_i
497: .ve

499:    Level: developer

501: .seealso: STMatSolve()
502: @*/
503: PetscErrorCode STMatSetUp(ST st,PetscScalar sigma,PetscScalar *coeffs)
504: {

511:   STCheckMatrices(st,1);

513:   PetscLogEventBegin(ST_MatSetUp,st,0,0,0);
514:   STMatMAXPY_Private(st,sigma,0.0,0,coeffs,PETSC_TRUE,&st->P);
515:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
516:   KSPSetOperators(st->ksp,st->P,st->P);
517:   KSPSetUp(st->ksp);
518:   PetscLogEventEnd(ST_MatSetUp,st,0,0,0);
519:   return(0);
520: }