Actual source code: ks-symm.c

  1: /*                       

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur

  7:    Algorithm:

  9:        Single-vector Krylov-Schur method for both symmetric and non-symmetric
 10:        problems.

 12:    References:

 14:        [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7, 
 15:            available at http://www.grycap.upv.es/slepc.

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

 20:    Last update: Feb 2009

 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 24:    Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain

 26:    This file is part of SLEPc.
 27:       
 28:    SLEPc is free software: you can redistribute it and/or modify it under  the
 29:    terms of version 3 of the GNU Lesser General Public License as published by
 30:    the Free Software Foundation.

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

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

 42:  #include private/epsimpl.h
 43:  #include slepcblaslapack.h

 47: /*
 48:    ArrowTridFlip - Solves the arrowhead-tridiagonal eigenproblem by flipping
 49:    the matrix and tridiagonalizing the bottom part.

 51:    On input:
 52:      l is the size of diagonal part
 53:      d contains diagonal elements (length n)
 54:      e contains offdiagonal elements (length n-1)

 56:    On output:
 57:      d contains the eigenvalues in ascending order
 58:      Q is the eigenvector matrix (order n)

 60:    Workspace:
 61:      S is workspace to store a copy of the full matrix (nxn reals)
 62: */
 63: PetscErrorCode ArrowTridFlip(PetscInt n_,PetscInt l,PetscReal *d,PetscReal *e,PetscReal *Q,PetscReal *S)
 64: {
 65: #if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR) || defined(SLEPC_MISSING_LAPACK_STEQR)
 67:   SETERRQ(PETSC_ERR_SUP,"SYTRD/ORGTR/STEQR - Lapack routine is unavailable.");
 68: #else
 70:   PetscInt       i,j;
 71:   PetscBLASInt   n,n1,n2,lwork,info;


 75:   n = PetscBLASIntCast(n_);
 76:   n1 = PetscBLASIntCast(l+1);    /* size of leading block, including residuals */
 77:   n2 = PetscBLASIntCast(n-l-1);  /* size of trailing block */
 78:   PetscMemzero(S,n*n*sizeof(PetscReal));

 80:   /* Flip matrix S, copying the values saved in Q */
 81:   for (i=0;i<n;i++)
 82:     S[(n-1-i)+(n-1-i)*n] = d[i];
 83:   for (i=0;i<l;i++)
 84:     S[(n-1-i)+(n-1-l)*n] = e[i];
 85:   for (i=l;i<n-1;i++)
 86:     S[(n-1-i)+(n-1-i-1)*n] = e[i];

 88:   /* Reduce (2,2)-block of flipped S to tridiagonal form */
 89:   lwork = PetscBLASIntCast(n_*n_-n_);
 90:   LAPACKsytrd_("L",&n1,S+n2*(n+1),&n,d,e,Q,Q+n,&lwork,&info);
 91:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);

 93:   /* Flip back diag and subdiag, put them in d and e */
 94:   for (i=0;i<n-1;i++) {
 95:     d[n-i-1] = S[i+i*n];
 96:     e[n-i-2] = S[i+1+i*n];
 97:   }
 98:   d[0] = S[n-1+(n-1)*n];

100:   /* Compute the orthogonal matrix used for tridiagonalization */
101:   LAPACKorgtr_("L",&n1,S+n2*(n+1),&n,Q,Q+n,&lwork,&info);
102:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);

104:   /* Create full-size Q, flipped back to original order */
105:   for (i=0;i<n;i++)
106:     for (j=0;j<n;j++)
107:       Q[i+j*n] = 0.0;
108:   for (i=n1;i<n;i++)
109:     Q[i+i*n] = 1.0;
110:   for (i=0;i<n1;i++)
111:     for (j=0;j<n1;j++)
112:       Q[i+j*n] = S[n-i-1+(n-j-1)*n];

114:   /* Solve the tridiagonal eigenproblem */
115:   LAPACKsteqr_("V",&n,d,e,Q,&n,S,&info);
116:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);

118:   return(0);
119: #endif
120: }

