Actual source code: ks-slice.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur with spectrum slicing for symmetric eigenproblems

  7:    References:

  9:        [1] R.G. Grimes et al., "A shifted block Lanczos algorithm for
 10:            solving sparse symmetric generalized eigenproblems", SIAM J.
 11:            Matrix Anal. Appl. 15(1):228-272, 1994.

 13:        [2] C. Campos and J.E. Roman, "Spectrum slicing strategies based
 14:            on restarted Lanczos methods", Numer. Algor. 60(2):279-295,
 15:            2012.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

 21:    This file is part of SLEPc.

 23:    SLEPc is free software: you can redistribute it and/or modify it under  the
 24:    terms of version 3 of the GNU Lesser General Public License as published by
 25:    the Free Software Foundation.

 27:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 28:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 29:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 30:    more details.

 32:    You  should have received a copy of the GNU Lesser General  Public  License
 33:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 34:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35: */

 37: #include <slepc-private/epsimpl.h>
 38:  #include krylovschur.h

 42: /*
 43:   EPSAllocateSolutionSlice - Allocate memory storage for common variables such
 44:   as eigenvalues and eigenvectors. The argument extra is used for methods
 45:   that require a working basis slightly larger than ncv.
 46: */
 47: PetscErrorCode EPSAllocateSolutionSlice(EPS eps,PetscInt extra)
 48: {
 50:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 51:   PetscInt        requested;
 52:   PetscReal       eta;
 53:   PetscLogDouble  cnt;
 54:   BVType          type;
 55:   BVOrthogType    orthog_type;
 56:   BVOrthogRefineType orthog_ref;
 57:   Mat             matrix;
 58:   Vec             t;
 59:   EPS_SR          sr = ctx->sr;

 62:   requested = ctx->ncv + extra;

 64:   /* allocate space for eigenvalues and friends */
 65:   PetscMalloc4(requested,&sr->eigr,requested,&sr->eigi,requested,&sr->errest,requested,&sr->perm);
 66:   cnt = 2*requested*sizeof(PetscScalar) + 2*requested*sizeof(PetscReal) + requested*sizeof(PetscInt);
 67:   PetscLogObjectMemory((PetscObject)eps,cnt);

 69:   /* allocate sr->V and transfer options from eps->V */
 70:   BVCreate(PetscObjectComm((PetscObject)eps),&sr->V);
 71:   PetscLogObjectParent((PetscObject)eps,(PetscObject)sr->V);
 72:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
 73:   if (!((PetscObject)(eps->V))->type_name) {
 74:     BVSetType(sr->V,BVSVEC);
 75:   } else {
 76:     BVGetType(eps->V,&type);
 77:     BVSetType(sr->V,type);
 78:   }
 79:   STMatGetVecs(eps->st,&t,NULL);
 80:   BVSetSizesFromVec(sr->V,t,requested);
 81:   VecDestroy(&t);
 82:   EPS_SetInnerProduct(eps);
 83:   BVGetMatrix(eps->V,&matrix,NULL);
 84:   BVSetMatrix(sr->V,matrix,PETSC_FALSE);
 85:   BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta);
 86:   BVSetOrthogonalization(sr->V,orthog_type,orthog_ref,eta);
 87:   return(0);
 88: }

 92: PetscErrorCode EPSSetUp_KrylovSchur_Slice(EPS eps)
 93: {
 94:   PetscErrorCode  ierr;
 95:   PetscBool       issinv;
 96:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 97:   EPS_SR          sr;
 98:   KSP             ksp;
 99:   PC              pc;
100:   Mat             F;

103:   if (eps->intb >= PETSC_MAX_REAL && eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
104:   if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Must define a computational interval when using EPS_ALL");
105:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
106:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
107:   if (!((PetscObject)(eps->st))->type_name) { /* default to shift-and-invert */
108:     STSetType(eps->st,STSINVERT);
109:   }
110:   PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
111:   if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
112:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
113:   if (!eps->max_it) eps->max_it = 100;
114:   if (ctx->nev==1) ctx->nev = 40;  /* nev not set, use default value */
115:   if (ctx->nev<10) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
116:   eps->ops->backtransform = NULL;

118:   /* create spectrum slicing context and initialize it */
119:   EPSReset_KrylovSchur(eps);
120:   PetscNewLog(eps,&sr);
121:   ctx->sr = sr;
122:   sr->itsKs = 0;
123:   sr->nleap = 0;
124:   sr->nMAXCompl = ctx->nev/4;
125:   sr->iterCompl = eps->max_it/4;
126:   sr->sPres = NULL;
127:   sr->nS = 0;

129:   /* check presence of ends and finding direction */
130:   if ((eps->inta > PETSC_MIN_REAL && eps->inta != 0.0) || eps->intb >= PETSC_MAX_REAL) {
131:     sr->int0 = eps->inta;
132:     sr->int1 = eps->intb;
133:     sr->dir = 1;
134:     if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
135:       sr->hasEnd = PETSC_FALSE;
136:       sr->inertia1 = eps->n;
137:     } else sr->hasEnd = PETSC_TRUE;
138:   } else {
139:     sr->int0 = eps->intb;
140:     sr->int1 = eps->inta;
141:     sr->dir = -1;
142:     if (eps->inta <= PETSC_MIN_REAL) { /* Left-open interval */
143:       sr->hasEnd = PETSC_FALSE;
144:       sr->inertia1 = 0;
145:     } else sr->hasEnd = PETSC_TRUE;
146:   }

148:   STGetKSP(eps->st,&ksp);
149:   KSPGetPC(ksp,&pc);
150:   /* compute inertia1 if necessary */
151:   if (sr->hasEnd) {
152:     STSetShift(eps->st,sr->int1);
153:     STSetUp(eps->st);
154:     PCFactorGetMatrix(pc,&F);
155:     MatGetInertia(F,&sr->inertia1,NULL,NULL);
156:   }

158:   /* compute inertia0 */
159:   STSetShift(eps->st,sr->int0);
160:   STSetUp(eps->st);
161:   PCFactorGetMatrix(pc,&F);
162:   MatGetInertia(F,&sr->inertia0,NULL,NULL);

164:   /* number of eigenvalues in interval */
165:   sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
166:   eps->nev = sr->numEigs;
167:   eps->ncv = sr->numEigs;
168:   eps->mpd = sr->numEigs;
169:   EPSSetDimensions_Default(eps,ctx->nev,&ctx->ncv,&ctx->mpd);

171:   /* allocate solution for subsolves */
172:   if (sr->numEigs) {
173:     EPSAllocateSolutionSlice(eps,1);
174:   }
175:   return(0);
176: }

