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

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

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

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

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

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

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

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

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

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

556:   sr = ctx->sr;
557:   sPres = sr->sPres;

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

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

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

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

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

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

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