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
 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-2013, 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>                /*I "slepceps.h" I*/
 38: #include <slepcblaslapack.h>
 39:  #include krylovschur.h

 41: /*
 42:    Fills the fields of a shift structure

 44: */
 47: static PetscErrorCode EPSCreateShift(EPS eps,PetscReal val,shift neighb0,shift neighb1)
 48: {
 49:   PetscErrorCode   ierr;
 50:   shift            s,*pending2;
 51:   PetscInt         i;
 52:   SR               sr;
 53:   EPS_KRYLOVSCHUR  *ctx = (EPS_KRYLOVSCHUR*)eps->data;

 56:   sr = ctx->sr;
 57:   PetscMalloc(sizeof(struct _n_shift),&s);
 58:   PetscLogObjectMemory(eps,sizeof(struct _n_shift));
 59:   s->value = val;
 60:   s->neighb[0] = neighb0;
 61:   if (neighb0) neighb0->neighb[1] = s;
 62:   s->neighb[1] = neighb1;
 63:   if (neighb1) neighb1->neighb[0] = s;
 64:   s->comp[0] = PETSC_FALSE;
 65:   s->comp[1] = PETSC_FALSE;
 66:   s->index = -1;
 67:   s->neigs = 0;
 68:   s->nconv[0] = s->nconv[1] = 0;
 69:   s->nsch[0] = s->nsch[1]=0;
 70:   /* Inserts in the stack of pending shifts */
 71:   /* If needed, the array is resized */
 72:   if (sr->nPend >= sr->maxPend) {
 73:     sr->maxPend *= 2;
 74:     PetscMalloc((sr->maxPend)*sizeof(shift),&pending2);
 75:     PetscLogObjectMemory(eps,sizeof(shift));
 76:     for (i=0;i<sr->nPend;i++) pending2[i] = sr->pending[i];
 77:     PetscFree(sr->pending);
 78:     sr->pending = pending2;
 79:   }
 80:   sr->pending[sr->nPend++]=s;
 81:   return(0);
 82: }

 84: /* Provides next shift to be computed */
 87: static PetscErrorCode EPSExtractShift(EPS eps)
 88: {
 89:   PetscErrorCode   ierr;
 90:   PetscInt         iner,dir,i,k,ld;
 91:   PetscScalar      *A;
 92:   Mat              F;
 93:   PC               pc;
 94:   KSP              ksp;
 95:   EPS_KRYLOVSCHUR  *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 96:   SR               sr;

 99:   sr = ctx->sr;
100:   if (sr->nPend > 0) {
101:     sr->sPrev = sr->sPres;
102:     sr->sPres = sr->pending[--sr->nPend];
103:     STSetShift(eps->st,sr->sPres->value);
104:     STGetKSP(eps->st,&ksp);
105:     KSPGetPC(ksp,&pc);
106:     PCFactorGetMatrix(pc,&F);
107:     MatGetInertia(F,&iner,NULL,NULL);
108:     sr->sPres->inertia = iner;
109:     eps->target = sr->sPres->value;
110:     eps->reason = EPS_CONVERGED_ITERATING;
111:     eps->its = 0;
112:     /* For rational Krylov */
113:     if (sr->nS>0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1])) {
114:       DSGetLeadingDimension(eps->ds,&ld);
115:       dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
116:       dir*=sr->dir;
117:       k = 0;
118:       for (i=0;i<sr->nS;i++) {
119:         if (dir*PetscRealPart(sr->S[i])>0.0) {
120:           sr->S[k] = sr->S[i];
121:           sr->S[sr->nS+k] = sr->S[sr->nS+i];
122:           VecCopy(eps->V[eps->nconv+i],eps->V[k++]);
123:           if (k>=sr->nS/2)break;
124:         }
125:       }
126:       /* Copy to DS */
127:       DSGetArray(eps->ds,DS_MAT_A,&A);
128:       PetscMemzero(A,ld*ld*sizeof(PetscScalar));
129:       for (i=0;i<k;i++) {
130:         A[i*(1+ld)] = sr->S[i];
131:         A[k+i*ld] = sr->S[sr->nS+i];
132:       }
133:       sr->nS = k;
134:       DSRestoreArray(eps->ds,DS_MAT_A,&A);
135:       DSSetDimensions(eps->ds,0,0,0,k);
136:       /* Normalize u and append it to V */
137:       VecAXPBY(eps->V[sr->nS],1.0/sr->beta,0.0,eps->work[0]);
138:     }
139:     eps->nconv = 0;
140:   } else sr->sPres = NULL;
141:   return(0);
142: }

