| 1474 |
slepc |
1 |
/*
|
|
|
2 |
|
|
|
3 |
SLEPc eigensolver: "krylovschur"
|
|
|
4 |
|
| 2404 |
jroman |
5 |
Method: Krylov-Schur for symmetric eigenproblems
|
| 1474 |
slepc |
6 |
|
|
|
7 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
| 1672 |
slepc |
8 |
SLEPc - Scalable Library for Eigenvalue Problem Computations
|
| 2575 |
eromero |
9 |
Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
|
| 1474 |
slepc |
10 |
|
| 1672 |
slepc |
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.
|
|
|
16 |
|
|
|
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.
|
|
|
21 |
|
|
|
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/>.
|
| 1474 |
slepc |
24 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
25 |
*/
|
|
|
26 |
|
| 2283 |
jroman |
27 |
#include <private/epsimpl.h> /*I "slepceps.h" I*/
|
|
|
28 |
#include <slepcblaslapack.h>
|
| 1474 |
slepc |
29 |
|
|
|
30 |
#undef __FUNCT__
|
| 1542 |
slepc |
31 |
#define __FUNCT__ "ArrowTridFlip"
|
|
|
32 |
/*
|
|
|
33 |
ArrowTridFlip - Solves the arrowhead-tridiagonal eigenproblem by flipping
|
|
|
34 |
the matrix and tridiagonalizing the bottom part.
|
|
|
35 |
|
|
|
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)
|
|
|
40 |
|
|
|
41 |
On output:
|
|
|
42 |
d contains the eigenvalues in ascending order
|
|
|
43 |
Q is the eigenvector matrix (order n)
|
|
|
44 |
|
|
|
45 |
Workspace:
|
|
|
46 |
S is workspace to store a copy of the full matrix (nxn reals)
|
|
|
47 |
*/
|
| 1643 |
slepc |
48 |
PetscErrorCode ArrowTridFlip(PetscInt n_,PetscInt l,PetscReal *d,PetscReal *e,PetscReal *Q,PetscReal *S)
|
| 1542 |
slepc |
49 |
{
|
| 1548 |
slepc |
50 |
#if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR) || defined(SLEPC_MISSING_LAPACK_STEQR)
|
|
|
51 |
PetscFunctionBegin;
|
| 2214 |
jroman |
52 |
SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYTRD/ORGTR/STEQR - Lapack routine is unavailable.");
|
| 1548 |
slepc |
53 |
#else
|
| 1617 |
slepc |
54 |
PetscErrorCode ierr;
|
| 1542 |
slepc |
55 |
PetscInt i,j;
|
| 1643 |
slepc |
56 |
PetscBLASInt n,n1,n2,lwork,info;
|
| 1542 |
slepc |
57 |
|
|
|
58 |
PetscFunctionBegin;
|
| 1840 |
antodo |
59 |
ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
|
| 1643 |
slepc |
60 |
n = PetscBLASIntCast(n_);
|
| 1718 |
antodo |
61 |
/* quick return */
|
|
|
62 |
if (n == 1) {
|
|
|
63 |
Q[0] = 1.0;
|
|
|
64 |
PetscFunctionReturn(0);
|
|
|
65 |
}
|
| 1643 |
slepc |
66 |
n1 = PetscBLASIntCast(l+1); /* size of leading block, including residuals */
|
|
|
67 |
n2 = PetscBLASIntCast(n-l-1); /* size of trailing block */
|
| 1617 |
slepc |
68 |
ierr = PetscMemzero(S,n*n*sizeof(PetscReal));CHKERRQ(ierr);
|
| 1542 |
slepc |
69 |
|
|
|
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];
|
| 1625 |
slepc |
75 |
for (i=l;i<n-1;i++)
|
| 1542 |
slepc |
76 |
S[(n-1-i)+(n-1-i-1)*n] = e[i];
|
|
|
77 |
|
|
|
78 |
/* Reduce (2,2)-block of flipped S to tridiagonal form */
|
| 1643 |
slepc |
79 |
lwork = PetscBLASIntCast(n_*n_-n_);
|
| 1542 |
slepc |
80 |
LAPACKsytrd_("L",&n1,S+n2*(n+1),&n,d,e,Q,Q+n,&lwork,&info);
|
| 2214 |
jroman |
81 |
if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);
|
| 1542 |
slepc |
82 |
|
|
|
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];
|
|
|
89 |
|
|
|
90 |
/* Compute the orthogonal matrix used for tridiagonalization */
|
|
|
91 |
LAPACKorgtr_("L",&n1,S+n2*(n+1),&n,Q,Q+n,&lwork,&info);
|
| 2214 |
jroman |
92 |
if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);
|
| 1542 |
slepc |
93 |
|
|
|
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];
|
|
|
103 |
|
|
|
104 |
/* Solve the tridiagonal eigenproblem */
|
|
|
105 |
LAPACKsteqr_("V",&n,d,e,Q,&n,S,&info);
|
| 2214 |
jroman |
106 |
if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
|
| 1542 |
slepc |
107 |
|
| 1840 |
antodo |
108 |
ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
|
| 1542 |
slepc |
109 |
PetscFunctionReturn(0);
|
| 1548 |
slepc |
110 |
#endif
|
| 1542 |
slepc |
111 |
}
|
|
|
112 |
|
|
|
113 |
#undef __FUNCT__
|
| 1475 |
slepc |
114 |
#define __FUNCT__ "EPSProjectedKSSym"
|
| 1474 |
slepc |
115 |
/*
|
|
|
116 |
EPSProjectedKSSym - Solves the projected eigenproblem in the Krylov-Schur
|
|
|
117 |
method (symmetric case).
|
|
|
118 |
|
|
|
119 |
On input:
|
| 1655 |
slepc |
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)
|
| 1474 |
slepc |
124 |
|
|
|
125 |
On output:
|
| 1477 |
slepc |
126 |
eig is the sorted list of eigenvalues
|
| 1655 |
slepc |
127 |
Q is the eigenvector matrix (order n, leading dimension n)
|
| 1474 |
slepc |
128 |
|
|
|
129 |
Workspace:
|
| 1655 |
slepc |
130 |
work is workspace to store a real square matrix of order n
|
|
|
131 |
perm is workspace to store 2n integers
|
| 1474 |
slepc |
132 |
*/
|
| 1655 |
slepc |
133 |
PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscReal *a,PetscReal *b,PetscScalar *eig,PetscScalar *Q,PetscReal *work,PetscInt *perm)
|
| 1474 |
slepc |
134 |
{
|
|
|
135 |
PetscErrorCode ierr;
|
| 1655 |
slepc |
136 |
PetscInt i,j,k,p;
|
|
|
137 |
PetscReal rtmp,*Qreal = (PetscReal*)Q;
|
| 1474 |
slepc |
138 |
|
|
|
139 |
PetscFunctionBegin;
|
| 1655 |
slepc |
140 |
/* Compute eigendecomposition of projected matrix */
|
|
|
141 |
ierr = ArrowTridFlip(n,l,a,b,Qreal,work);CHKERRQ(ierr);
|
| 1477 |
slepc |
142 |
|
| 1655 |
slepc |
143 |
/* Sort eigendecomposition according to eps->which */
|
| 1782 |
antodo |
144 |
ierr = EPSSortEigenvaluesReal(eps,n,a,perm);CHKERRQ(ierr);
|
| 1616 |
slepc |
145 |
for (i=0;i<n;i++)
|
| 1655 |
slepc |
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++) {
|
| 1669 |
slepc |
155 |
rtmp = Qreal[k+p*n]; Qreal[k+p*n] = Qreal[k+i*n]; Qreal[k+i*n] = rtmp;
|
| 1655 |
slepc |
156 |
}
|
|
|
157 |
}
|
|
|
158 |
}
|
|
|
159 |
|
| 1542 |
slepc |
160 |
#if defined(PETSC_USE_COMPLEX)
|
| 1616 |
slepc |
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];
|
| 1542 |
slepc |
164 |
#endif
|
| 1474 |
slepc |
165 |
PetscFunctionReturn(0);
|
|
|
166 |
}
|
|
|
167 |
|
|
|
168 |
#undef __FUNCT__
|
| 2324 |
jroman |
169 |
#define __FUNCT__ "EPSSolve_KrylovSchur_Symm"
|
|
|
170 |
PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS eps)
|
| 1474 |
slepc |
171 |
{
|
|
|
172 |
PetscErrorCode ierr;
|
| 2317 |
jroman |
173 |
PetscInt i,k,l,lds,lt,nv,m,*iwork;
|
| 1755 |
antodo |
174 |
Vec u=eps->work[0];
|
| 2059 |
jroman |
175 |
PetscScalar *Q;
|
|
|
176 |
PetscReal *a,*b,*work,beta;
|
| 2216 |
jroman |
177 |
PetscBool breakdown;
|
| 1474 |
slepc |
178 |
|
|
|
179 |
PetscFunctionBegin;
|
| 1691 |
slepc |
180 |
lds = PetscMin(eps->mpd,eps->ncv);
|
| 1655 |
slepc |
181 |
ierr = PetscMalloc(lds*lds*sizeof(PetscReal),&work);CHKERRQ(ierr);
|
| 2062 |
jroman |
182 |
ierr = PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
|
| 1655 |
slepc |
183 |
ierr = PetscMalloc(2*lds*sizeof(PetscInt),&iwork);CHKERRQ(ierr);
|
| 1691 |
slepc |
184 |
lt = PetscMin(eps->nev+eps->mpd,eps->ncv);
|
|
|
185 |
ierr = PetscMalloc(lt*sizeof(PetscReal),&a);CHKERRQ(ierr);
|
|
|
186 |
ierr = PetscMalloc(lt*sizeof(PetscReal),&b);CHKERRQ(ierr);
|
| 1474 |
slepc |
187 |
|
| 1610 |
slepc |
188 |
/* Get the starting Lanczos vector */
|
| 1474 |
slepc |
189 |
ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
|
|
|
190 |
l = 0;
|
|
|
191 |
|
|
|
192 |
/* Restart loop */
|
|
|
193 |
while (eps->reason == EPS_CONVERGED_ITERATING) {
|
|
|
194 |
eps->its++;
|
|
|
195 |
|
| 1610 |
slepc |
196 |
/* Compute an nv-step Lanczos factorization */
|
| 1622 |
slepc |
197 |
m = PetscMin(eps->nconv+eps->mpd,eps->ncv);
|
| 1655 |
slepc |
198 |
ierr = EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);CHKERRQ(ierr);
|
| 1830 |
antodo |
199 |
nv = m - eps->nconv;
|
|
|
200 |
beta = b[nv-1];
|
| 1474 |
slepc |
201 |
|
|
|
202 |
/* Solve projected problem and compute residual norm estimates */
|
| 1830 |
antodo |
203 |
ierr = EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);CHKERRQ(ierr);
|
| 1474 |
slepc |
204 |
|
|
|
205 |
/* Check convergence */
|
| 2583 |
jroman |
206 |
ierr = EPSKrylovConvergence(eps,PETSC_TRUE,eps->trackall,eps->nconv,nv,PETSC_NULL,nv,Q,eps->V+eps->nconv,nv,beta,1.0,&k,PETSC_NULL);CHKERRQ(ierr);
|
| 1474 |
slepc |
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;
|
| 1830 |
antodo |
212 |
else l = (eps->nconv+nv-k)/2;
|
| 1574 |
slepc |
213 |
|
| 1474 |
slepc |
214 |
if (eps->reason == EPS_CONVERGED_ITERATING) {
|
|
|
215 |
if (breakdown) {
|
| 1610 |
slepc |
216 |
/* Start a new Lanczos factorization */
|
| 2499 |
jroman |
217 |
ierr = PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);CHKERRQ(ierr);
|
| 1474 |
slepc |
218 |
ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
|
|
|
219 |
if (breakdown) {
|
|
|
220 |
eps->reason = EPS_DIVERGED_BREAKDOWN;
|
| 2499 |
jroman |
221 |
ierr = PetscInfo(eps,"Unable to generate more start vectors\n");CHKERRQ(ierr);
|
| 1474 |
slepc |
222 |
}
|
|
|
223 |
} else {
|
|
|
224 |
/* Prepare the Rayleigh quotient for restart */
|
| 1612 |
slepc |
225 |
for (i=0;i<l;i++) {
|
| 1682 |
slepc |
226 |
a[i] = PetscRealPart(eps->eigr[i+k]);
|
| 1830 |
antodo |
227 |
b[i] = PetscRealPart(Q[nv-1+(i+k-eps->nconv)*nv]*beta);
|
| 1474 |
slepc |
228 |
}
|
|
|
229 |
}
|
|
|
230 |
}
|
|
|
231 |
/* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
|
| 1830 |
antodo |
232 |
ierr = SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
|
| 1475 |
slepc |
233 |
/* Normalize u and append it to V */
|
| 1474 |
slepc |
234 |
if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
|
| 1475 |
slepc |
235 |
ierr = VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);CHKERRQ(ierr);
|
| 1474 |
slepc |
236 |
}
|
| 1610 |
slepc |
237 |
|
| 2313 |
jroman |
238 |
ierr = EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv+eps->nconv);CHKERRQ(ierr);
|
| 1474 |
slepc |
239 |
eps->nconv = k;
|
|
|
240 |
}
|
| 2062 |
jroman |
241 |
ierr = PetscFree(Q);CHKERRQ(ierr);
|
| 1622 |
slepc |
242 |
ierr = PetscFree(a);CHKERRQ(ierr);
|
|
|
243 |
ierr = PetscFree(b);CHKERRQ(ierr);
|
| 1474 |
slepc |
244 |
ierr = PetscFree(work);CHKERRQ(ierr);
|
| 1655 |
slepc |
245 |
ierr = PetscFree(iwork);CHKERRQ(ierr);
|
| 1474 |
slepc |
246 |
PetscFunctionReturn(0);
|
|
|
247 |
}
|
|
|
248 |
|