Actual source code: epsdefault.c

slepc-3.11.0 2019-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    This file contains some simple default routines for common operations
 12: */

 14: #include <slepc/private/epsimpl.h>   /*I "slepceps.h" I*/
 15: #include <slepcvec.h>

 17: PetscErrorCode EPSBackTransform_Default(EPS eps)
 18: {

 22:   STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
 23:   return(0);
 24: }

 26: /*
 27:   EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
 28:   using purification for generalized eigenproblems.
 29:  */
 30: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
 31: {
 33:   PetscInt       i;
 34:   PetscScalar    dot;
 35:   PetscReal      norm;
 36:   PetscBool      iscayley;
 37:   Mat            B;
 38:   Vec            w,x;

 41:   if (eps->purify) {
 42:     EPS_Purify(eps,eps->nconv);
 43:     for (i=0;i<eps->nconv;i++) {
 44:       BVNormColumn(eps->V,i,NORM_2,&norm);
 45:       BVScaleColumn(eps->V,i,1.0/norm);
 46:     }
 47:   } else {
 48:     /* In the case of Cayley transform, eigenvectors need to be B-normalized */
 49:     PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
 50:     if (iscayley && eps->isgeneralized) {
 51:       STGetMatrix(eps->st,1,&B);
 52:       MatCreateVecs(B,NULL,&w);
 53:       for (i=0;i<eps->nconv;i++) {
 54:         BVGetColumn(eps->V,i,&x);
 55:         MatMult(B,x,w);
 56:         VecDot(w,x,&dot);
 57:         VecScale(x,1.0/PetscSqrtScalar(dot));
 58:         BVRestoreColumn(eps->V,i,&x);
 59:       }
 60:       VecDestroy(&w);
 61:     }
 62:   }
 63:   return(0);
 64: }

 66: /*
 67:   EPSComputeVectors_Indefinite - similar to the Schur version but
 68:   for indefinite problems
 69:  */
 70: PetscErrorCode EPSComputeVectors_Indefinite(EPS eps)
 71: {
 73:   PetscInt       n,i;
 74:   Mat            X;
 75:   Vec            v;
 76: #if !defined(PETSC_USE_COMPLEX)
 77:   Vec            v1;
 78:   PetscScalar    tmp;
 79:   PetscReal      norm,normi;
 80: #endif

 83:   DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
 84:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
 85:   DSGetMat(eps->ds,DS_MAT_X,&X);
 86:   BVSetActiveColumns(eps->V,0,n);
 87:   BVMultInPlace(eps->V,X,0,n);
 88:   MatDestroy(&X);

 90:   /* purification */
 91:   if (eps->purify) { EPS_Purify(eps,eps->nconv); }

 93:   /* normalization */
 94:   for (i=0;i<n;i++) {
 95: #if !defined(PETSC_USE_COMPLEX)
 96:     if (eps->eigi[i] != 0.0) {
 97:       BVGetColumn(eps->V,i,&v);
 98:       BVGetColumn(eps->V,i+1,&v1);
 99:       VecNorm(v,NORM_2,&norm);
100:       VecNorm(v1,NORM_2,&normi);
101:       tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
102:       VecScale(v,tmp);
103:       VecScale(v1,tmp);
104:       BVRestoreColumn(eps->V,i,&v);
105:       BVRestoreColumn(eps->V,i+1,&v1);
106:       i++;
107:     } else
108: #endif
109:     {
110:       BVGetColumn(eps->V,i,&v);
111:       VecNormalize(v,NULL);
112:       BVRestoreColumn(eps->V,i,&v);
113:     }
114:   }
115:   return(0);
116: }

118: /*
119:   EPSComputeVectors_Twosided - Adjust left eigenvectors in generalized problems: y = B^-* y.
120:  */
121: PetscErrorCode EPSComputeVectors_Twosided(EPS eps)
122: {
124:   PetscReal      norm;
125:   PetscInt       i;
126:   Vec            w,y;
127: #if !defined(PETSC_USE_COMPLEX)
128:   Vec            z;
129: #endif

132:   if (!eps->twosided || !eps->isgeneralized) return(0);
133:   EPSSetWorkVecs(eps,1);
134:   w = eps->work[0];
135:   for (i=0;i<eps->nconv;i++) {
136:     BVCopyVec(eps->W,i,w);
137:     VecConjugate(w);
138:     BVGetColumn(eps->W,i,&y);
139:     STMatSolveTranspose(eps->st,w,y);
140:     VecConjugate(y);
141:     BVRestoreColumn(eps->W,i,&y);
142:   }
143:   /* normalize */
144:   for (i=0;i<eps->nconv;i++) {
145: #if !defined(PETSC_USE_COMPLEX)
146:     if (eps->eigi[i]!=0.0) {   /* first eigenvalue of a complex conjugate pair */
147:       BVGetColumn(eps->W,i,&y);
148:       BVGetColumn(eps->W,i+1,&z);
149:       VecNormalizeComplex(y,z,PETSC_TRUE,NULL);
150:       BVRestoreColumn(eps->W,i,&y);
151:       BVRestoreColumn(eps->W,i+1,&z);
152:       i++;
153:     } else   /* real eigenvalue */
154: #endif
155:     {
156:       BVNormColumn(eps->W,i,NORM_2,&norm);
157:       BVScaleColumn(eps->W,i,1.0/norm);
158:     }
159:   }
160:   return(0);
161: }