144: /*
145:    Symmetric KrylovSchur adapted to spectrum slicing:
146:    Allows searching an specific amount of eigenvalues in the subintervals left and right.
147:    Returns whether the search has succeeded
148: */
151: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
152: {
153:   PetscErrorCode  ierr;
154:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
155:   PetscInt        i,conv,k,l,ld,nv,*iwork,j,p;
156:   Vec             u=eps->work[0];
157:   PetscScalar     *Q,*A,nu,rtmp;
158:   PetscReal       *a,*b,beta;
159:   PetscBool       breakdown;
160:   PetscInt        count0,count1;
161:   PetscReal       lambda;
162:   shift           sPres;
163:   PetscBool       complIterating,iscayley;
164:   PetscBool       sch0,sch1;
165:   PetscInt        iterCompl=0,n0,n1,aux,auxc;
166:   SR              sr;

169:   /* Spectrum slicing data */
170:   sr = ctx->sr;
171:   sPres = sr->sPres;
172:   complIterating =PETSC_FALSE;
173:   sch1 = sch0 = PETSC_TRUE;
174:   DSGetLeadingDimension(eps->ds,&ld);
175:   PetscMalloc(2*ld*sizeof(PetscInt),&iwork);
176:   count0=0;count1=0; /* Found on both sides */
177:   /* filling in values for the monitor */
178:   if (eps->numbermonitors >0) {
179:     PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
180:     if (iscayley) {
181:       STCayleyGetAntishift(eps->st,&nu);
182:       for (i=0;i<sr->indexEig;i++) {
183:         sr->monit[i]=(nu + sr->eig[i])/(sr->eig[i] - sPres->value);
184:       }
185:     } else {
186:       for (i=0;i<sr->indexEig;i++) {
187:         sr->monit[i]=1.0/(sr->eig[i] - sPres->value);
188:       }
189:     }
190:   }
191:   if (sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev)) {
192:     /* Rational Krylov */
193:     DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
194:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
195:     DSGetDimensions(eps->ds,NULL,NULL,NULL,&l,NULL);
196:     SlepcUpdateVectors(l+1,eps->V,0,l+1,Q,ld,PETSC_FALSE);
197:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
198:   } else {
199:     /* Get the starting Lanczos vector */
200:     EPSGetStartVector(eps,0,eps->V[0],NULL);
201:     l = 0;
202:   }
203:   /* Restart loop */
204:   while (eps->reason == EPS_CONVERGED_ITERATING) {
205:     eps->its++; sr->itsKs++;
206:     /* Compute an nv-step Lanczos factorization */
207:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
208:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
209:     b = a + ld;
210:     EPSFullLanczos(eps,a,b,eps->V,eps->nconv+l,&nv,u,&breakdown);
211:     beta = b[nv-1];
212:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
213:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
214:     if (l==0) {
215:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
216:     } else {
217:       DSSetState(eps->ds,DS_STATE_RAW);
218:     }

220:     /* Solve projected problem and compute residual norm estimates */
221:     if (eps->its == 1 && l > 0) {/* After rational update */
222:       DSGetArray(eps->ds,DS_MAT_A,&A);
223:       DSGetArrayReal(eps->ds,DS_MAT_T,&a);
224:       b = a + ld;
225:       k = eps->nconv+l;
226:       A[k*ld+k-1] = A[(k-1)*ld+k];
227:       A[k*ld+k] = a[k];
228:       for (j=k+1; j< nv; j++) {
229:         A[j*ld+j] = a[j];
230:         A[j*ld+j-1] = b[j-1] ;
231:         A[(j-1)*ld+j] = b[j-1];
232:       }
233:       DSRestoreArray(eps->ds,DS_MAT_A,&A);
234:       DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
235:       DSSolve(eps->ds,eps->eigr,NULL);
236:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
237:       DSSetCompact(eps->ds,PETSC_TRUE);
238:     } else { /* Restart */
239:       DSSolve(eps->ds,eps->eigr,NULL);
240:       DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
241:     }
242:     /* Residual */
243:     EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,1.0,&k);

