Actual source code: ks-symm.c
1: /*
3: SLEPc eigensolver: "krylovschur"
5: Method: Krylov-Schur for symmetric eigenproblems
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
12:
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
27: #include <private/epsimpl.h> /*I "slepceps.h" I*/
28: #include <slepcblaslapack.h>
32: /*
33: ArrowTridFlip - Solves the arrowhead-tridiagonal eigenproblem by flipping
34: the matrix and tridiagonalizing the bottom part.
36: On input:
37: l is the size of diagonal part
38: d contains diagonal elements (length n)
39: e contains offdiagonal elements (length n-1)
41: On output:
42: d contains the eigenvalues in ascending order
43: Q is the eigenvector matrix (order n)
45: Workspace:
46: S is workspace to store a copy of the full matrix (nxn reals)
47: */
48: PetscErrorCode ArrowTridFlip(PetscInt n_,PetscInt l,PetscReal *d,PetscReal *e,PetscReal *Q,PetscReal *S)
49: {
50: #if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR) || defined(SLEPC_MISSING_LAPACK_STEQR)
52: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYTRD/ORGTR/STEQR - Lapack routine is unavailable.");
53: #else
55: PetscInt i,j;
56: PetscBLASInt n,n1,n2,lwork,info;
59: PetscLogEventBegin(EPS_Dense,0,0,0,0);
60: n = PetscBLASIntCast(n_);
61: /* quick return */
62: if (n == 1) {
63: Q[0] = 1.0;
64: return(0);
65: }
66: n1 = PetscBLASIntCast(l+1); /* size of leading block, including residuals */
67: n2 = PetscBLASIntCast(n-l-1); /* size of trailing block */
68: PetscMemzero(S,n*n*sizeof(PetscReal));
70: /* Flip matrix S, copying the values saved in Q */
71: for (i=0;i<n;i++)
72: S[(n-1-i)+(n-1-i)*n] = d[i];
73: for (i=0;i<l;i++)
74: S[(n-1-i)+(n-1-l)*n] = e[i];
75: for (i=l;i<n-1;i++)
76: S[(n-1-i)+(n-1-i-1)*n] = e[i];
78: /* Reduce (2,2)-block of flipped S to tridiagonal form */
79: lwork = PetscBLASIntCast(n_*n_-n_);
80: LAPACKsytrd_("L",&n1,S+n2*(n+1),&n,d,e,Q,Q+n,&lwork,&info);
81: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);
83: /* Flip back diag and subdiag, put them in d and e */
84: for (i=0;i<n-1;i++) {
85: d[n-i-1] = S[i+i*n];
86: e[n-i-2] = S[i+1+i*n];
87: }
88: d[0] = S[n-1+(n-1)*n];
90: /* Compute the orthogonal matrix used for tridiagonalization */
91: LAPACKorgtr_("L",&n1,S+n2*(n+1),&n,Q,Q+n,&lwork,&info);
92: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);
94: /* Create full-size Q, flipped back to original order */
95: for (i=0;i<n;i++)
96: for (j=0;j<n;j++)
97: Q[i+j*n] = 0.0;
98: for (i=n1;i<n;i++)
99: Q[i+i*n] = 1.0;
100: for (i=0;i<n1;i++)
101: for (j=0;j<n1;j++)
102: Q[i+j*n] = S[n-i-1+(n-j-1)*n];
104: /* Solve the tridiagonal eigenproblem */
105: LAPACKsteqr_("V",&n,d,e,Q,&n,S,&info);
106: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
108: PetscLogEventEnd(EPS_Dense,0,0,0,0);
109: return(0);
110: #endif
111: }
115: /*
116: EPSProjectedKSSym - Solves the projected eigenproblem in the Krylov-Schur
117: method (symmetric case).
119: On input:
120: n is the matrix dimension
121: l is the number of vectors kept in previous restart
122: a contains diagonal elements (length n)
123: b contains offdiagonal elements (length n-1)
125: On output:
126: eig is the sorted list of eigenvalues
127: Q is the eigenvector matrix (order n, leading dimension n)
129: Workspace:
130: work is workspace to store a real square matrix of order n
131: perm is workspace to store 2n integers
132: */
133: PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscReal *a,PetscReal *b,PetscScalar *eig,PetscScalar *Q,PetscReal *work,PetscInt *perm)
134: {
136: PetscInt i,j,k,p;
137: PetscReal rtmp,*Qreal = (PetscReal*)Q;
140: /* Compute eigendecomposition of projected matrix */
141: ArrowTridFlip(n,l,a,b,Qreal,work);
143: /* Sort eigendecomposition according to eps->which */
144: EPSSortEigenvaluesReal(eps,n,a,perm);
145: for (i=0;i<n;i++)
146: eig[i] = a[perm[i]];
147: for (i=0;i<n;i++) {
148: p = perm[i];
149: if (p != i) {
150: j = i + 1;
151: while (perm[j] != i) j++;
152: perm[j] = p; perm[i] = i;
153: /* swap eigenvectors i and j */
154: for (k=0;k<n;k++) {
155: rtmp = Qreal[k+p*n]; Qreal[k+p*n] = Qreal[k+i*n]; Qreal[k+i*n] = rtmp;
156: }
157: }
158: }
160: #if defined(PETSC_USE_COMPLEX)
161: for (j=n-1;j>=0;j--)
162: for (i=n-1;i>=0;i--)
163: Q[i+j*n] = Qreal[i+j*n];
164: #endif
165: return(0);
166: }
170: PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS eps)
171: {
173: PetscInt i,k,l,lds,lt,nv,m,*iwork;
174: Vec u=eps->work[0];
175: PetscScalar *Q;
176: PetscReal *a,*b,*work,beta;
177: PetscBool breakdown;
180: lds = PetscMin(eps->mpd,eps->ncv);
181: PetscMalloc(lds*lds*sizeof(PetscReal),&work);
182: PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);
183: PetscMalloc(2*lds*sizeof(PetscInt),&iwork);
184: lt = PetscMin(eps->nev+eps->mpd,eps->ncv);
185: PetscMalloc(lt*sizeof(PetscReal),&a);
186: PetscMalloc(lt*sizeof(PetscReal),&b);
188: /* Get the starting Lanczos vector */
189: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
190: l = 0;
191:
192: /* Restart loop */
193: while (eps->reason == EPS_CONVERGED_ITERATING) {
194: eps->its++;
196: /* Compute an nv-step Lanczos factorization */
197: m = PetscMin(eps->nconv+eps->mpd,eps->ncv);
198: EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);
199: nv = m - eps->nconv;
200: beta = b[nv-1];
202: /* Solve projected problem and compute residual norm estimates */
203: EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);
205: /* Check convergence */
206: EPSKrylovConvergence(eps,PETSC_TRUE,eps->trackall,eps->nconv,nv,PETSC_NULL,nv,Q,eps->V+eps->nconv,nv,beta,1.0,&k,PETSC_NULL);
207: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
208: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
209:
210: /* Update l */
211: if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
212: else l = (eps->nconv+nv-k)/2;
214: if (eps->reason == EPS_CONVERGED_ITERATING) {
215: if (breakdown) {
216: /* Start a new Lanczos factorization */
217: PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
218: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
219: if (breakdown) {
220: eps->reason = EPS_DIVERGED_BREAKDOWN;
221: PetscInfo(eps,"Unable to generate more start vectors\n");
222: }
223: } else {
224: /* Prepare the Rayleigh quotient for restart */
225: for (i=0;i<l;i++) {
226: a[i] = PetscRealPart(eps->eigr[i+k]);
227: b[i] = PetscRealPart(Q[nv-1+(i+k-eps->nconv)*nv]*beta);
228: }
229: }
230: }
231: /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
232: SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);
233: /* Normalize u and append it to V */
234: if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
235: VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);
236: }
238: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv+eps->nconv);
239: eps->nconv = k;
240: }
241: PetscFree(Q);
242: PetscFree(a);
243: PetscFree(b);
244: PetscFree(work);
245: PetscFree(iwork);
246: return(0);
247: }