163: /*
164:   EPSComputeVectors_Schur - Compute eigenvectors from the vectors
165:   provided by the eigensolver. This version is intended for solvers
166:   that provide Schur vectors. Given the partial Schur decomposition
167:   OP*V=V*T, the following steps are performed:
168:       1) compute eigenvectors of T: T*Z=Z*D
169:       2) compute eigenvectors of OP: X=V*Z
170:  */
171: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
172: {
174:   PetscInt       n,i;
175:   Mat            Z;
176:   Vec            z,v;
177: #if !defined(PETSC_USE_COMPLEX)
178:   Vec            v1;
179:   PetscScalar    tmp;
180:   PetscReal      norm,normi;
181: #endif

184:   if (eps->ishermitian) {
185:     if (eps->isgeneralized && !eps->ispositive) {
186:        EPSComputeVectors_Indefinite(eps);
187:     } else {
188:       EPSComputeVectors_Hermitian(eps);
189:     }
190:     return(0);
191:   }
192:   DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);

194:   /* right eigenvectors */
195:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);

197:   /* V = V * Z */
198:   DSGetMat(eps->ds,DS_MAT_X,&Z);
199:   BVSetActiveColumns(eps->V,0,n);
200:   BVMultInPlace(eps->V,Z,0,n);
201:   MatDestroy(&Z);

203:   /* Purify eigenvectors */
204:   if (eps->purify) { EPS_Purify(eps,eps->nconv); }

206:   /* Fix eigenvectors if balancing was used */
207:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
208:     for (i=0;i<n;i++) {
209:       BVGetColumn(eps->V,i,&z);
210:       VecPointwiseDivide(z,z,eps->D);
211:       BVRestoreColumn(eps->V,i,&z);
212:     }
213:   }

215:   /* normalize eigenvectors (when using purification or balancing) */
216:   if (eps->purify || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
217:     for (i=0;i<n;i++) {
218: #if !defined(PETSC_USE_COMPLEX)
219:       if (eps->eigi[i] != 0.0) {
220:         BVGetColumn(eps->V,i,&v);
221:         BVGetColumn(eps->V,i+1,&v1);
222:         VecNorm(v,NORM_2,&norm);
223:         VecNorm(v1,NORM_2,&normi);
224:         tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
225:         VecScale(v,tmp);
226:         VecScale(v1,tmp);
227:         BVRestoreColumn(eps->V,i,&v);
228:         BVRestoreColumn(eps->V,i+1,&v1);
229:         i++;
230:       } else
231: #endif
232:       {
233:         BVGetColumn(eps->V,i,&v);
234:         VecNormalize(v,NULL);
235:         BVRestoreColumn(eps->V,i,&v);
236:       }
237:     }
238:   }

240:   /* left eigenvectors from eps->dsts */
241:   if (eps->twosided) {
242:     DSVectors(eps->dsts,DS_MAT_X,NULL,NULL);
243:     /* W = W * Z */
244:     DSGetMat(eps->dsts,DS_MAT_X,&Z);
245:     BVSetActiveColumns(eps->W,0,n);
246:     BVMultInPlace(eps->W,Z,0,n);
247:     MatDestroy(&Z);
248:     EPSComputeVectors_Twosided(eps);
249: #if !defined(PETSC_USE_COMPLEX)
250:     for (i=0;i<eps->nconv-1;i++) {
251:       if (eps->eigi[i] != 0.0) {
252:         if (eps->eigi[i] > 0.0) { BVScaleColumn(eps->W,i+1,-1.0); }
253:         i++;
254:       }
255:     }
256: #endif
257:   }
258:   return(0);
259: }