245:     /* Check convergence */
246:     DSGetArrayReal(eps->ds,DS_MAT_T,&a);
247:     b = a + ld;
248:     conv = 0;
249:     j = k = eps->nconv;
250:     for (i=eps->nconv;i<nv;i++) if (eps->errest[i] < eps->tol) conv++;
251:     for (i=eps->nconv;i<nv;i++) {
252:       if (eps->errest[i] < eps->tol) {
253:         iwork[j++]=i;
254:       } else iwork[conv+k++]=i;
255:     }
256:     for (i=eps->nconv;i<nv;i++) {
257:       a[i]=PetscRealPart(eps->eigr[i]);
258:       b[i]=eps->errest[i];
259:     }
260:     for (i=eps->nconv;i<nv;i++) {
261:       eps->eigr[i] = a[iwork[i]];
262:       eps->errest[i] = b[iwork[i]];
263:     }
264:     for (i=eps->nconv;i<nv;i++) {
265:       a[i]=PetscRealPart(eps->eigr[i]);
266:       b[i]=eps->errest[i];
267:     }
268:     DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
269:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
270:     for (i=eps->nconv;i<nv;i++) {
271:       p=iwork[i];
272:       if (p!=i) {
273:         j=i+1;
274:         while (iwork[j]!=i) j++;
275:         iwork[j]=p;iwork[i]=i;
276:         for (k=0;k<nv;k++) {
277:           rtmp=Q[k+p*ld];Q[k+p*ld]=Q[k+i*ld];Q[k+i*ld]=rtmp;
278:         }
279:       }
280:     }
281:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
282:     k=eps->nconv+conv;

284:     /* Checking values obtained for completing */
285:     for (i=0;i<k;i++) {
286:       sr->back[i]=eps->eigr[i];
287:     }
288:     STBackTransform(eps->st,k,sr->back,eps->eigi);
289:     count0=count1=0;
290:     for (i=0;i<k;i++) {
291:       lambda = PetscRealPart(sr->back[i]);
292:       if (((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0)) count0++;
293:       if (((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0)) count1++;
294:     }

296:     /* Checks completion */
297:     if ((!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1])) {
298:       eps->reason = EPS_CONVERGED_TOL;
299:     } else {
300:       if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
301:       if (complIterating) {
302:         if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
303:       } else if (k >= eps->nev) {
304:         n0 = sPres->nsch[0]-count0;
305:         n1 = sPres->nsch[1]-count1;
306:         if (sr->iterCompl>0 && ((n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl))) {
307:           /* Iterating for completion*/
308:           complIterating = PETSC_TRUE;
309:           if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
310:           if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
311:           iterCompl = sr->iterCompl;
312:         } else eps->reason = EPS_CONVERGED_TOL;
313:       }
314:     }
315:     /* Update l */
316:     if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
317:     else l = nv-k;
318:     if (breakdown) l=0;

320:     if (eps->reason == EPS_CONVERGED_ITERATING) {
321:       if (breakdown) {
322:         /* Start a new Lanczos factorization */
323:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
324:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
325:         if (breakdown) {
326:           eps->reason = EPS_DIVERGED_BREAKDOWN;
327:           PetscInfo(eps,"Unable to generate more start vectors\n");
328:         }
329:       } else {
330:         /* Prepare the Rayleigh quotient for restart */
331:         DSGetArrayReal(eps->ds,DS_MAT_T,&a);
332:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
333:         b = a + ld;
334:         for (i=k;i<k+l;i++) {
335:           a[i] = PetscRealPart(eps->eigr[i]);
336:           b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
337:         }
338:         DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
339:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
340:       }
341:     }
342:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
343:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
344:     SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
345:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
346:     /* Normalize u and append it to V */
347:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
348:       VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
349:     }
350:     /* Monitor */
351:     if (eps->numbermonitors >0) {
352:       aux = auxc = 0;
353:       for (i=0;i<nv;i++) {
354:         sr->back[i]=eps->eigr[i];
355:       }
356:       STBackTransform(eps->st,nv,sr->back,eps->eigi);
357:       for (i=0;i<nv;i++) {
358:         lambda = PetscRealPart(sr->back[i]);
359:         if (((sr->dir)*(lambda - sPres->ext[0]) > 0)&& ((sr->dir)*(sPres->ext[1] - lambda) > 0)) {
360:           sr->monit[sr->indexEig+aux]=eps->eigr[i];
361:           sr->errest[sr->indexEig+aux]=eps->errest[i];
362:           aux++;
363:           if (eps->errest[i] < eps->tol)auxc++;
364:         }
365:       }
366:       EPSMonitor(eps,eps->its,auxc+sr->indexEig,sr->monit,sr->eigi,sr->errest,sr->indexEig+aux);
367:     }
368:     eps->nconv = k;
369:     if (eps->reason != EPS_CONVERGED_ITERATING) {
370:       /* Store approximated values for next shift */
371:       DSGetArray(eps->ds,DS_MAT_Q,&Q);
372:       sr->nS = l;
373:       for (i=0;i<l;i++) {
374:         sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
375:         sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
376:       }
377:       sr->beta = beta;
378:       DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
379:     }
380:   }
381:   /* Check for completion */
382:   for (i=0;i< eps->nconv; i++) {
383:     if ((sr->dir)*PetscRealPart(eps->eigr[i])>0) sPres->nconv[1]++;
384:     else sPres->nconv[0]++;
385:   }
386:   sPres->comp[0] = (count0 >= sPres->nsch[0])?PETSC_TRUE:PETSC_FALSE;
387:   sPres->comp[1] = (count1 >= sPres->nsch[1])?PETSC_TRUE:PETSC_FALSE;
388:   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");
389:   PetscFree(iwork);
390:   return(0);
391: }

