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