Actual source code: ks-slice.c

  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 solving
 10:            sparse symmetric generalized eigenproblems", SIAM J. Matrix Analysis
 11:            and App., 15(1), pp. 228–272, 1994.

 13:        [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
 14:            SIAM J. Matrix Analysis and App., 23(3), pp. 601-614, 2001. 

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

 20:    This file is part of SLEPc.

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

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

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

 36: #include <private/epsimpl.h>                /*I "slepceps.h" I*/
 37: #include <slepcblaslapack.h>



 42: /* Type of data characterizing a shift (place from where an eps is applied) */
 43: typedef struct _n_shift *shift;
 44: struct _n_shift{
 45:   PetscReal        value;
 46:   PetscInt        inertia;
 47:   PetscBool        comp[2]; /* Shows completion of subintervals (left and right) */
 48:   shift          neighb[2];/* Adjacent shifts */
 49:   PetscInt        index;/* Index in eig where found values are stored */
 50:   PetscInt        neigs; /* Number of values found */
 51:   PetscReal     ext[2];   /* Limits for accepted values */
 52:   PetscInt      nsch[2];  /* Number of missing values for each subinterval */
 53:   PetscInt      nconv[2]; /* Converged on each side (accepted or not)*/
 54: };

 56: /* Type of data  for storing the state of spectrum slicing*/
 57: struct _n_SR{
 58:   PetscReal       int0,int1; /* Extremes of the interval */
 59:   PetscInt        dir; /* Determines the order of values in eig (+1 incr, -1 decr) */
 60:   PetscBool       hasEnd; /* Tells whether the interval has an end */
 61:   PetscInt        inertia0,inertia1;
 62:   Vec             *V;
 63:   PetscScalar     *eig,*eigi,*monit,*back;
 64:   PetscReal       *errest;
 65:   PetscInt        *perm;/* Permutation for keeping the eigenvalues in order */
 66:   PetscInt        numEigs; /* Number of eigenvalues in the interval */
 67:   PetscInt        indexEig;
 68:   shift           sPres; /* Present shift */
 69:   shift           *pending;/* Pending shifts array */
 70:   PetscInt        nPend;/* Number of pending shifts */
 71:   PetscInt        maxPend;/* Size of "pending" array */
 72:   Vec             *VDef; /* Vector for deflation */
 73:   PetscInt        *idxDef;/* For deflation */
 74:   PetscInt        nMAXCompl;
 75:   PetscInt        iterCompl;
 76:   PetscInt        itsKs; /* Krylovschur restarts */
 77:   PetscInt        nleap;
 78:   shift           s0;/* Initial shift */
 79: };
 80: typedef struct _n_SR  *SR;

 82: /* 
 83:    Fills the fields of a shift structure

 85: */
 88: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val, shift neighb0,shift neighb1)
 89: {
 90:   PetscErrorCode   ierr;
 91:   shift            s,*pending2;
 92:   PetscInt         i;
 93:   SR               sr;

 96:   sr = (SR)(eps->data);
 97:   PetscMalloc(sizeof(struct _n_shift),&s);
 98:   s->value = val;
 99:   s->neighb[0] = neighb0;
100:   if(neighb0) neighb0->neighb[1] = s;
101:   s->neighb[1] = neighb1;
102:   if(neighb1) neighb1->neighb[0] = s;
103:   s->comp[0] = PETSC_FALSE;
104:   s->comp[1] = PETSC_FALSE;
105:   s->index = -1;
106:   s->neigs = 0;
107:   s->nconv[0] = s->nconv[1] = 0;
108:   s->nsch[0] = s->nsch[1]=0;
109:   /* Inserts in the stack of pending shifts */
110:   /* If needed, the array is resized */
111:   if(sr->nPend >= sr->maxPend){
112:     sr->maxPend *= 2;
113:     PetscMalloc((sr->maxPend)*sizeof(shift),&pending2);
114:     for(i=0;i < sr->nPend; i++)pending2[i] = sr->pending[i];
115:     PetscFree(sr->pending);
116:     sr->pending = pending2;
117:   }
118:   sr->pending[sr->nPend++]=s;
119:   return(0);
120: }

