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] C. Campos and J. E. Roman, "Spectrum slicing strategies based on
14: restarted Lanczos methods", Numer. Algor., 2012.
15: DOI: 10.1007/s11075-012-9564-z
17: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
18: SLEPc - Scalable Library for Eigenvalue Problem Computations
19: Copyright (c) 2002-2012, 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: s->value = val;
59: s->neighb[0] = neighb0;
60: if (neighb0) neighb0->neighb[1] = s;
61: s->neighb[1] = neighb1;
62: if (neighb1) neighb1->neighb[0] = s;
63: s->comp[0] = PETSC_FALSE;
64: s->comp[1] = PETSC_FALSE;
65: s->index = -1;
66: s->neigs = 0;
67: s->nconv[0] = s->nconv[1] = 0;
68: s->nsch[0] = s->nsch[1]=0;
69: /* Inserts in the stack of pending shifts */
70: /* If needed, the array is resized */
71: if (sr->nPend >= sr->maxPend) {
72: sr->maxPend *= 2;
73: PetscMalloc((sr->maxPend)*sizeof(shift),&pending2);
74: for (i=0;i < sr->nPend; i++)pending2[i] = sr->pending[i];
75: PetscFree(sr->pending);
76: sr->pending = pending2;
77: }
78: sr->pending[sr->nPend++]=s;
79: return(0);
80: }
82: /* Provides next shift to be computed */
85: static PetscErrorCode EPSExtractShift(EPS eps) {
86: PetscErrorCode ierr;
87: PetscInt iner,dir,i,k,ld;
88: PetscScalar *A;
89: Mat F;
90: PC pc;
91: KSP ksp;
92: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
93: SR sr;
96: sr = ctx->sr;
97: if (sr->nPend > 0) {
98: sr->sPrev = sr->sPres;
99: sr->sPres = sr->pending[--sr->nPend];
100: STSetShift(eps->OP, sr->sPres->value);
101: STGetKSP(eps->OP, &ksp);
102: KSPGetPC(ksp, &pc);
103: PCFactorGetMatrix(pc,&F);
104: MatGetInertia(F,&iner,PETSC_NULL,PETSC_NULL);
105: sr->sPres->inertia = iner;
106: eps->target = sr->sPres->value;
107: eps->reason = EPS_CONVERGED_ITERATING;
108: eps->its = 0;
109: /* For rational Krylov */
110: if (sr->nS >0 && (sr->sPrev == sr->sPres->neighb[0] || sr->sPrev == sr->sPres->neighb[1]) ) {
111: DSGetLeadingDimension(eps->ds,&ld);
112: dir = (sr->sPres->neighb[0] == sr->sPrev)?1:-1;
113: dir*=sr->dir;
114: k = 0;
115: for (i=0;i<sr->nS;i++) {
116: if (dir*PetscRealPart(sr->S[i])>0.0) {
117: sr->S[k] = sr->S[i];
118: sr->S[sr->nS+k] = sr->S[sr->nS+i];
119: VecCopy(eps->V[eps->nconv+i],eps->V[k++]);
120: if (k>=sr->nS/2)break;
121: }
122: }
123: /* Copy to DS */
124: DSGetArray(eps->ds,DS_MAT_A,&A);
125: PetscMemzero(A,ld*ld*sizeof(PetscScalar));
126: for (i=0;i<k;i++) {
127: A[i*(1+ld)] = sr->S[i];
128: A[k+i*ld] = sr->S[sr->nS+i];
129: }
130: sr->nS = k;
131: DSRestoreArray(eps->ds,DS_MAT_A,&A);
132: DSSetDimensions(eps->ds,PETSC_IGNORE,PETSC_IGNORE,0,k);
133: /* Normalize u and append it to V */
134: VecAXPBY(eps->V[sr->nS],1.0/sr->beta,0.0,eps->work[0]);
135: }
136: eps->nconv = 0;
137: }else sr->sPres = PETSC_NULL;
138: return(0);
139: }
141: /*
142: Symmetric KrylovSchur adapted to spectrum slicing:
143: Allows searching an specific amount of eigenvalues in the subintervals left and right.
144: Returns whether the search has succeeded
145: */
148: static PetscErrorCode EPSKrylovSchur_Slice(EPS eps)
149: {
150: PetscErrorCode ierr;
151: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
152: PetscInt i,conv,k,l,ld,nv,*iwork,j,p;
153: Vec u=eps->work[0];
154: PetscScalar *Q,*A,nu,rtmp;
155: PetscReal *a,*b,beta;
156: PetscBool breakdown;
157: PetscInt count0,count1;
158: PetscReal lambda;
159: shift sPres;
160: PetscBool complIterating,iscayley;
161: PetscBool sch0,sch1;
162: PetscInt iterCompl=0,n0,n1,aux,auxc;
163: SR sr;
166: /* Spectrum slicing data */
167: sr = ctx->sr;
168: sPres = sr->sPres;
169: complIterating =PETSC_FALSE;
170: sch1 = sch0 = PETSC_TRUE;
171: DSGetLeadingDimension(eps->ds,&ld);
172: PetscMalloc(2*ld*sizeof(PetscInt),&iwork);
173: count0=0;count1=0; /* Found on both sides */
174: /* filling in values for the monitor */
175: if (eps->numbermonitors >0) {
176: PetscObjectTypeCompare((PetscObject)eps->OP,STCAYLEY,&iscayley);
177: if (iscayley) {
178: STCayleyGetAntishift(eps->OP,&nu);
179: for (i=0;i<sr->indexEig;i++) {
180: sr->monit[i]=(nu + sr->eig[i])/(sr->eig[i] - sPres->value);
181: }
182: } else {
183: for (i=0;i<sr->indexEig;i++) {
184: sr->monit[i]=1.0/(sr->eig[i] - sPres->value);
185: }
186: }
187: }
188: if (sr->nS > 0 && (sPres->neighb[0] == sr->sPrev || sPres->neighb[1] == sr->sPrev) ) {
189: /* Rational Krylov */
190: DSTranslateRKS(eps->ds,sr->sPrev->value-sPres->value);
191: DSGetArray(eps->ds,DS_MAT_Q,&Q);
192: DSGetDimensions(eps->ds,PETSC_NULL,PETSC_NULL,PETSC_NULL,&l);
193: SlepcUpdateVectors(l+1,eps->V,0,l+1,Q,ld,PETSC_FALSE);
194: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
195: } else {
196: /* Get the starting Lanczos vector */
197: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
198: l = 0;
199: }
200: /* Restart loop */
201: while (eps->reason == EPS_CONVERGED_ITERATING) {
202: eps->its++; sr->itsKs++;
203: /* Compute an nv-step Lanczos factorization */
204: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
205: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
206: b = a + ld;
207: EPSFullLanczos(eps,a,b,eps->V,eps->nconv+l,&nv,u,&breakdown);
208: beta = b[nv-1];
209: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
210: DSSetDimensions(eps->ds,nv,PETSC_IGNORE,eps->nconv,eps->nconv+l);
211: if (l==0) {
212: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
213: } else {
214: DSSetState(eps->ds,DS_STATE_RAW);
215: }
217: /* Solve projected problem and compute residual norm estimates */
218: if (eps->its == 1 && l > 0) {/* After rational update */
219: DSGetArray(eps->ds,DS_MAT_A,&A);
220: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
221: b = a + ld;
222: k = eps->nconv+l;
223: A[k*ld+k-1] = A[(k-1)*ld+k];
224: A[k*ld+k] = a[k];
225: for (j=k+1; j< nv; j++) {
226: A[j*ld+j] = a[j];
227: A[j*ld+j-1] = b[j-1] ;
228: A[(j-1)*ld+j] = b[j-1];
229: }
230: DSRestoreArray(eps->ds,DS_MAT_A,&A);
231: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
232: DSSolve(eps->ds,eps->eigr,PETSC_NULL);
233: DSSort(eps->ds,eps->eigr,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
234: DSSetCompact(eps->ds,PETSC_TRUE);
235: } else {/* Restart */
236: DSSolve(eps->ds,eps->eigr,PETSC_NULL);
237: DSSort(eps->ds,eps->eigr,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
238: }
239: /* Residual */
240: EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,1.0,&k);
242: /* Check convergence */
243: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
244: b = a + ld;
245: conv = 0;
246: j = k = eps->nconv;
247: for (i=eps->nconv;i<nv;i++)if (eps->errest[i] < eps->tol)conv++;
248: for (i=eps->nconv;i<nv;i++) {
249: if (eps->errest[i] < eps->tol) {
250: iwork[j++]=i;
251: }else iwork[conv+k++]=i;
252: }
253: for (i=eps->nconv;i<nv;i++) {
254: a[i]=PetscRealPart(eps->eigr[i]);
255: b[i]=eps->errest[i];
256: }
257: for (i=eps->nconv;i<nv;i++) {
258: eps->eigr[i] = a[iwork[i]];
259: eps->errest[i] = b[iwork[i]];
260: }
261: for (i=eps->nconv;i<nv;i++) {
262: a[i]=PetscRealPart(eps->eigr[i]);
263: b[i]=eps->errest[i];
264: }
265: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
266: DSGetArray(eps->ds,DS_MAT_Q,&Q);
267: for ( i=eps->nconv;i<nv;i++) {
268: p=iwork[i];
269: if (p!=i) {
270: j=i+1;
271: while (iwork[j]!=i) j++;
272: iwork[j]=p;iwork[i]=i;
273: for (k=0;k<nv;k++) {
274: rtmp=Q[k+p*ld];Q[k+p*ld]=Q[k+i*ld];Q[k+i*ld]=rtmp;
275: }
276: }
277: }
278: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
279: k=eps->nconv+conv;
281: /* Checking values obtained for completing */
282: for (i=0;i<k;i++) {
283: sr->back[i]=eps->eigr[i];
284: }
285: STBackTransform(eps->OP,k,sr->back,eps->eigi);
286: count0=count1=0;
287: for (i=0;i<k;i++) {
288: lambda = PetscRealPart(sr->back[i]);
289: if ( ((sr->dir)*(sPres->value - lambda) > 0) && ((sr->dir)*(lambda - sPres->ext[0]) > 0))count0++;
290: if ( ((sr->dir)*(lambda - sPres->value) > 0) && ((sr->dir)*(sPres->ext[1] - lambda) > 0))count1++;
291: }
292:
293: /* Checks completion */
294: if ( (!sch0||count0 >= sPres->nsch[0]) && (!sch1 ||count1 >= sPres->nsch[1]) ) {
295: eps->reason = EPS_CONVERGED_TOL;
296: }else {
297: if (!complIterating && eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
298: if (complIterating) {
299: if (--iterCompl <= 0) eps->reason = EPS_DIVERGED_ITS;
300: }else if (k >= eps->nev) {
301: n0 = sPres->nsch[0]-count0;
302: n1 = sPres->nsch[1]-count1;
303: if ( sr->iterCompl>0 && ( (n0>0 && n0<= sr->nMAXCompl)||(n1>0&&n1<=sr->nMAXCompl) )) {
304: /* Iterating for completion*/
305: complIterating = PETSC_TRUE;
306: if (n0 >sr->nMAXCompl)sch0 = PETSC_FALSE;
307: if (n1 >sr->nMAXCompl)sch1 = PETSC_FALSE;
308: iterCompl = sr->iterCompl;
309: }else eps->reason = EPS_CONVERGED_TOL;
310: }
311: }
312: /* Update l */
313: if (eps->reason == EPS_CONVERGED_ITERATING) l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
314: else l = nv-k;
315: if (breakdown)l=0;
317: if (eps->reason == EPS_CONVERGED_ITERATING) {
318: if (breakdown) {
319: /* Start a new Lanczos factorization */
320: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
321: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
322: if (breakdown) {
323: eps->reason = EPS_DIVERGED_BREAKDOWN;
324: PetscInfo(eps,"Unable to generate more start vectors\n");
325: }
326: } else {
327: /* Prepare the Rayleigh quotient for restart */
328: DSGetArrayReal(eps->ds,DS_MAT_T,&a);
329: DSGetArray(eps->ds,DS_MAT_Q,&Q);
330: b = a + ld;
331: for (i=k;i<k+l;i++) {
332: a[i] = PetscRealPart(eps->eigr[i]);
333: b[i] = PetscRealPart(Q[nv-1+i*ld]*beta);
334: }
335: DSRestoreArrayReal(eps->ds,DS_MAT_T,&a);
336: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
337: }
338: }
339: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
340: DSGetArray(eps->ds,DS_MAT_Q,&Q);
341: SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
342: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
343: /* Normalize u and append it to V */
344: if ( eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
345: VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
346: }
347: /* Monitor */
348: if (eps->numbermonitors >0) {
349: aux = auxc = 0;
350: for (i=0;i<nv;i++) {
351: sr->back[i]=eps->eigr[i];
352: }
353: STBackTransform(eps->OP,nv,sr->back,eps->eigi);
354: for (i=0;i<nv;i++) {
355: lambda = PetscRealPart(sr->back[i]);
356: if ( ((sr->dir)*(lambda - sPres->ext[0]) > 0)&& ((sr->dir)*(sPres->ext[1] - lambda) > 0)) {
357: sr->monit[sr->indexEig+aux]=eps->eigr[i];
358: sr->errest[sr->indexEig+aux]=eps->errest[i];
359: aux++;
360: if (eps->errest[i] < eps->tol)auxc++;
361: }
362: }
363: EPSMonitor(eps,eps->its,auxc+sr->indexEig,sr->monit,sr->eigi,sr->errest,sr->indexEig+aux);
364: }
365: eps->nconv = k;
366: if (eps->reason != EPS_CONVERGED_ITERATING) {
367: /* Store approximated values for next shift */
368: DSGetArray(eps->ds,DS_MAT_Q,&Q);
369: sr->nS = l;
370: for (i=0;i<l;i++) {
371: sr->S[i] = eps->eigr[i+k];/* Diagonal elements */
372: sr->S[i+l] = Q[nv-1+(i+k)*ld]*beta; /* Out of diagonal elements */
373: }
374: sr->beta = beta;
375: DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
376: }
377: }
378: /* Check for completion */
379: for (i=0;i< eps->nconv; i++) {
380: if ( (sr->dir)*PetscRealPart(eps->eigr[i])>0 )sPres->nconv[1]++;
381: else sPres->nconv[0]++;
382: }
383: sPres->comp[0] = (count0 >= sPres->nsch[0])?PETSC_TRUE:PETSC_FALSE;
384: sPres->comp[1] = (count1 >= sPres->nsch[1])?PETSC_TRUE:PETSC_FALSE;
385: if (count0 > sPres->nsch[0] || count1 > sPres->nsch[1])SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error in Spectrum Slicing!\nMismatch between number of values found and information from inertia");
386: PetscFree(iwork);
387: return(0);
388: }
390: /*
391: Obtains value of subsequent shift
392: */
395: static PetscErrorCode EPSGetNewShiftValue(EPS eps,PetscInt side,PetscReal *newS)
396: {
397: PetscReal lambda,d_prev;
398: PetscInt i,idxP;
399: SR sr;
400: shift sPres,s;
401: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
404: sr = ctx->sr;
405: sPres = sr->sPres;
406: if ( sPres->neighb[side]) {
407: /* Completing a previous interval */
408: if (!sPres->neighb[side]->neighb[side] && sPres->neighb[side]->nconv[side]==0) { /* One of the ends might be too far from eigenvalues */
409: if (side) *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[sr->indexEig-1]]))/2;
410: else *newS = (sPres->value + PetscRealPart(sr->eig[sr->perm[0]]))/2;
411: }else *newS=(sPres->value + sPres->neighb[side]->value)/2;
412: } else { /* (Only for side=1). Creating a new interval. */
413: if (sPres->neigs==0) {/* No value has been accepted*/
414: if (sPres->neighb[0]) {
415: /* Multiplying by 10 the previous distance */
416: *newS = sPres->value + 10*(sr->dir)*PetscAbsReal(sPres->value - sPres->neighb[0]->value);
417: sr->nleap++;
418: /* Stops when the interval is open and no values are found in the last 5 shifts (there might be infinite eigenvalues) */
419: if ( !sr->hasEnd && sr->nleap > 5)SETERRQ(((PetscObject)eps)->comm,1,"Unable to compute the wanted eigenvalues with open interval");
420: }else {/* First shift */
421: if (eps->nconv != 0) {
422: /* Unaccepted values give information for next shift */
423: idxP=0;/* Number of values left from shift */
424: for (i=0;i<eps->nconv;i++) {
425: lambda = PetscRealPart(eps->eigr[i]);
426: if ( (sr->dir)*(lambda - sPres->value) <0)idxP++;
427: else break;
428: }
429: /* Avoiding subtraction of eigenvalues (might be the same).*/
430: if (idxP>0) {
431: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[0]))/(idxP+0.3);
432: }else {
433: d_prev = PetscAbsReal(sPres->value - PetscRealPart(eps->eigr[eps->nconv-1]))/(eps->nconv+0.3);
434: }
435: *newS = sPres->value + ((sr->dir)*d_prev*eps->nev)/2;
436: } else {/* No values found, no information for next shift */
437: SETERRQ(((PetscObject)eps)->comm,1,"First shift renders no information");
438: }
439: }
440: } else {/* Accepted values found */
441: sr->nleap = 0;
442: /* Average distance of values in previous subinterval */
443: s = sPres->neighb[0];
444: while (s && PetscAbs(s->inertia - sPres->inertia)==0) {
445: s = s->neighb[0];/* Looking for previous shifts with eigenvalues within */
446: }
447: if (s) {
448: d_prev = PetscAbsReal( (sPres->value - s->value)/(sPres->inertia - s->inertia));
449: } else {/* First shift. Average distance obtained with values in this shift */
450: /* first shift might be too far from first wanted eigenvalue (no values found outside the interval)*/
451: 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) ) {
452: d_prev = PetscAbsReal( (PetscRealPart(sr->eig[sr->indexEig-1]) - PetscRealPart(sr->eig[0])))/(sPres->neigs+0.3);
453: } else {
454: d_prev = PetscAbsReal( PetscRealPart(sr->eig[sr->indexEig-1]) - sPres->value)/(sPres->neigs+0.3);
455: }
456: }
457: /* Average distance is used for next shift by adding it to value on the right or to shift */
458: if ( (sr->dir)*(PetscRealPart(sr->eig[sPres->index + sPres->neigs -1]) - sPres->value) >0) {
459: *newS = PetscRealPart(sr->eig[sPres->index + sPres->neigs -1])+ ((sr->dir)*d_prev*(eps->nev))/2;
460: } else {/* Last accepted value is on the left of shift. Adding to shift */
461: *newS = sPres->value + ((sr->dir)*d_prev*(eps->nev))/2;
462: }
463: }
464: /* End of interval can not be surpassed */
465: if ((sr->dir)*( sr->int1 - *newS) < 0) *newS = sr->int1;
466: }/* of neighb[side]==null */
467: return(0);
468: }
470: /*
471: Function for sorting an array of real values
472: */
475: static PetscErrorCode sortRealEigenvalues(PetscScalar *r,PetscInt *perm,PetscInt nr,PetscBool prev,PetscInt dir)
476: {
477: PetscReal re;
478: PetscInt i,j,tmp;
479:
481: if (!prev) for (i=0; i < nr; i++) { perm[i] = i; }
482: /* Insertion sort */
483: for (i=1; i < nr; i++) {
484: re = PetscRealPart(r[perm[i]]);
485: j = i-1;
486: while ( j>=0 && dir*(re - PetscRealPart(r[perm[j]])) <= 0 ) {
487: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
488: }
489: }
490: return(0);
491: }
493: /* Stores the pairs obtained since the last shift in the global arrays */
496: PetscErrorCode EPSStoreEigenpairs(EPS eps)
497: {
498: PetscErrorCode ierr;
499: PetscReal lambda,err,norm;
500: PetscInt i,count;
501: PetscBool iscayley;
502: SR sr;
503: shift sPres;
504: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
507: sr = ctx->sr;
508: sPres = sr->sPres;
509: sPres->index = sr->indexEig;
510: count = sr->indexEig;
511: /* Back-transform */
512: EPSBackTransform_Default(eps);
513: PetscObjectTypeCompare((PetscObject)eps->OP,STCAYLEY,&iscayley);
514: /* Sort eigenvalues */
515: sortRealEigenvalues(eps->eigr,eps->perm,eps->nconv,PETSC_FALSE,sr->dir);
516: /* Values stored in global array */
517: for ( i=0; i < eps->nconv ;i++ ) {
518: lambda = PetscRealPart(eps->eigr[eps->perm[i]]);
519: err = eps->errest[eps->perm[i]];
521: if ( (sr->dir)*(lambda - sPres->ext[0]) > 0 && (sr->dir)*(sPres->ext[1] - lambda) > 0 ) {/* Valid value */
522: if (count>=sr->numEigs) {/* Error found */
523: SETERRQ(((PetscObject)eps)->comm,1,"Unexpected error in Spectrum Slicing");
524: }
525: sr->eig[count] = lambda;
526: sr->errest[count] = err;
527: /* Explicit purification */
528: STApply(eps->OP,eps->V[eps->perm[i]],sr->V[count]);
529: IPNorm(eps->ip,sr->V[count],&norm);
530: VecScale(sr->V[count],1.0/norm);
531: count++;
532: }
533: }
534: sPres->neigs = count - sr->indexEig;
535: sr->indexEig = count;
536: /* Global ordering array updating */
537: sortRealEigenvalues(sr->eig,sr->perm,count,PETSC_TRUE,sr->dir);
538: return(0);
539: }
543: PetscErrorCode EPSLookForDeflation(EPS eps)
544: {
545: PetscReal val;
546: PetscInt i,count0=0,count1=0;
547: shift sPres;
548: PetscInt ini,fin,k,idx0,idx1;
549: SR sr;
550: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
553: sr = ctx->sr;
554: sPres = sr->sPres;
556: if (sPres->neighb[0]) ini = (sr->dir)*(sPres->neighb[0]->inertia - sr->inertia0);
557: else ini = 0;
558: fin = sr->indexEig;
559: /* Selection of ends for searching new values */
560: if (!sPres->neighb[0]) sPres->ext[0] = sr->int0;/* First shift */
561: else sPres->ext[0] = sPres->neighb[0]->value;
562: if (!sPres->neighb[1]) {
563: if (sr->hasEnd) sPres->ext[1] = sr->int1;
564: else sPres->ext[1] = (sr->dir > 0)?PETSC_MAX_REAL:PETSC_MIN_REAL;
565: }else sPres->ext[1] = sPres->neighb[1]->value;
566: /* Selection of values between right and left ends */
567: for (i=ini;i<fin;i++) {
568: val=PetscRealPart(sr->eig[sr->perm[i]]);
569: /* Values to the right of left shift */
570: if ( (sr->dir)*(val - sPres->ext[1]) < 0 ) {
571: if ((sr->dir)*(val - sPres->value) < 0)count0++;
572: else count1++;
573: }else break;
574: }
575: /* The number of values on each side are found */
576: if (sPres->neighb[0]) {
577: sPres->nsch[0] = (sr->dir)*(sPres->inertia - sPres->neighb[0]->inertia)-count0;
578: 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");
579: }else sPres->nsch[0] = 0;
581: if (sPres->neighb[1]) {
582: sPres->nsch[1] = (sr->dir)*(sPres->neighb[1]->inertia - sPres->inertia) - count1;
583: 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");
584: }else sPres->nsch[1] = (sr->dir)*(sr->inertia1 - sPres->inertia);
586: /* Completing vector of indexes for deflation */
587: idx0 = ini;
588: idx1 = ini+count0+count1;
589: k=0;
590: for (i=idx0;i<idx1;i++)sr->idxDef[k++]=sr->perm[i];
591: for (i=0;i<k;i++)sr->VDef[i]=sr->V[sr->idxDef[i]];
592: eps->defl = sr->VDef;
593: eps->nds = k;
594: return(0);
595: }
599: PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS eps)
600: {
601: PetscErrorCode ierr;
602: PetscInt i,lds;
603: PetscReal newS;
604: KSP ksp;
605: PC pc;
606: Mat F;
607: PetscReal *errest_left;
608: Vec t;
609: SR sr;
610: shift s;
611: EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
612:
614: PetscMalloc(sizeof(struct _n_SR),&sr);
615: ctx->sr = sr;
616: sr->itsKs = 0;
617: sr->nleap = 0;
618: sr->nMAXCompl = eps->nev/4;
619: sr->iterCompl = eps->max_it/4;
620: sr->sPres = PETSC_NULL;
621: sr->nS = 0;
622: lds = PetscMin(eps->mpd,eps->ncv);
623: /* Checking presence of ends and finding direction */
624: if ( eps->inta > PETSC_MIN_REAL) {
625: sr->int0 = eps->inta;
626: sr->int1 = eps->intb;
627: sr->dir = 1;
628: if (eps->intb >= PETSC_MAX_REAL) { /* Right-open interval */
629: sr->hasEnd = PETSC_FALSE;
630: sr->inertia1 = eps->n;
631: }else sr->hasEnd = PETSC_TRUE;
632: } else { /* Left-open interval */
633: sr->int0 = eps->intb;
634: sr->int1 = eps->inta;
635: sr->dir = -1;
636: sr->hasEnd = PETSC_FALSE;
637: sr->inertia1 = 0;
638: }
639: /* Array of pending shifts */
640: sr->maxPend = 100;/* Initial size */
641: PetscMalloc((sr->maxPend)*sizeof(shift),&sr->pending);
642: if (sr->hasEnd) {
643: STGetKSP(eps->OP, &ksp);
644: KSPGetPC(ksp, &pc);
645: PCFactorGetMatrix(pc,&F);
646: /* Not looking for values in b (just inertia).*/
647: MatGetInertia(F,&sr->inertia1,PETSC_NULL,PETSC_NULL);
648: PCReset(pc); /* avoiding memory leak */
649: }
650: sr->nPend = 0;
651: EPSCreateShift(eps,sr->int0,PETSC_NULL,PETSC_NULL);
652: EPSExtractShift(eps);
653: sr->s0 = sr->sPres;
654: sr->inertia0 = sr->s0->inertia;
655: sr->numEigs = (sr->dir)*(sr->inertia1 - sr->inertia0);
656: sr->indexEig = 0;
657: /* Only with eigenvalues present in the interval ...*/
658: if (sr->numEigs==0) {
659: eps->reason = EPS_CONVERGED_TOL;
660: PetscFree(sr->s0);
661: PetscFree(sr->pending);
662: PetscFree(sr);
663: return(0);
664: }
665: /* Memory reservation for eig, V and perm */
666: PetscMalloc(lds*lds*sizeof(PetscScalar),&sr->S);
667: PetscMemzero(sr->S,lds*lds*sizeof(PetscScalar));
668: PetscMalloc((sr->numEigs)*sizeof(PetscScalar),&sr->eig);
669: PetscMalloc((sr->numEigs)*sizeof(PetscScalar),&sr->eigi);
670: PetscMalloc((sr->numEigs+eps->ncv) *sizeof(PetscReal),&sr->errest);
671: PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscReal),&errest_left);
672: PetscMalloc((sr->numEigs+eps->ncv)*sizeof(PetscScalar),&sr->monit);
673: PetscMalloc((eps->ncv)*sizeof(PetscScalar),&sr->back);
674: for (i=0;i<sr->numEigs;i++) {sr->eigi[i]=0;sr->eig[i] = 0;}
675: for (i=0;i<sr->numEigs+eps->ncv;i++) {errest_left[i]=0;sr->errest[i]=0;sr->monit[i]=0;}
676: VecCreateMPI(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,&t);
677: VecDuplicateVecs(t,sr->numEigs,&sr->V);
678: VecDestroy(&t);
679: /* Vector for maintaining order of eigenvalues */
680: PetscMalloc((sr->numEigs)*sizeof(PetscInt),&sr->perm);
681: for (i=0;i< sr->numEigs;i++)sr->perm[i]=i;
682: /* Vectors for deflation */
683: PetscMalloc((sr->numEigs)*sizeof(PetscInt),&sr->idxDef);
684: PetscMalloc((sr->numEigs)*sizeof(Vec),&sr->VDef);
685: sr->indexEig = 0;
686: /* Main loop */
687: while (sr->sPres) {
688: /* Search for deflation */
689: EPSLookForDeflation(eps);
690: /* KrylovSchur */
691: EPSKrylovSchur_Slice(eps);
692:
693: EPSStoreEigenpairs(eps);
694: /* Select new shift */
695: if (!sr->sPres->comp[1]) {
696: EPSGetNewShiftValue(eps,1,&newS);
697: EPSCreateShift(eps,newS,sr->sPres,sr->sPres->neighb[1]);
698: }
699: if (!sr->sPres->comp[0]) {
700: /* Completing earlier interval */
701: EPSGetNewShiftValue(eps,0,&newS);
702: EPSCreateShift(eps,newS,sr->sPres->neighb[0],sr->sPres);
703: }
704: /* Preparing for a new search of values */
705: EPSExtractShift(eps);
706: }
708: /* Updating eps values prior to exit */
709:
710: VecDestroyVecs(eps->allocated_ncv,&eps->V);
711: eps->V = sr->V;
712: PetscFree(sr->S);
713: PetscFree(eps->eigr);
714: PetscFree(eps->eigi);
715: PetscFree(eps->errest);
716: PetscFree(eps->errest_left);
717: PetscFree(eps->perm);
718: eps->eigr = sr->eig;
719: eps->eigi = sr->eigi;
720: eps->errest = sr->errest;
721: eps->errest_left = errest_left;
722: eps->perm = sr->perm;
723: eps->ncv = eps->allocated_ncv = sr->numEigs;
724: eps->nconv = sr->indexEig;
725: eps->reason = EPS_CONVERGED_TOL;
726: eps->its = sr->itsKs;
727: eps->nds = 0;
728: eps->defl = PETSC_NULL;
729: eps->evecsavailable = PETSC_TRUE;
730: PetscFree(sr->VDef);
731: PetscFree(sr->idxDef);
732: PetscFree(sr->pending);
733: PetscFree(sr->monit);
734: PetscFree(sr->back);
735: /* Reviewing list of shifts to free memory */
736: s = sr->s0;
737: if (s) {
738: while (s->neighb[1]) {
739: s = s->neighb[1];
740: PetscFree(s->neighb[0]);
741: }
742: PetscFree(s);
743: }
744: PetscFree(sr);
745: return(0);
746: }