261: /*@
262:    EPSSetWorkVecs - Sets a number of work vectors into an EPS object.

264:    Collective on EPS

266:    Input Parameters:
267: +  eps - eigensolver context
268: -  nw  - number of work vectors to allocate

270:    Developers Note:
271:    This is SLEPC_EXTERN because it may be required by user plugin EPS
272:    implementations.

274:    Level: developer
275: @*/
276: PetscErrorCode EPSSetWorkVecs(EPS eps,PetscInt nw)
277: {
279:   Vec            t;

284:   if (nw <= 0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"nw must be > 0: nw = %D",nw);
285:   if (eps->nwork < nw) {
286:     VecDestroyVecs(eps->nwork,&eps->work);
287:     eps->nwork = nw;
288:     BVGetColumn(eps->V,0,&t);
289:     VecDuplicateVecs(t,nw,&eps->work);
290:     BVRestoreColumn(eps->V,0,&t);
291:     PetscLogObjectParents(eps,nw,eps->work);
292:   }
293:   return(0);
294: }

296: /*
297:   EPSSetWhichEigenpairs_Default - Sets the default value for which,
298:   depending on the ST.
299:  */
300: PetscErrorCode EPSSetWhichEigenpairs_Default(EPS eps)
301: {
302:   PetscBool      target;

306:   PetscObjectTypeCompareAny((PetscObject)eps->st,&target,STSINVERT,STCAYLEY,"");
307:   if (target) eps->which = EPS_TARGET_MAGNITUDE;
308:   else eps->which = EPS_LARGEST_MAGNITUDE;
309:   return(0);
310: }

312: /*
313:   EPSConvergedRelative - Checks convergence relative to the eigenvalue.
314: */
315: PetscErrorCode EPSConvergedRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
316: {
317:   PetscReal w;

320:   w = SlepcAbsEigenvalue(eigr,eigi);
321:   *errest = res/w;
322:   return(0);
323: }

325: /*
326:   EPSConvergedAbsolute - Checks convergence absolutely.
327: */
328: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
329: {
331:   *errest = res;
332:   return(0);
333: }

335: /*
336:   EPSConvergedNorm - Checks convergence relative to the eigenvalue and
337:   the matrix norms.
338: */
339: PetscErrorCode EPSConvergedNorm(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
340: {
341:   PetscReal w;

344:   w = SlepcAbsEigenvalue(eigr,eigi);
345:   *errest = res / (eps->nrma + w*eps->nrmb);
346:   return(0);
347: }

349: /*@C
350:    EPSStoppingBasic - Default routine to determine whether the outer eigensolver
351:    iteration must be stopped.

353:    Collective on EPS

355:    Input Parameters:
356: +  eps    - eigensolver context obtained from EPSCreate()
357: .  its    - current number of iterations
358: .  max_it - maximum number of iterations
359: .  nconv  - number of currently converged eigenpairs
360: .  nev    - number of requested eigenpairs
361: -  ctx    - context (not used here)

363:    Output Parameter:
364: .  reason - result of the stopping test

366:    Notes:
367:    A positive value of reason indicates that the iteration has finished successfully
368:    (converged), and a negative value indicates an error condition (diverged). If
369:    the iteration needs to be continued, reason must be set to EPS_CONVERGED_ITERATING
370:    (zero).

372:    EPSStoppingBasic() will stop if all requested eigenvalues are converged, or if
373:    the maximum number of iterations has been reached.

375:    Use EPSSetStoppingTest() to provide your own test instead of using this one.

377:    Level: advanced

379: .seealso: EPSSetStoppingTest(), EPSConvergedReason, EPSGetConvergedReason()
380: @*/
381: PetscErrorCode EPSStoppingBasic(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
382: {

386:   *reason = EPS_CONVERGED_ITERATING;
387:   if (nconv >= nev) {
388:     PetscInfo2(eps,"Linear eigensolver finished successfully: %D eigenpairs converged at iteration %D\n",nconv,its);
389:     *reason = EPS_CONVERGED_TOL;
390:   } else if (its >= max_it) {
391:     *reason = EPS_DIVERGED_ITS;
392:     PetscInfo1(eps,"Linear eigensolver iteration reached maximum number of iterations (%D)\n",its);
393:   }
394:   return(0);
395: }

397: /*
398:   EPSComputeRitzVector - Computes the current Ritz vector.

400:   Simple case (complex scalars or real scalars with Zi=NULL):
401:     x = V*Zr  (V is a basis of nv vectors, Zr has length nv)

403:   Split case:
404:     x = V*Zr  y = V*Zi  (Zr and Zi have length nv)
405: */
406: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,BV V,Vec x,Vec y)
407: {
409:   PetscInt       l,k;
410:   PetscReal      norm;
411: #if !defined(PETSC_USE_COMPLEX)
412:   Vec            z;
413:   PetscReal      normi;
414:   PetscScalar    tmp;
415: #endif

418:   /* compute eigenvector */
419:   BVGetActiveColumns(V,&l,&k);
420:   BVSetActiveColumns(V,0,k);
421:   BVMultVec(V,1.0,0.0,x,Zr);

423:   /* purify eigenvector if necessary */
424:   if (eps->purify) {
425:     STApply(eps->st,x,y);
426:     if (eps->ishermitian) {
427:       BVNormVec(eps->V,y,NORM_2,&norm);
428:     } else {
429:       VecNorm(y,NORM_2,&norm);
430:     }
431:     VecScale(y,1.0/norm);
432:     VecCopy(y,x);
433:   }
434:   /* fix eigenvector if balancing is used */
435:   if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
436:     VecPointwiseDivide(x,x,eps->D);
437:   }
438: #if !defined(PETSC_USE_COMPLEX)
439:   /* compute imaginary part of eigenvector */
440:   if (Zi) {
441:     BVMultVec(V,1.0,0.0,y,Zi);
442:     if (eps->ispositive) {
443:       BVCreateVec(V,&z);
444:       STApply(eps->st,y,z);
445:       VecNorm(z,NORM_2,&norm);
446:       VecScale(z,1.0/norm);
447:       VecCopy(z,y);
448:       VecDestroy(&z);
449:     }
450:     if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
451:       VecPointwiseDivide(y,y,eps->D);
452:     }
453:   } else
454: #endif
455:   { VecSet(y,0.0); }