122: /* Provides next shift to be computed */
125: static PetscErrorCode EPSExtractShift(EPS eps){
126:   PetscErrorCode   ierr;
127:   PetscInt         iner;
128:   Mat              F;
129:   PC               pc;
130:   KSP              ksp;
131:   SR               sr;

134:   sr = (SR)(eps->data);
135:   if(sr->nPend > 0){
136:     sr->sPres = sr->pending[--sr->nPend];
137:     STSetShift(eps->OP, sr->sPres->value);
138:     STGetKSP(eps->OP, &ksp);
139:     KSPGetPC(ksp, &pc);
140:     PCFactorGetMatrix(pc,&F);
141:     MatGetInertia(F,&iner,PETSC_NULL,PETSC_NULL);
142:     sr->sPres->inertia = iner;
143:     eps->target = sr->sPres->value;
144:     eps->nconv = 0;
145:     eps->reason = EPS_CONVERGED_ITERATING;
146:     eps->its = 0;
147:    }else sr->sPres = PETSC_NULL;
148:    return(0);
149: }
150: 
151: /*
152:    Symmetric KrylovSchur adapted to spectrum slicing:
153:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
154:    Returns whether the search has succeeded
155: */
158: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
159: {
161:   PetscInt       i,conv,k,l,lds,lt,nv,m,*iwork,p,j;
162:   Vec            u=eps->work[0];
163:   PetscScalar    *Q,nu;
164:   PetscReal      *a,*b,*work,beta,rtmp;
165:   PetscBool      breakdown;
166:   PetscInt       count0,count1;
167:   PetscReal      theta,lambda;
168:   shift          sPres;
169:   PetscBool      complIterating,iscayley;/* Shows whether iterations are made for completion */
170:   PetscBool      sch0,sch1;/* Shows whether values are looked after on each side */
171:   PetscInt       iterCompl=0,n0,n1,aux,auxc;
172:   SR             sr;

175:   /* Spectrum slicing data */
176:   sr = (SR)eps->data;
177:   sPres = sr->sPres;
178:   complIterating =PETSC_FALSE;
179:   sch1 = sch0 = PETSC_TRUE;
180:   lds = PetscMin(eps->mpd,eps->ncv);
181:   PetscMalloc(lds*lds*sizeof(PetscReal),&work);
182:   PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);
183:   PetscMalloc(2*lds*sizeof(PetscInt),&iwork);
184:   lt = PetscMin(eps->nev+eps->mpd,eps->ncv);
185:   PetscMalloc(lt*sizeof(PetscReal),&a);
186:   PetscMalloc(lt*sizeof(PetscReal),&b);
187:   count0=0;count1=0; /* Found on both sides */

189:    /* filling in values for the monitor */
190:   PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&iscayley);
191:   if(iscayley){
192:     STCayleyGetAntishift(eps->OP,&nu);
193:     for(i=0;i<sr->indexEig;i++){
194:       sr->monit[i]=(nu + sr->eig[i])/(sr->eig[i] - sPres->value);
195:     }
196:   }else{
197:     for(i=0;i<sr->indexEig;i++){
198:       sr->monit[i]=1.0/(sr->eig[i] - sPres->value);
199:     }
200:   }
201: 
202: 
203:   /* Get the starting Lanczos vector */
204:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
205:   l = 0;

207:   /* Restart loop */
208:   while (eps->reason == EPS_CONVERGED_ITERATING) {
209:     eps->its++; sr->itsKs++;
210:     /* Compute an nv-step Lanczos factorization */
211:     m = PetscMin(eps->nconv+eps->mpd,eps->ncv);

213:     EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);
214:     nv = m - eps->nconv;
215:     beta = b[nv-1];