124: /*
125:    EPSProjectedKSSym - Solves the projected eigenproblem in the Krylov-Schur
126:    method (symmetric case).

128:    On input:
129:      n is the matrix dimension
130:      l is the number of vectors kept in previous restart
131:      a contains diagonal elements (length n)
132:      b contains offdiagonal elements (length n-1)

134:    On output:
135:      eig is the sorted list of eigenvalues
136:      Q is the eigenvector matrix (order n, leading dimension n)

138:    Workspace:
139:      work is workspace to store a real square matrix of order n
140:      perm is workspace to store 2n integers
141: */
142: PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscReal *a,PetscReal *b,PetscScalar *eig,PetscScalar *Q,PetscReal *work,PetscInt *perm)
143: {
145:   PetscInt       i,j,k,p;
146:   PetscReal      rtmp,*Qreal = (PetscReal*)Q;

149:   /* Compute eigendecomposition of projected matrix */
150:   ArrowTridFlip(n,l,a,b,Qreal,work);

152:   /* Sort eigendecomposition according to eps->which */
153:   EPSSortEigenvaluesReal(n,a,eps->which,n,perm,work);
154:   for (i=0;i<n;i++)
155:     eig[i] = a[perm[i]];
156:   for (i=0;i<n;i++) {
157:     p = perm[i];
158:     if (p != i) {
159:       j = i + 1;
160:       while (perm[j] != i) j++;
161:       perm[j] = p; perm[i] = i;
162:       /* swap eigenvectors i and j */
163:       for (k=0;k<n;k++) {
164:         rtmp = Qreal[k+p*n]; Qreal[k+p*n] = Qreal[k+i*n]; Qreal[k+i*n] = rtmp;
165:       }
166:     }
167:   }

169: #if defined(PETSC_USE_COMPLEX)
170:   for (j=n-1;j>=0;j--)
171:     for (i=n-1;i>=0;i--)
172:       Q[i+j*n] = Qreal[i+j*n];
173: #endif
174:   return(0);
175: }

179: PetscErrorCode EPSSolve_KRYLOVSCHUR_SYMM(EPS eps)
180: {
182:   PetscInt       i,k,l,lds,lt,nv,m;
183:   Vec            u=eps->work[1];
184:   PetscScalar    *Q;
185:   PetscReal      *a,*b,*work,beta;
186:   PetscInt       *iwork;
187:   PetscTruth     breakdown;

190:   lds = PetscMin(eps->mpd,eps->ncv);
191:   PetscMalloc(lds*lds*sizeof(PetscReal),&work);
192:   PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);
193:   PetscMalloc(2*lds*sizeof(PetscInt),&iwork);
194:   lt = PetscMin(eps->nev+eps->mpd,eps->ncv);
195:   PetscMalloc(lt*sizeof(PetscReal),&a);
196:   PetscMalloc(lt*sizeof(PetscReal),&b);

198:   /* Get the starting Lanczos vector */
199:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
200:   l = 0;
201: 
202:   /* Restart loop */
203:   while (eps->reason == EPS_CONVERGED_ITERATING) {
204:     eps->its++;

206:     /* Compute an nv-step Lanczos factorization */
207:     m = PetscMin(eps->nconv+eps->mpd,eps->ncv);
208:     EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);
209:     nv = m - eps->nconv;
210:     beta = b[nv-1];

212:     /* Solve projected problem and compute residual norm estimates */
213:     EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);
214:     for (i=0;i<nv;i++)
215:       eps->errest[i+eps->nconv] = beta*PetscAbsScalar(Q[(i+1)*nv-1]) / PetscAbsScalar(eps->eigr[i+eps->nconv]);

217:     /* Check convergence */
218:     k = eps->nconv;
219:     while (k<eps->nconv+nv && eps->errest[k]<eps->tol) k++;
220:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
221:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
222: 
223:     /* Update l */
224:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
225:     else l = (eps->nconv+nv-k)/2;

227:     if (eps->reason == EPS_CONVERGED_ITERATING) {
228:       if (breakdown) {
229:         /* Start a new Lanczos factorization */
230:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
231:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
232:         if (breakdown) {
233:           eps->reason = EPS_DIVERGED_BREAKDOWN;
234:           PetscInfo(eps,"Unable to generate more start vectors\n");
235:         }
236:       } else {
237:         /* Prepare the Rayleigh quotient for restart */
238:         for (i=0;i<l;i++) {
239:           a[i] = PetscRealPart(eps->eigr[i+k]);
240:           b[i] = PetscRealPart(Q[nv-1+(i+k-eps->nconv)*nv]*beta);
241:         }
242:       }
243:     }
244:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
245:     SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);
246:     /* Normalize u and append it to V */
247:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
248:       VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
249:     }

251:     EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv+eps->nconv);
252:     eps->nconv = k;
253: 
254:   }

256:   PetscFree(Q);
257:   PetscFree(a);
258:   PetscFree(b);
259:   PetscFree(work);
260:   PetscFree(iwork);
261:   return(0);
262: }