178: /*
179:    Fills the fields of a shift structure
180: */
183: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,EPS_shift neighb0,EPS_shift neighb1)
184: {
185:   PetscErrorCode  ierr;
186:   EPS_shift       s,*pending2;
187:   PetscInt        i;
188:   EPS_SR          sr;
189:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

192:   sr = ctx->sr;
193:   PetscNewLog(eps,&s);
194:   s->value = val;
195:   s->neighb[0] = neighb0;
196:   if (neighb0) neighb0->neighb[1] = s;
197:   s->neighb[1] = neighb1;
198:   if (neighb1) neighb1->neighb[0] = s;
199:   s->comp[0] = PETSC_FALSE;
200:   s->comp[1] = PETSC_FALSE;
201:   s->index = -1;
202:   s->neigs = 0;
203:   s->nconv[0] = s->nconv[1] = 0;
204:   s->nsch[0] = s->nsch[1]=0;
205:   /* Inserts in the stack of pending shifts */
206:   /* If needed, the array is resized */
207:   if (sr->nPend >= sr->maxPend) {
208:     sr->maxPend *= 2;
209:     PetscMalloc1(sr->maxPend,&pending2);
210:     PetscLogObjectMemory((PetscObject)eps,sizeof(EPS_shift));
211:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
212:     PetscFree(sr->pending);
213:     sr->pending = pending2;
214:   }
215:   sr->pending[sr->nPend++]=s;
216:   return(0);
217: }