393: /*
394:   Obtains value of subsequent shift
395: */
398: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
399: {
400:   PetscReal       lambda,d_prev;
401:   PetscInt        i,idxP;
402:   SR              sr;
403:   shift           sPres,s;
404:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

407:   sr = ctx->sr;
408:   sPres = sr->sPres;
409:   if (sPres->neighb[side]) {
410:   /* Completing a previous interval */
411:     if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
412:       if (side) *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[sr->indexEig-1]]))/2;
413:       else *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[0]]))/2;
414:     } else *newS=(sPres->value + sPres->neighb[side]->value)/2;
415:   } else { /* (Only for side=1). Creating a new interval. */
416:     if (sPres->neigs==0) {/* No value has been accepted*/
417:       if (sPres->neighb[0]) {
418:         /* Multiplying by 10 the previous distance */
419:         *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
420:         sr->nleap++;
421:         /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
422:         if (!sr->hasEnd && sr->nleap > 5) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to compute the wanted eigenvalues with open interval");
423:       } else { /* First shift */
424:         if (eps->nconv != 0) {
425:           /* Unaccepted values give information for next shift */
426:           idxP=0;/* Number of values left from shift */
427:           for (i=0;i<eps->nconv;i++) {
428:             lambda = PetscRealPart(eps->eigr[i]);
429:             if ((sr->dir)*(lambda - sPres->value) <0) idxP++;
430:             else break;
431:           }
432:           /* Avoiding subtraction of eigenvalues (might be the same).*/
433:           if (idxP>0) {
434:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
435:           } else {
436:             d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
437:           }
438:           *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
439:         } else { /* No values found, no information for next shift */
440:           SETERRQ(PetscObjectComm((PetscObject)eps),1,"First shift renders no information");
441:         }
442:       }
443:     } else { /* Accepted values found */
444:       sr->nleap = 0;
445:       /* Average distance of values in previous subinterval */
446:       s = sPres->neighb[0];
447:       while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
448:         s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
449:       }
450:       if (s) {
451:         d_prev = PetscAbsReal((sPres->value - s->value)/(sPres->inertia - s->inertia));
452:       } else { /* First shift. Average distance obtained with values in this shift */
453:         /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
454:         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)) {
455:           d_prev =  PetscAbsReal((PetscRealPart(sr->eig[sr->indexEig-1]) - PetscRealPart(sr->eig[0])))/(sPres->neigs+0.3);
456:         } else {
457:           d_prev = PetscAbsReal(PetscRealPart(sr->eig[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
458:         }
459:       }
460:       /* Average distance is used for next shift by adding it to value on the right or to shift */
461:       if ((sr->dir)*(PetscRealPart(sr->eig[sPres->index + sPres->neigs -1]) - sPres->value)>0) {
462:         *newS = PetscRealPart(sr->eig[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
463:       } else { /* Last accepted value is on the left of shift. Adding to shift */
464:         *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
465:       }
466:     }
467:     /* End of interval can not be surpassed */
468:     if ((sr->dir)*(sr->int1 - *newS) < 0) *newS = sr->int1;
469:   }/* of neighb[side]==null */
470:   return(0);
471: }

473: /*
474:   Function for sorting an array of real values
475: */
478: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
479: {
480:   PetscReal      re;
481:   PetscInt       i,j,tmp;

484:   if (!prev) for (i=0;i<nr;i++) perm[i] = i;
485:   /* Insertion sort */
486:   for (i=1;i<nr;i++) {
487:     re = PetscRealPart(r[perm[i]]);
488:     j = i-1;
489:     while (j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0) {
490:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
491:     }
492:   }
493:   return(0);
494: }

496: /* Stores the pairs obtained since the last shift in the global arrays */
499: static PetscErrorCode EPSStoreEigenpairs(EPS eps)
500: {
501:   PetscErrorCode  ierr;
502:   PetscReal       lambda,err,norm;
503:   PetscInt        i,count;
504:   PetscBool       iscayley;
505:   SR              sr;
506:   shift           sPres;
507:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

510:   sr = ctx->sr;
511:   sPres = sr->sPres;
512:   sPres->index = sr->indexEig;
513:   count = sr->indexEig;
514:   /* Back-transform */
515:   EPSBackTransform_Default(eps);
516:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
517:   /* Sort eigenvalues */
518:   sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
519:   /* Values stored in global array */
520:   for (i=0;i<eps->nconv;i++) {
521:     lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
522:     err = eps->errest[eps->perm[i]];

524:     if ((sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0) {/* Valid value */
525:       if (count>=sr->numEigs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error in Spectrum Slicing");
526:       sr->eig[count] = lambda;
527:       sr->errest[count] = err;
528:       /* Explicit purification */
529:       STApply(eps->st,eps->V[eps->perm[i]],sr->V[count]);
530:       IPNorm(eps->ip,sr->V[count],&norm);
531:       VecScale(sr->V[count],1.0/norm);
532:       count++;
533:     }
534:   }
535:   sPres->neigs = count - sr->indexEig;
536:   sr->indexEig = count;
537:   /* Global ordering array updating */
538:   sortRealEigenvalues(sr->eig,sr->perm,count,PETSC_TRUE,sr->dir);
539:   return(0);
540: }

544: static PetscErrorCode EPSLookForDeflation(EPS eps)
545: {
546:   PetscReal       val;
547:   PetscInt        i,count0=0,count1=0;
548:   shift           sPres;
549:   PetscInt        ini,fin,k,idx0,idx1;
550:   SR              sr;
551:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

554:   sr = ctx->sr;
555:   sPres = sr->sPres;

557:   if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
558:   else ini = 0;
559:   fin = sr->indexEig;
560:   /* Selection of ends for searching new values */
561:   if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
562:   else sPres->ext[0] = sPres->neighb[0]->value;
563:   if (!sPres->neighb[1]) {
564:     if (sr->hasEnd) sPres->ext[1] = sr->int1;
565:     else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
566:   } else sPres->ext[1] = sPres->neighb[1]->value;
567:   /* Selection of values between right and left ends */
568:   for (i=ini;i<fin;i++) {
569:     val=PetscRealPart(sr->eig[sr->perm[i]]);
570:     /* Values to the right of left shift */
571:     if ((sr->dir)*(val - sPres->ext[1]) < 0) {
572:       if ((sr->dir)*(val - sPres->value) < 0) count0++;
573:       else count1++;
574:     } else break;
575:   }
576:   /* The number of values on each side are found */
577:   if (sPres->neighb[0]) {
578:     sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
579:     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");
580:   } else sPres->nsch[0] = 0;

582:   if (sPres->neighb[1]) {
583:     sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
584:     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");
585:   } else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);

587:   /* Completing vector of indexes for deflation */
588:   idx0 = ini;
589:   idx1 = ini+count0+count1;
590:   k=0;
591:   for (i=idx0;i<idx1;i++) sr->idxDef[k++]=sr->perm[i];
592:   for (i=0;i<k;i++) sr->VDef[i]=sr->V[sr->idxDef[i]];
593:   eps->defl = sr->VDef;
594:   eps->nds = k;
595:   return(0);
596: }

600: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
601: {
602:   PetscErrorCode  ierr;
603:   PetscInt        i,lds;
604:   PetscReal       newS;
605:   KSP             ksp;
606:   PC              pc;
607:   Mat             F;
608:   PetscReal       *errest_left;
609:   Vec             t;
610:   SR              sr;
611:   shift           s;
612:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

615:   PetscMalloc(sizeof(struct _n_SR),&sr);
616:   PetscLogObjectMemory(eps,sizeof(struct _n_SR));
617:   ctx->sr = sr;
618:   sr->itsKs = 0;
619:   sr->nleap = 0;
620:   sr->nMAXCompl = eps->nev/4;
621:   sr->iterCompl = eps->max_it/4;
622:   sr->sPres = NULL;
623:   sr->nS = 0;
624:   lds = PetscMin(eps->mpd,eps->ncv);
625:   /* Checking presence of ends and finding direction */
626:   if (eps->inta > PETSC_MIN_REAL) {
627:     sr->int0 = eps->inta;
628:     sr->int1 = eps->intb;
629:     sr->dir = 1;
630:     if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
631:       sr->hasEnd = PETSC_FALSE;
632:       sr->inertia1 = eps->n;
633:     } else sr->hasEnd = PETSC_TRUE;
634:   } else { /* Left-open interval */
635:     sr->int0 = eps->intb;
636:     sr->int1 = eps->inta;
637:     sr->dir = -1;
638:     sr->hasEnd = PETSC_FALSE;
639:     sr->inertia1 = 0;
640:   }
641:   /* Array of pending shifts */
642:   sr->maxPend = 100;/* Initial size */
643:   PetscMalloc((sr->maxPend)*sizeof(shift),&sr->pending);
644:   PetscLogObjectMemory(eps,(sr->maxPend)*sizeof(shift));
645:   if (sr->hasEnd) {
646:     STGetKSP(eps->st,&ksp);
647:     KSPGetPC(ksp,&pc);
648:     PCFactorGetMatrix(pc,&F);
649:     /* Not looking for values in b (just inertia).*/
650:     MatGetInertia(F,&sr->inertia1,NULL,NULL);
651:     PCReset(pc); /* avoiding memory leak */
652:   }
653:   sr->nPend = 0;
654:   EPSCreateShift(eps,sr->int0,NULL,NULL);
655:   EPSExtractShift(eps);
656:   sr->s0 = sr->sPres;
657:   sr->inertia0 = sr->s0->inertia;
658:   sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
659:   sr->indexEig = 0;
660:   /* Only with eigenvalues present in the interval ...*/
661:   if (sr->numEigs==0) {
662:     eps->reason = EPS_CONVERGED_TOL;
663:     PetscFree(sr->s0);
664:     PetscFree(sr->pending);
665:     PetscFree(sr);
666:     return(0);
667:   }
668:   /* Memory reservation for eig, V and perm */
669:   PetscMalloc(lds*lds*sizeof(PetscScalar),&sr->S);
670:   PetscMemzero(sr->S,lds*lds*sizeof(PetscScalar));
671:   PetscMalloc(sr->numEigs*sizeof(PetscScalar),&sr->eig);
672:   PetscMalloc(sr->numEigs*sizeof(PetscScalar),&sr->eigi);
673:   PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscReal),&sr->errest);
674:   PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscReal),&errest_left);
675:   PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscScalar),&sr->monit);
676:   PetscMalloc((eps->ncv)*sizeof(PetscScalar),&sr->back);
677:   PetscLogObjectMemory(eps,(lds*lds+3*sr->numEigs+eps->ncv)*sizeof(PetscScalar)+2*(sr->numEigs+eps->ncv)*sizeof(PetscReal));
678:   for (i=0;i<sr->numEigs;i++) {
679:     sr->eig[i]  = 0.0;
680:     sr->eigi[i] = 0.0;
681:   }
682:   for (i=0;i<sr->numEigs+eps->ncv;i++) {
683:     errest_left[i] = 0.0;
684:     sr->errest[i]  = 0.0;
685:     sr->monit[i]   = 0.0;
686:   }
687:   VecCreateMPI(PetscObjectComm((PetscObject)eps),eps->nloc,PETSC_DECIDE,&t);
688:   VecDuplicateVecs(t,sr->numEigs,&sr->V);
689:   PetscLogObjectParents(eps,sr->numEigs,sr->V);
690:   VecDestroy(&t);
691:   /* Vector for maintaining order of eigenvalues */
692:   PetscMalloc(sr->numEigs*sizeof(PetscInt),&sr->perm);
693:   PetscLogObjectMemory(eps,sr->numEigs*sizeof(PetscInt));
694:   for (i=0;i< sr->numEigs;i++) sr->perm[i]=i;
695:   /* Vectors for deflation */
696:   PetscMalloc(sr->numEigs*sizeof(PetscInt),&sr->idxDef);
697:   PetscLogObjectMemory(eps,sr->numEigs*sizeof(PetscInt));
698:   PetscMalloc(sr->numEigs*sizeof(Vec),&sr->VDef);
699:   PetscLogObjectMemory(eps,sr->numEigs*sizeof(Vec));
700:   sr->indexEig = 0;
701:   /* Main loop */
702:   while (sr->sPres) {
703:     /* Search for deflation */
704:     EPSLookForDeflation(eps);
705:     /* KrylovSchur */
706:     EPSKrylovSchur_Slice(eps);

708:     EPSStoreEigenpairs(eps);
709:     /* Select new shift */
710:     if (!sr->sPres->comp[1]) {
711:       EPSGetNewShiftValue(eps,1,&newS);
712:       EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
713:     }
714:     if (!sr->sPres->comp[0]) {
715:       /* Completing earlier interval */
716:       EPSGetNewShiftValue(eps,0,&newS);
717:       EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
718:     }
719:     /* Preparing for a new search of values */
720:     EPSExtractShift(eps);
721:   }