217:     /* Solve projected problem and compute residual norm estimates */
218:     EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);
219:     /* Residual */
220:     EPSKrylovConvergence(eps,PETSC_TRUE,PETSC_TRUE,eps->nconv,nv,PETSC_NULL,nv,Q,eps->V+eps->nconv,nv,beta,1.0,&k,PETSC_NULL);
221:     /* Check convergence */
222:     conv=k=j=0;
223:     for(i=0;i<nv;i++)if(eps->errest[eps->nconv+i] < eps->tol)conv++;
224:     for(i=0;i<nv;i++){
225:       if(eps->errest[eps->nconv+i] < eps->tol){
226:         iwork[j++]=i;
227:       }else iwork[conv+k++]=i;
228:     }
229:     for(i=0;i<nv;i++){
230:       a[i]=PetscRealPart(eps->eigr[eps->nconv+i]);
231:       b[i]=eps->errest[eps->nconv+i];
232:     }
233:     for(i=0;i<nv;i++){
234:       eps->eigr[eps->nconv+i] = a[iwork[i]];
235:       eps->errest[eps->nconv+i] = b[iwork[i]];
236:     }
237:     for( i=0;i<nv;i++){
238:       p=iwork[i];
239:       if(p!=i){
240:         j=i+1;
241:         while(iwork[j]!=i)j++;
242:         iwork[j]=p;iwork[i]=i;
243:         for(k=0;k<nv;k++){
244:           rtmp=Q[k+p*nv];Q[k+p*nv]=Q[k+i*nv];Q[k+i*nv]=rtmp;
245:         }
246:       }
247:     }
248:     k=eps->nconv+conv;
249:     /* Checking values obtained for completing */
250:     for(i=0;i<k;i++){
251:       sr->back[i]=eps->eigr[i];
252:     }
253:     STBackTransform(eps->OP,k,sr->back,eps->eigi);
254:     count0=count1=0;
255:     for(i=0;i<k;i++){
256:       theta = PetscRealPart(eps->eigr[i]);
257:       lambda = sr->back[i];
258:       if( ((sr->dir)*theta < 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0))count0++;
259:       if( ((sr->dir)*theta > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0))count1++;
260:     }
261: 
262:     /* Checks completion */
263:     if( (!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1]) ) {
264:       eps->reason = EPS_CONVERGED_TOL;
265:     }else {
266:       if(!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
267:       if(complIterating){
268:         if(--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
269:       }else if (k >= eps->nev) {
270:         n0 = sPres->nsch[0]-count0;
271:         n1 = sPres->nsch[1]-count1;
272:         if( sr->iterCompl>0 && ( (n0>0&& n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl) )){
273:           /* Iterating for completion*/
274:           complIterating = PETSC_TRUE;
275:           if(n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
276:           if(n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
277:           iterCompl = sr->iterCompl;
278:         }else eps->reason = EPS_CONVERGED_TOL;
279:       }
280:     }
281: 
282:     /* Update l */
283:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
284:     else l = (eps->nconv+nv-k)/2;

286:     if (eps->reason == EPS_CONVERGED_ITERATING) {
287:       if (breakdown) {
288:         /* Start a new Lanczos factorization */
289:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
290:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
291:         if (breakdown) {
292:           eps->reason = EPS_DIVERGED_BREAKDOWN;
293:           PetscInfo(eps,"Unable to generate more start vectors\n");
294:         }
295:       } else {
296:         /* Prepare the Rayleigh quotient for restart */
297:         for (i=0;i<l;i++) {
298:           a[i] = PetscRealPart(eps->eigr[i+k]);
299:           b[i] = PetscRealPart(Q[nv-1+(i+k-eps->nconv)*nv]*beta);
300:         }
301:       }
302:     }
303:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
304:     SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);
305:     /* Normalize u and append it to V */
306:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
307:       VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
308:     }
309:     if(eps->numbermonitors >0){
310:       aux = auxc = 0;
311:       for(i=0;i<nv+eps->nconv;i++){
312:         sr->back[i]=eps->eigr[i];
313:       }
314:       STBackTransform(eps->OP,nv+eps->nconv,sr->back,eps->eigi);
315:       for(i=0;i<nv+eps->nconv;i++){
316:         lambda = sr->back[i];
317:         if( ((sr->dir)*(lambda - sPres->ext[0]) > 0)&& ((sr->dir)*(sPres->ext[1] - lambda) > 0)){
318:           sr->monit[sr->indexEig+aux]=eps->eigr[i];
319:           sr->errest[sr->indexEig+aux]=eps->errest[i];
320:           aux++;
321:           if(eps->errest[i] < eps->tol)auxc++;
322:         }
323:       }
324:       EPSMonitor(eps,eps->its,auxc+sr->indexEig,sr->monit,sr->eigi,sr->errest,sr->indexEig+aux);
325:     }
326:     eps->nconv = k;
327:   }
328:   /* Check for completion */
329:   for(i=0;i< eps->nconv; i++){
330:     if( (sr->dir)*PetscRealPart(eps->eigr[i])>0 )sPres->nconv[1]++;
331:     else sPres->nconv[0]++;
332:   }
333:   sPres->comp[0] = (count0 >= sPres->nsch[0])?PETSC_TRUE:PETSC_FALSE;
334:   sPres->comp[1] = (count1 >= sPres->nsch[1])?PETSC_TRUE:PETSC_FALSE;
335:   PetscFree(Q);
336:   PetscFree(a);
337:   PetscFree(b);
338:   PetscFree(work);
339:   PetscFree(iwork);
340:   return(0);
341: }