219: /* Prepare for Rational Krylov update */
222: static PetscErrorCode EPSPrepareRational(EPS eps)
223: {
224:   EPS_KRYLOVSCHUR  *ctx = (EPS_KRYLOVSCHUR*)eps->data;
225:   PetscErrorCode   ierr;
226:   PetscInt         dir,i,k,ld,nv;
227:   PetscScalar      *A;
228:   EPS_SR           sr = ctx->sr;
229:   Vec              v;

232:   DSGetLeadingDimension(eps->ds,&ld);
233:   dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
234:   dir*=sr->dir;
235:   k = 0;
236:   for (i=0;i<sr->nS;i++) {
237:     if (dir*PetscRealPart(sr->S[i])>0.0) {
238:       sr->S[k] = sr->S[i];
239:       sr->S[sr->nS+k] = sr->S[sr->nS+i];
240:       BVGetColumn(sr->Vnext,k,&v);
241:       BVCopyVec(sr->V,eps->nconv+i,v);
242:       BVRestoreColumn(sr->Vnext,k,&v);
243:       k++;
244:       if (k>=sr->nS/2)break;
245:     }
246:   }
247:   /* Copy to DS */
248:   DSGetArray(eps->ds,DS_MAT_A,&A);
249:   PetscMemzero(A,ld*ld*sizeof(PetscScalar));
250:   for (i=0;i<k;i++) {
251:     A[i*(1+ld)] = sr->S[i];
252:     A[k+i*ld] = sr->S[sr->nS+i];
253:   }
254:   sr->nS = k;
255:   DSRestoreArray(eps->ds,DS_MAT_A,&A);
256:   DSGetDimensions(eps->ds,&nv,NULL,NULL,NULL,NULL);
257:   DSSetDimensions(eps->ds,nv,0,0,k);
258:   /* Append u to V */
259:   BVGetColumn(sr->Vnext,sr->nS,&v);
260:   BVCopyVec(sr->V,sr->nv,v);
261:   BVRestoreColumn(sr->Vnext,sr->nS,&v);
262:   return(0);
263: }

265: /* Provides next shift to be computed */
268: static PetscErrorCode EPSExtractShift(EPS eps)
269: {
270:   PetscErrorCode   ierr;
271:   PetscInt         iner;
272:   Mat              F;
273:   PC               pc;
274:   KSP              ksp;
275:   EPS_KRYLOVSCHUR  *ctx = (EPS_KRYLOVSCHUR*)eps->data;
276:   EPS_SR           sr;

279:   sr = ctx->sr;
280:   if (sr->nPend > 0) {
281:     sr->sPrev = sr->sPres;
282:     sr->sPres = sr->pending[--sr->nPend];
283:     STSetShift(eps->st,sr->sPres->value);
284:     STGetKSP(eps->st,&ksp);
285:     KSPGetPC(ksp,&pc);
286:     PCFactorGetMatrix(pc,&F);
287:     MatGetInertia(F,&iner,NULL,NULL);
288:     sr->sPres->inertia = iner;
289:     eps->target = sr->sPres->value;
290:     eps->reason = EPS_CONVERGED_ITERATING;
291:     eps->its = 0;
292:   } else sr->sPres = NULL;
293:   return(0);
294: }