723:   /* Updating eps values prior to exit */
724:   VecDestroyVecs(eps->allocated_ncv,&eps->V);
725:   eps->V = sr->V;
726:   PetscFree(sr->S);
727:   PetscFree(eps->eigr);
728:   PetscFree(eps->eigi);
729:   PetscFree(eps->errest);
730:   PetscFree(eps->errest_left);
731:   PetscFree(eps->perm);
732:   eps->eigr = sr->eig;
733:   eps->eigi = sr->eigi;
734:   eps->errest = sr->errest;
735:   eps->errest_left = errest_left;
736:   eps->perm = sr->perm;
737:   eps->ncv = eps->allocated_ncv = sr->numEigs;
738:   eps->nconv = sr->indexEig;
739:   eps->reason = EPS_CONVERGED_TOL;
740:   eps->its = sr->itsKs;
741:   eps->nds = 0;
742:   eps->defl = NULL;
743:   eps->evecsavailable = PETSC_TRUE;
744:   PetscFree(sr->VDef);
745:   PetscFree(sr->idxDef);
746:   PetscFree(sr->pending);
747:   PetscFree(sr->monit);
748:   PetscFree(sr->back);
749:   /* Reviewing list of shifts to free memory */
750:   s = sr->s0;
751:   if (s) {
752:     while (s->neighb[1]) {
753:       s = s->neighb[1];
754:       PetscFree(s->neighb[0]);
755:     }
756:     PetscFree(s);
757:   }
758:   PetscFree(sr);
759:   return(0);
760: }