343: /* 
344:   Obtains value of subsequent shift
345: */
348: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS){
349:   PetscReal   lambda,d_prev;
350:   PetscInt    i,idxP;
351:   SR          sr;
352:   shift       sPres;

355:   sr = (SR)eps->data;
356:   sPres = sr->sPres;
357:   if( sPres->neighb[side]){
358:   /* Completing a previous interval */
359:     if(!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0){ /* One of the ends might be too far from eigenvalues */
360:       if(side) *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[sr->indexEig-1]]))/2;
361:       else *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[0]]))/2;
362:     }else *newS=(sPres->value + sPres->neighb[side]->value)/2;
363:   }else{ /* (Only for side=1). Creating a new interval. */
364:     if(sPres->neigs==0){/* No value has been accepted*/
365:       if(sPres->neighb[0]){
366:         /* Multiplying by 10 the previous distance */
367:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
368:         sr->nleap++;
369:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
370:         if( !sr->hasEnd && sr->nleap > 5)SETERRQ(((PetscObject)eps)->comm,1,"Unable to compute the wanted eigenvalues with open interval");
371:       }else {/* First shift */
372:         if(eps->nconv != 0){
373:            /* Unaccepted values give information for next shift */
374:            idxP=0;/* Number of values left from shift */
375:            for(i=0;i<eps->nconv;i++){
376:              lambda = PetscRealPart(eps->eigr[i]);
377:              if( (sr->dir)*(lambda - sPres->value) <0)idxP++;
378:              else break;
379:            }
380:            /* Avoiding subtraction of eigenvalues (might be the same).*/
381:            if(idxP>0){
382:              d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
383:            }else {
384:              d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
385:            }
386:            *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
387:         }else{/* No values found, no information for next shift */
388:           SETERRQ(((PetscObject)eps)->comm,1,"First shift renders no information");
389:         }
390:       }
391:     }else{/* Accepted values found */
392:       sr->nleap = 0;
393:       /* Average distance of values in previous subinterval */
394:       shift s = sPres->neighb[0];
395:       while(s && PetscAbs(s->inertia - sPres->inertia)==0){
396:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
397:       }
398:       if(s){
399:         d_prev = PetscAbsReal( (sPres->value - s->value)/(sPres->inertia - s->inertia));
400:       }else{/* First shift. Average distance obtained with values in this shift */
401:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
402:         if( (sr->dir)*(PetscRealPart(sr->eig[0])-sPres->value)>0 && PetscAbsReal( (PetscRealPart(sr->eig[sr->indexEig-1]) - PetscRealPart(sr->eig[0]))/PetscRealPart(sr->eig[0])) > PetscSqrtReal(eps->tol) ){
403:           d_prev =  PetscAbsReal( (PetscRealPart(sr->eig[sr->indexEig-1]) - PetscRealPart(sr->eig[0])))/(sPres->neigs+0.3);
404:         }else{
405:           d_prev = PetscAbsReal( PetscRealPart(sr->eig[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
406:         }
407:       }
408:       /* Average distance is used for next shift by adding it to value on the right or to shift */
409:       if( (sr->dir)*(PetscRealPart(sr->eig[sPres->index + sPres->neigs -1]) - sPres->value) >0){
410:         *newS = PetscRealPart(sr->eig[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
411:       }else{/* Last accepted value is on the left of shift. Adding to shift */
412:         *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
413:       }
414:     }
415:     /* End of interval can not be surpassed */
416:     if((sr->dir)*( sr->int1 - *newS) < 0) *newS = sr->int1;
417:   }/* of neighb[side]==null */
418:   return(0);
419: }

421: /* 
422:   Function for sorting an array of real values
423: */
426: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
427: {
428:   PetscReal      re;
429:   PetscInt       i,j,tmp;
430: 
432:   if(!prev) for (i=0; i < nr; i++) { perm[i] = i; }
433:   /* Insertion sort */
434:   for (i=1; i < nr; i++) {
435:     re = PetscRealPart(r[perm[i]]);
436:     j = i-1;
437:     while ( j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0 ) {
438:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
439:     }
440:   }
441:   return(0);
442: }

444: /* Stores the pairs obtained since the last shift in the global arrays */
447: PetscErrorCode EPSStoreEigenpairs(EPS eps)
448: {
450:   PetscReal      lambda,err,norm;
451:   PetscInt       i,count;
452:   PetscBool      cond,iscayley;
453:   SR             sr;
454:   shift          sPres;

457:   sr = (SR)(eps->data);
458:   sPres = sr->sPres;
459:   sPres->index = sr->indexEig;
460:   count = sr->indexEig;
461:   /* Backtransform */
462:   EPSBackTransform_Default(eps);
463:   PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&iscayley);
464:   /* Sort eigenvalues */
465:   sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
466:   /* Values stored in global array */
467:   /* Condition for avoiding comparing with a non-existing end */
468:   cond = (!sPres->neighb[1] && !sr->hasEnd)?PETSC_TRUE:PETSC_FALSE;
469:   for( i=0; i < eps->nconv ;i++ ){
470:     lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
471:     err = eps->errest[eps->perm[i]];
472:     if(  (sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0  ){/* Valid value */
473:       if(count>=sr->numEigs){/* Error found */
474:          SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error in Spectrum Slicing!");
475:       }
476:       sr->eig[count] = lambda;
477:       sr->errest[count] = err;
478:       /* Purification */
479:       if (eps->isgeneralized && !iscayley){
480:         STApply(eps->OP,eps->V[eps->perm[i]],sr->V[count]);
481:         IPNorm(eps->ip,sr->V[count],&norm);
482:         VecScale(sr->V[count],1.0/norm);
483:       }else{
484:         VecCopy(eps->V[eps->perm[i]], sr->V[count]);
485:       }
486:       count++;
487:     }
488:   }
489:   sPres->neigs = count - sr->indexEig;
490:   sr->indexEig = count;
491:   /* Global ordering array updating */
492:   sortRealEigenvalues(sr->eig,sr->perm,count,PETSC_TRUE,sr->dir);
493:   return(0);
494: }

498: PetscErrorCode EPSLookForDeflation(EPS eps)
499: {
500:   PetscReal       val;
501:   PetscInt        i,count0=0,count1=0;
502:   shift           sPres;
503:   PetscInt        ini,fin,k,idx0,idx1;
504:   SR              sr;

507:   sr = (SR)(eps->data);
508:   sPres = sr->sPres;

510:   if(sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
511:   else ini = 0;
512:   fin = sr->indexEig;
513:   /* Selection of ends for searching new values */
514:   if(!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
515:   else sPres->ext[0] = sPres->neighb[0]->value;
516:   if(!sPres->neighb[1]) {
517:     if(sr->hasEnd) sPres->ext[1] = sr->int1;
518:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
519:   }else sPres->ext[1] = sPres->neighb[1]->value;
520:   /* Selection of values between right and left ends */
521:   for(i=ini;i<fin;i++){
522:     val=PetscRealPart(sr->eig[sr->perm[i]]);
523:     /* Values to the right of left shift */
524:     if( (sr->dir)*(val - sPres->ext[1]) < 0 ){
525:       if((sr->dir)*(val - sPres->value) < 0)count0++;
526:       else count1++;
527:     }else break;
528:   }
529:   /* The number of values on each side are found */
530:   if(sPres->neighb[0]){
531:      sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
532:      if(sPres->nsch[0]<0)SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
533:   }else sPres->nsch[0] = 0;

535:   if(sPres->neighb[1]){
536:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
537:     if(sPres->nsch[1]<0)SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
538:   }else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
539: 
540:   /* Completing vector of indexes for deflation */
541:   idx0 = ini;
542:   idx1 = ini+count0+count1;
543:   k=0;
544:   for(i=idx0;i<idx1;i++)sr->idxDef[k++]=sr->perm[i];
545:   for(i=0;i<k;i++)sr->VDef[i]=sr->V[sr->idxDef[i]];
546:   eps->DS = sr->VDef;
547:   eps->nds = k;
548:   return(0);
549: }

553: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
554: {
556:   PetscInt       i;
557:   PetscReal      newS;
558:   KSP            ksp;
559:   PC             pc;
560:   Mat            F;
561:   PetscReal      *errest_left;
562:   Vec            t;
563:   SR             sr;
564: 
566: #if defined(PETSC_USE_COMPLEX)
567:   SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Spectrum slicing not supported in complex scalars");
568: #endif
569:   PetscMalloc(sizeof(struct _n_SR),&sr);
570:   eps->data = sr;
571:   sr->itsKs = 0;
572:   sr->nleap = 0;
573:   sr->nMAXCompl = eps->nev/4;
574:   sr->iterCompl = eps->max_it/4;
575:   /* Checking presence of ends and finding direction */
576:   if( eps->inta > PETSC_MIN_REAL){
577:     sr->int0 = eps->inta;
578:     sr->int1 = eps->intb;
579:     sr->dir = 1;
580:     if(eps->intb >= PETSC_MAX_REAL){ /* Right-open interval */
581:       sr->hasEnd = PETSC_FALSE;
582:       sr->inertia1 = eps->n;
583:     }else sr->hasEnd = PETSC_TRUE;
584:   }else{ /* Left-open interval */
585:     sr->int0 = eps->intb;
586:     sr->int1 = eps->inta;
587:     sr->dir = -1;
588:     sr->hasEnd = PETSC_FALSE;
589:     sr->inertia1 = 0;
590:   }
591:   /* Array of pending shifts */
592:   sr->maxPend = 100;/* Initial size */
593:   PetscMalloc((sr->maxPend)*sizeof(shift),&sr->pending);
594:   if(sr->hasEnd){
595:     STGetKSP(eps->OP, &ksp);
596:     KSPGetPC(ksp, &pc);
597:     PCFactorGetMatrix(pc,&F);
598:     /* Not looking for values in b (just inertia).*/
599:     MatGetInertia(F,&sr->inertia1,PETSC_NULL,PETSC_NULL);
600:     PCReset(pc); /* avoiding memory leak */
601:   }
602:   sr->nPend = 0;
603:   EPSCreateShift(eps,sr->int0,PETSC_NULL,PETSC_NULL);
604:   EPSExtractShift(eps);
605:   sr->s0 = sr->sPres;
606:   sr->inertia0 = sr->s0->inertia;
607:   sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
608:   sr->indexEig = 0;
609:   /* Only with eigenvalues present in the interval ...*/
610:   if(sr->numEigs==0){
611:     eps->reason = EPS_CONVERGED_TOL;
612:     PetscFree(sr->s0);
613:     PetscFree(sr->pending);
614:     PetscFree(sr);
615:     return(0);
616:   }
617:   /* Memory reservation for eig, V and perm */
618:   PetscMalloc((sr->numEigs)*sizeof(PetscScalar),&sr->eig);
619:   PetscMalloc((sr->numEigs)*sizeof(PetscScalar),&sr->eigi);
620:   PetscMalloc((sr->numEigs+eps->ncv) *sizeof(PetscReal),&sr->errest);
621:   PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscReal),&errest_left);
622:   PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscScalar),&sr->monit);
623:   PetscMalloc((eps->ncv)*sizeof(PetscScalar),&sr->back);
624:   for(i=0;i<sr->numEigs;i++){sr->eigi[i]=0;sr->eig[i] = 0;}
625:   for(i=0;i<sr->numEigs+eps->ncv;i++){errest_left[i]=0;sr->errest[i]=0;sr->monit[i]=0;}
626:   VecCreateMPI(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,&t);
627:   VecDuplicateVecs(t,sr->numEigs,&sr->V);
628:   VecDestroy(&t);
629:   /* Vector for maintaining order of eigenvalues */
630:   PetscMalloc((sr->numEigs)*sizeof(PetscInt),&sr->perm);
631:   for(i=0;i< sr->numEigs;i++)sr->perm[i]=i;
632:   /* Vectors for deflation */
633:   PetscMalloc((sr->numEigs)*sizeof(PetscInt),&sr->idxDef);
634:   PetscMalloc((sr->numEigs)*sizeof(Vec),&sr->VDef);
635:   sr->indexEig = 0;

637:   while(sr->sPres){
638:     /* Search for deflation */
639:     EPSLookForDeflation(eps);
640:     /* KrylovSchur */
641:     EPSKrylovSchur_Slice(eps);
642: 
643:     EPSStoreEigenpairs(eps);
644:     /* Select new shift */
645:     if(!sr->sPres->comp[1]){
646:       EPSGetNewShiftValue(eps,1,&newS);
647:       EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
648:     }
649:     if(!sr->sPres->comp[0]){
650:       /* Completing earlier interval */
651:       EPSGetNewShiftValue(eps,0,&newS);
652:       EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
653:     }
654:     /* Preparing for a new search of values */
655:     EPSExtractShift(eps);
656:   }

658:   /* Updating eps values prior to exit */
659: 
660:   VecDestroyVecs(eps->allocated_ncv,&eps->V);
661:   eps->V = sr->V;
662:   PetscFree(eps->eigr);
663:   PetscFree(eps->eigi);
664:   PetscFree(eps->errest);
665:   PetscFree(eps->errest_left);
666:   PetscFree(eps->perm);
667:   eps->eigr = sr->eig;
668:   eps->eigi = sr->eigi;
669:   eps->errest = sr->errest;
670:   eps->errest_left = errest_left;
671:   eps->perm = sr->perm;
672:   eps->ncv = eps->allocated_ncv = sr->numEigs;
673:   eps->nconv = sr->indexEig;
674:   eps->reason = EPS_CONVERGED_TOL;
675:   eps->its = sr->itsKs;
676:   eps->nds = 0;
677:   eps->DS = PETSC_NULL;
678:   eps->evecsavailable = PETSC_TRUE;
679:   PetscFree(sr->VDef);
680:   PetscFree(sr->idxDef);
681:   PetscFree(sr->pending);
682:   PetscFree(sr->monit);
683:   PetscFree(sr->back);
684:   /* Reviewing list of shifts to free memory */
685:   shift s = sr->s0;
686:   if(s){
687:     while(s->neighb[1]){
688:       s = s->neighb[1];
689:       PetscFree(s->neighb[0]);
690:     }
691:     PetscFree(s);
692:   }
693:   PetscFree(sr);
694:   return(0);
695: }