296: /*
297:    Symmetric KrylovSchur adapted to spectrum slicing:
298:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
299:    Returns whether the search has succeeded
300: */
303: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
304: {
305:   PetscErrorCode  ierr;
306:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
307:   PetscInt        i,conv,k,l,ld,nv,*iwork,j,p;
308:   Mat             U;
309:   PetscScalar     *Q,*A,rtmp,*eigrsave,*eigisave;
310:   PetscReal       *a,*b,beta,*errestsave;
311:   PetscBool       breakdown;
312:   PetscInt        count0,count1;
313:   PetscReal       lambda;
314:   EPS_shift       sPres;
315:   PetscBool       complIterating;
316:   PetscBool       sch0,sch1;
317:   PetscInt        iterCompl=0,n0,n1;
318:   EPS_SR          sr = ctx->sr;
319:   BV              bvsave;

322:   bvsave = eps->V;  /* temporarily swap basis vectors */
323:   eps->V = sr->V;
324:   eigrsave = eps->eigr;
325:   eps->eigr = sr->eigr;
326:   eigisave = eps->eigi;
327:   eps->eigi = sr->eigi;
328:   errestsave = eps->errest;
329:   eps->errest = sr->errest;
330:   /* Spectrum slicing data */
331:   sPres = sr->sPres;
332:   complIterating =PETSC_FALSE;
333:   sch1 = sch0 = PETSC_TRUE;
334:   DSGetLeadingDimension(eps->ds,&ld);
335:   PetscMalloc1(2*ld,&iwork);
336:   count0=0;count1=0; /* Found on both sides */
337:   if (sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
338:     /* Rational Krylov */
339:     DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
340:     DSGetDimensions(eps->ds,NULL,NULL,NULL,&l,NULL);
341:     DSSetDimensions(eps->ds,l+1,0,0,0);
342:     BVSetActiveColumns(sr->V,0,l+1);
343:     DSGetMat(eps->ds,DS_MAT_Q,&U);
344:     BVMultInPlace(sr->V,U,0,l+1);
345:     MatDestroy(&U);
346:   } else {
347:     /* Get the starting Lanczos vector */
348:     EPSGetStartVector(eps,0,NULL);
349:     l = 0;
350:   }
351:   /* Restart loop */
352:   while (eps->reason == EPS_CONVERGED_ITERATING) {
353:     eps->its++; sr->itsKs++;
354:     /* Compute an nv-step Lanczos factorization */
355:     nv = PetscMin(eps->nconv+ctx->mpd,ctx->ncv);
356:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
357:     b = a + ld;
358:     EPSFullLanczos(eps,a,b,eps->nconv+l,&nv,&breakdown);
359:     sr->nv = nv;
360:     beta = b[nv-1];
361:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
362:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
363:     if (l==0) {
364:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
365:     } else {
366:       DSSetState(eps->ds,DS_STATE_RAW);
367:     }
368:     BVSetActiveColumns(sr->V,eps->nconv,nv);

370:     /* Solve projected problem and compute residual norm estimates */
371:     if (eps->its == 1 && l > 0) {/* After rational update */
372:       DSGetArray(eps->ds,DS_MAT_A,&A);
373:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
374:       b = a + ld;
375:       k = eps->nconv+l;
376:       A[k*ld+k-1] = A[(k-1)*ld+k];
377:       A[k*ld+k] = a[k];
378:       for (j=k+1; j< nv; j++) {
379:         A[j*ld+j] = a[j];
380:         A[j*ld+j-1] = b[j-1] ;
381:         A[(j-1)*ld+j] = b[j-1];
382:       }
383:       DSRestoreArray(eps->ds,DS_MAT_A,&A);
384:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
385:       DSSolve(eps->ds,sr->eigr,NULL);
386:       DSSort(eps->ds,sr->eigr,NULL,NULL,NULL,NULL);
387:       DSSetCompact(eps->ds,PETSC_TRUE);
388:     } else { /* Restart */
389:       DSSolve(eps->ds,sr->eigr,NULL);
390:       DSSort(eps->ds,sr->eigr,NULL,NULL,NULL,NULL);
391:     }
392:     /* Residual */
393:     EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,beta,1.0,&k);

395:     /* Check convergence */
396:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
397:     b = a + ld;
398:     conv = 0;
399:     j = k = eps->nconv;
400:     for (i=eps->nconv;i<nv;i++) if (sr->errest[i] < eps->tol) conv++;
401:     for (i=eps->nconv;i<nv;i++) {
402:       if (sr->errest[i] < eps->tol) {
403:         iwork[j++]=i;
404:       } else iwork[conv+k++]=i;
405:     }
406:     for (i=eps->nconv;i<nv;i++) {
407:       a[i]=PetscRealPart(sr->eigr[i]);
408:       b[i]=sr->errest[i];
409:     }
410:     for (i=eps->nconv;i<nv;i++) {
411:       sr->eigr[i] = a[iwork[i]];
412:       sr->errest[i] = b[iwork[i]];
413:     }
414:     for (i=eps->nconv;i<nv;i++) {
415:       a[i]=PetscRealPart(sr->eigr[i]);
416:       b[i]=sr->errest[i];
417:     }
418:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
419:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
420:     for (i=eps->nconv;i<nv;i++) {
421:       p=iwork[i];
422:       if (p!=i) {
423:         j=i+1;
424:         while (iwork[j]!=i) j++;
425:         iwork[j]=p;iwork[i]=i;
426:         for (k=0;k<nv;k++) {
427:           rtmp=Q[k+p*ld];Q[k+p*ld]=Q[k+i*ld];Q[k+i*ld]=rtmp;
428:         }
429:       }
430:     }
431:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
432:     k=eps->nconv+conv;

434:     /* Checking values obtained for completing */
435:     for (i=0;i<k;i++) {
436:       sr->back[i]=sr->eigr[i];
437:     }
438:     STBackTransform(eps->st,k,sr->back,sr->eigi);
439:     count0=count1=0;
440:     for (i=0;i<k;i++) {
441:       lambda = PetscRealPart(sr->back[i]);
442:       if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
443:       if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
444:     }
445:     if (k>ctx->nev && ctx->ncv-k<5) eps->reason = EPS_CONVERGED_TOL;
446:     else {
447:       /* Checks completion */
448:       if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
449:         eps->reason = EPS_CONVERGED_TOL;
450:       } else {
451:         if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
452:         if (complIterating) {
453:           if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
454:         } else if (k >= ctx->nev) {
455:           n0 = sPres->nsch[0]-count0;
456:           n1 = sPres->nsch[1]-count1;
457:           if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
458:             /* Iterating for completion*/
459:             complIterating = PETSC_TRUE;
460:             if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
461:             if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
462:             iterCompl = sr->iterCompl;
463:           } else eps->reason = EPS_CONVERGED_TOL;
464:         }
465:       }
466:     }
467:     /* Update l */
468:     if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
469:     else l = nv-k;
470:     if (breakdown) l=0;