457:   /* normalize eigenvectors (when using balancing) */
458:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
459: #if !defined(PETSC_USE_COMPLEX)
460:     if (Zi) {
461:       VecNorm(x,NORM_2,&norm);
462:       VecNorm(y,NORM_2,&normi);
463:       tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
464:       VecScale(x,tmp);
465:       VecScale(y,tmp);
466:     } else
467: #endif
468:     {
469:       VecNormalize(x,NULL);
470:     }
471:   }
472:   BVSetActiveColumns(V,l,k);
473:   return(0);
474: }

476: /*
477:   EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
478:   diagonal matrix to be applied for balancing in non-Hermitian problems.
479: */
480: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
481: {
482:   Vec               z,p,r;
483:   PetscInt          i,j;
484:   PetscReal         norma;
485:   PetscScalar       *pz,*pD;
486:   const PetscScalar *pr,*pp;
487:   PetscRandom       rand;
488:   PetscErrorCode    ierr;

491:   EPSSetWorkVecs(eps,3);
492:   BVGetRandomContext(eps->V,&rand);
493:   r = eps->work[0];
494:   p = eps->work[1];
495:   z = eps->work[2];
496:   VecSet(eps->D,1.0);

498:   for (j=0;j<eps->balance_its;j++) {

500:     /* Build a random vector of +-1's */
501:     VecSetRandom(z,rand);
502:     VecGetArray(z,&pz);
503:     for (i=0;i<eps->nloc;i++) {
504:       if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
505:       else pz[i]=1.0;
506:     }
507:     VecRestoreArray(z,&pz);

509:     /* Compute p=DA(D\z) */
510:     VecPointwiseDivide(r,z,eps->D);
511:     STApply(eps->st,r,p);
512:     VecPointwiseMult(p,p,eps->D);
513:     if (eps->balance == EPS_BALANCE_TWOSIDE) {
514:       if (j==0) {
515:         /* Estimate the matrix inf-norm */
516:         VecAbs(p);
517:         VecMax(p,NULL,&norma);
518:       }
519:       /* Compute r=D\(A'Dz) */
520:       VecPointwiseMult(z,z,eps->D);
521:       STApplyTranspose(eps->st,z,r);
522:       VecPointwiseDivide(r,r,eps->D);
523:     }

525:     /* Adjust values of D */
526:     VecGetArrayRead(r,&pr);
527:     VecGetArrayRead(p,&pp);
528:     VecGetArray(eps->D,&pD);
529:     for (i=0;i<eps->nloc;i++) {
530:       if (eps->balance == EPS_BALANCE_TWOSIDE) {
531:         if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
532:           pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
533:       } else {
534:         if (pp[i]!=0.0) pD[i] /= PetscAbsScalar(pp[i]);
535:       }
536:     }
537:     VecRestoreArrayRead(r,&pr);
538:     VecRestoreArrayRead(p,&pp);
539:     VecRestoreArray(eps->D,&pD);
540:   }
541:   return(0);
542: }