Actual source code: ks-slice.c

slepc-3.5.1 2014-09-01
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 (!((PetscObject)(eps->V))->type_name) {
 73:     BVSetType(sr->V,BVSVEC);
 74:   } else {
 75:     BVGetType(eps->V,&type);
 76:     BVSetType(sr->V,type);
 77:   }
 78:   STMatGetVecs(eps->st,&t,NULL);
 79:   BVSetSizesFromVec(sr->V,t,requested);
 80:   VecDestroy(&t);
 81:   EPS_SetInnerProduct(eps);
 82:   BVGetMatrix(eps->V,&matrix,NULL);
 83:   BVSetMatrix(sr->V,matrix,PETSC_FALSE);
 84:   BVGetOrthogonalization(eps->V,&orthog_type,&orthog_ref,&eta);
 85:   BVSetOrthogonalization(sr->V,orthog_type,orthog_ref,eta);
 86:   return(0);
 87: }

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

102:   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");
103:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
104:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
105:   if (!((PetscObject)(eps->st))->type_name) { /* default to shift-and-invert */
106:     STSetType(eps->st,STSINVERT);
107:   }
108:   PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
109:   if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
110:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL*1e-2;  /* use tighter tolerance */
111:   if (!eps->max_it) eps->max_it = 100;
112:   if (ctx->nev==1) ctx->nev = 40;  /* nev not set, use default value */
113:   if (ctx->nev<10) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
114:   eps->ops->backtransform = NULL;

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

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

146:   if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
147:     if (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");
148:     STSetShift(eps->st,eps->inta);
149:   } else {
150:     STSetShift(eps->st,eps->intb);
151:   }
152:   STSetUp(eps->st);
153:   STGetKSP(eps->st,&ksp);
154:   KSPGetPC(ksp,&pc);

156:   /* compute inertia1 if necessary */
157:   if (sr->hasEnd) {
158:     PCFactorGetMatrix(pc,&F);
159:     MatGetInertia(F,&sr->inertia1,NULL,NULL);
160:   }

162:   /* compute inertia0 */
163:   STSetShift(eps->st,sr->int0);
164:   PCFactorGetMatrix(pc,&F);
165:   MatGetInertia(F,&sr->inertia0,NULL,NULL);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

616: /*
617:   Function for sorting an array of real values
618: */
621: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
622: {
623:   PetscReal      re;
624:   PetscInt       i,j,tmp;

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

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

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

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

703:   sr = ctx->sr;
704:   sPres = sr->sPres;

706:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
707:   else ini = 0;
708:   fin = sr->indexEig;
709:   /* Selection of ends for searching new values */
710:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
711:   else sPres->ext[0] = sPres->neighb[0]->value;
712:   if (!sPres->neighb[1]) {
713:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
714:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
715:   } else sPres->ext[1] = sPres->neighb[1]->value;
716:   /* Selection of values between right and left ends */
717:   for (i=ini;i<fin;i++) {
718:     val=PetscRealPart(eps->eigr[eps->perm[i]]);
719:     /* Values to the right of left shift */
720:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
721:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
722:       else count1++;
723:     } else break;
724:   }
725:   /* The number of values on each side are found */
726:   if (sPres->neighb[0]) {
727:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
728:     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");
729:   } else sPres->nsch[0] = 0;

731:   if (sPres->neighb[1]) {
732:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
733:     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");
734:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

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

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

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

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

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

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