472:     if (eps->reason == EPS_CONVERGED_ITERATING) {
473:       if (breakdown) {
474:         /* Start a new Lanczos factorization */
475:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%g)\n",eps->its,(double)beta);
476:         EPSGetStartVector(eps,k,&breakdown);
477:         if (breakdown) {
478:           eps->reason = EPS_DIVERGED_BREAKDOWN;
479:           PetscInfo(eps,"Unable to generate more start vectors\n");
480:         }
481:       } else {
482:         /* Prepare the Rayleigh quotient for restart */
483:         DSGetArrayReal(eps->ds,DS_MAT_T,&a);
484:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
485:         b = a + ld;
486:         for (i=k;i<k+l;i++) {
487:           a[i] = PetscRealPart(sr->eigr[i]);
488:           b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
489:         }
490:         DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
491:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
492:       }
493:     }
494:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
495:     DSGetMat(eps->ds,DS_MAT_Q,&U);
496:     BVMultInPlace(eps->V,U,eps->nconv,k+l);
497:     MatDestroy(&U);

499:     /* Normalize u and append it to V */
500:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
501:       BVCopyColumn(sr->V,nv,k+l);
502:     }
503:     /* Monitor */
504:     eps->nconv = k;
505:     EPSMonitor(eps,ctx->sr->itsKs,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
506:     if (eps->reason != EPS_CONVERGED_ITERATING) {
507:       /* Store approximated values for next shift */
508:       DSGetArray(eps->ds,DS_MAT_Q,&Q);
509:       sr->nS = l;
510:       for (i=0;i<l;i++) {
511:         sr->S[i] = sr->eigr[i+k];/* Diagonal elements */
512:         sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
513:       }
514:       DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
515:     }
516:   }
517:   /* Check for completion */
518:   for (i=0;i< eps->nconv; i++) {
519:     if ((sr->dir)*PetscRealPart(sr->eigr[i])>0) sPres->nconv[1]++;
520:     else sPres->nconv[0]++;
521:   }
522:   sPres->comp[0] = (count0 >= sPres->nsch[0])?PETSC_TRUE:PETSC_FALSE;
523:   sPres->comp[1] = (count1 >= sPres->nsch[1])?PETSC_TRUE:PETSC_FALSE;
524:   if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1])SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
525:   PetscFree(iwork);
526:   eps->V = bvsave;   /* restore basis */
527:   eps->eigr = eigrsave;
528:   eps->eigi = eigisave;
529:   eps->errest = errestsave;
530:   return(0);
531: }

533: /*
534:   Obtains value of subsequent shift
535: */
538: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
539: {
540:   PetscReal       lambda,d_prev;
541:   PetscInt        i,idxP;
542:   EPS_SR          sr;
543:   EPS_shift       sPres,s;
544:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

547:   sr = ctx->sr;
548:   sPres = sr->sPres;
549:   if (sPres->neighb[side]) {
550:   /* Completing a previous interval */
551:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
552:       if (side) *newS = (sPres->value + PetscRealPart(eps->eigr[eps->perm[sr->indexEig-1]]))/2;
553:       else *newS = (sPres->value + PetscRealPart(eps->eigr[eps->perm[0]]))/2;
554:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
555:   } else { /* (Only for side=1). Creating a new interval. */
556:     if (sPres->neigs==0) {/* No value has been accepted*/
557:       if (sPres->neighb[0]) {
558:         /* Multiplying by 10 the previous distance */
559:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
560:         sr->nleap++;
561:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
562:         if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
563:       } else { /* First shift */
564:         if (eps->nconv != 0) {
565:           /* Unaccepted values give information for next shift */
566:           idxP=0;/* Number of values left from shift */
567:           for (i=0;i<eps->nconv;i++) {
568:             lambda = PetscRealPart(sr->eigr[i]);
569:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
570:             else break;
571:           }
572:           /* Avoiding subtraction of eigenvalues (might be the same).*/
573:           if (idxP>0) {
574:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(sr->eigr[0]))/(idxP+0.3);
575:           } else {
576:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(sr->eigr[eps->nconv-1]))/(eps->nconv+0.3);
577:           }
578:           *newS = sPres->value + ((sr->dir)*d_prev*ctx->nev)/2;
579:         } else { /* No values found, no information for next shift */
580:           SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
581:         }
582:       }
583:     } else { /* Accepted values found */
584:       sr->nleap = 0;
585:       /* Average distance of values in previous subinterval */
586:       s = sPres->neighb[0];
587:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
588:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
589:       }
590:       if (s) {
591:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
592:       } else { /* First shift. Average distance obtained with values in this shift */
593:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
594:         if ((sr->dir)*(PetscRealPart(eps->eigr[0])-sPres->value)>0 && PetscAbsReal((PetscRealPart(eps->eigr[sr->indexEig-1]) - PetscRealPart(eps->eigr[0]))/PetscRealPart(eps->eigr[0])) > PetscSqrtReal(eps->tol)) {
595:           d_prev =  PetscAbsReal((PetscRealPart(eps->eigr[sr->indexEig-1]) - PetscRealPart(eps->eigr[0])))/(sPres->neigs+0.3);
596:         } else {
597:           d_prev = PetscAbsReal(PetscRealPart(eps->eigr[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
598:         }
599:       }
600:       /* Average distance is used for next shift by adding it to value on the right or to shift */
601:       if ((sr->dir)*(PetscRealPart(eps->eigr[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
602:         *newS = PetscRealPart(eps->eigr[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(ctx->nev))/2;
603:       } else { /* Last accepted value is on the left of shift. Adding to shift */
604:         *newS = sPres->value + ((sr->dir)*d_prev*(ctx->nev))/2;
605:       }
606:     }
607:     /* End of interval can not be surpassed */
608:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
609:   }/* of neighb[side]==null */
610:   return(0);
611: }

613: /*
614:   Function for sorting an array of real values
615: */
618: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
619: {
620:   PetscReal      re;
621:   PetscInt       i,j,tmp;

624:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
625:   /* Insertion sort */
626:   for (i=1;i<nr;i++) {
627:     re = PetscRealPart(r[perm[i]]);
628:     j = i-1;
629:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
630:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
631:     }
632:   }
633:   return(0);
634: }

636: /* Stores the pairs obtained since the last shift in the global arrays */
639: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
640: {
641:   PetscErrorCode  ierr;
642:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
643:   PetscReal       lambda,err,norm;
644:   PetscInt        i,count;
645:   PetscBool       iscayley;
646:   EPS_SR          sr = ctx->sr;
647:   EPS_shift       sPres;
648:   Vec             v,w;
649:  
651:   sPres = sr->sPres;
652:   sPres->index = sr->indexEig;
653:   count = sr->indexEig;
654:   /* Back-transform */
655:   STBackTransform(eps->st,eps->nconv,sr->eigr,sr->eigi);
656:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
657:   /* Sort eigenvalues */
658:   sortRealEigenvalues(sr->eigr,sr->perm,eps->nconv,PETSC_FALSE,sr->dir);
659:   /* Values stored in global array */
660:   for (i=0;i<eps->nconv;i++) {
661:     lambda = PetscRealPart(sr->eigr[sr->perm[i]]);
662:     err = sr->errest[sr->perm[i]];

664:     if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
665:       if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
666:       eps->eigr[count] = lambda;
667:       eps->errest[count] = err;
668:       /* Explicit purification */
669:       BVGetColumn(eps->V,count,&v);
670:       BVGetColumn(sr->V,sr->perm[i],&w);
671:       STApply(eps->st,w,v);
672:       BVRestoreColumn(eps->V,count,&v);
673:       BVRestoreColumn(sr->V,sr->perm[i],&w);
674:       BVNormColumn(eps->V,count,NORM_2,&norm);
675:       BVScaleColumn(eps->V,count,1.0/norm);
676:       count++;
677:     }
678:   }
679:   sPres->neigs = count - sr->indexEig;
680:   sr->indexEig = count;
681:   /* Global ordering array updating */
682:   sortRealEigenvalues(eps->eigr,eps->perm,count,PETSC_TRUE,sr->dir);
683:   return(0);
684: }

688: static PetscErrorCode EPSLookForDeflation(EPS eps)
689: {
690:   PetscErrorCode  ierr;
691:   PetscReal       val;
692:   PetscInt        i,count0=0,count1=0;
693:   EPS_shift       sPres;
694:   PetscInt        ini,fin,k,idx0,idx1;
695:   EPS_SR          sr;
696:   Vec             v;
697:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

700:   sr = ctx->sr;
701:   sPres = sr->sPres;

703:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
704:   else ini = 0;
705:   fin = sr->indexEig;
706:   /* Selection of ends for searching new values */
707:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
708:   else sPres->ext[0] = sPres->neighb[0]->value;
709:   if (!sPres->neighb[1]) {
710:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
711:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
712:   } else sPres->ext[1] = sPres->neighb[1]->value;
713:   /* Selection of values between right and left ends */
714:   for (i=ini;i<fin;i++) {
715:     val=PetscRealPart(eps->eigr[eps->perm[i]]);
716:     /* Values to the right of left shift */
717:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
718:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
719:       else count1++;
720:     } else break;
721:   }
722:   /* The number of values on each side are found */
723:   if (sPres->neighb[0]) {
724:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
725:     if (sPres->nsch[0]<0)SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
726:   } else sPres->nsch[0] = 0;

728:   if (sPres->neighb[1]) {
729:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
730:     if (sPres->nsch[1]<0)SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
731:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

733:   /* Completing vector of indexes for deflation */
734:   idx0 = ini;
735:   idx1 = ini+count0+count1;
736:   k=0;
737:   for (i=idx0;i<idx1;i++) sr->idxDef[k++]=eps->perm[i];
738:   BVDuplicateResize(sr->V,k+ctx->ncv+1,&sr->Vnext);
739:   BVSetNumConstraints(sr->Vnext,k);
740:   for (i=0;i<k;i++) {
741:     BVGetColumn(sr->Vnext,-i-1,&v);
742:     BVCopyVec(eps->V,sr->idxDef[i],v);
743:     BVRestoreColumn(sr->Vnext,-i-1,&v);
744:   }

746:   /* For rational Krylov */
747:   if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
748:     EPSPrepareRational(eps);
749:   }
750:   eps->nconv = 0;
751:   /* Get rid of temporary Vnext */
752:   BVDestroy(&sr->V);
753:   sr->V = sr->Vnext;
754:   sr->Vnext = NULL;
755:   return(0);
756: }

760: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
761: {
762:   PetscErrorCode  ierr;
763:   PetscInt        i,lds;
764:   PetscReal       newS;
765:   EPS_SR          sr;
766:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

769:   sr = ctx->sr;
770:   /* Only with eigenvalues present in the interval ...*/
771:   if (sr->numEigs==0) {
772:     eps->reason = EPS_CONVERGED_TOL;
773:     return(0);
774:   }
775:   /* Array of pending shifts */
776:   sr->maxPend = 100; /* Initial size */
777:   sr->nPend = 0;
778:   PetscMalloc1(sr->maxPend,&sr->pending);
779:   PetscLogObjectMemory((PetscObject)eps,(sr->maxPend)*sizeof(EPS_shift));
780:   EPSCreateShift(eps,sr->int0,NULL,NULL);
781:   /* extract first shift */
782:   sr->sPrev = NULL;
783:   sr->sPres = sr->pending[--sr->nPend];
784:   sr->sPres->inertia = sr->inertia0;
785:   eps->target = sr->sPres->value;
786:   sr->s0 = sr->sPres;
787:   sr->indexEig = 0;
788:   /* Memory reservation for auxiliary variables */
789:   lds = PetscMin(ctx->mpd,ctx->ncv);
790:   PetscCalloc1(lds*lds,&sr->S);
791:   PetscMalloc1(ctx->ncv,&sr->back);
792:   PetscLogObjectMemory((PetscObject)eps,(sr->numEigs+2*ctx->ncv)*sizeof(PetscScalar));
793:   for (i=0;i<ctx->ncv;i++) {
794:     sr->eigr[i]    = 0.0;
795:     sr->eigi[i]   = 0.0;
796:     sr->errest[i] = 0.0;
797:   }
798:   for (i=0;i<sr->numEigs;i++) eps->perm[i] = i;
799:   /* Vectors for deflation */
800:   PetscMalloc1(sr->numEigs,&sr->idxDef);
801:   PetscLogObjectMemory((PetscObject)eps,sr->numEigs*sizeof(PetscInt));
802:   sr->indexEig = 0;
803:   /* Main loop */
804:   while (sr->sPres) {
805:     /* Search for deflation */
806:     EPSLookForDeflation(eps);
807:     /* KrylovSchur */
808:     EPSKrylovSchur_Slice(eps);

810:     EPSStoreEigenpairs(eps);
811:     /* Select new shift */
812:     if (!sr->sPres->comp[1]) {
813:       EPSGetNewShiftValue(eps,1,&newS);
814:       EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
815:     }
816:     if (!sr->sPres->comp[0]) {
817:       /* Completing earlier interval */
818:       EPSGetNewShiftValue(eps,0,&newS);
819:       EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
820:     }
821:     /* Preparing for a new search of values */
822:     EPSExtractShift(eps);
823:   }

825:   /* Updating eps values prior to exit */
826:   BVDestroy(&sr->V);
827:   PetscFree4(sr->eigr,sr->eigi,sr->errest,sr->perm);
828:   PetscFree(sr->S);
829:   PetscFree(sr->idxDef);
830:   PetscFree(sr->pending);
831:   PetscFree(sr->back);
832:   eps->nconv  = sr->indexEig;
833:   eps->reason = EPS_CONVERGED_TOL;
834:   eps->its    = sr->itsKs;
835:   eps->nds    = 0;
836:   return(0);
837: }