Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1474 slepc 1
/*                      
2
 
3
   SLEPc eigensolver: "krylovschur"
4
 
5
   Method: Krylov-Schur
6
 
7
   Algorithm:
8
 
9
       Single-vector Krylov-Schur method for both symmetric and non-symmetric
10
       problems.
11
 
12
   References:
13
 
14
       [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
15
           available at http://www.grycap.upv.es/slepc.
16
 
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.
19
 
20
   Last update: Oct 2006
21
 
22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23
      SLEPc - Scalable Library for Eigenvalue Problem Computations
24
      Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
25
 
26
      This file is part of SLEPc. See the README file for conditions of use
27
      and additional information.
28
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29
*/
30
 
1521 slepc 31
#include "private/epsimpl.h"                /*I "slepceps.h" I*/
1474 slepc 32
#include "slepcblaslapack.h"
33
 
34
#undef __FUNCT__  
1542 slepc 35
#define __FUNCT__ "ArrowTridFlip"
36
/*
37
   ArrowTridFlip - Solves the arrowhead-tridiagonal eigenproblem by flipping
38
   the matrix and tridiagonalizing the bottom part.
39
 
40
   On input:
41
     l is the size of diagonal part
42
     d contains diagonal elements (length n)
43
     e contains offdiagonal elements (length n-1)
44
 
45
   On output:
46
     d contains the eigenvalues in ascending order
47
     Q is the eigenvector matrix (order n)
48
 
49
   Workspace:
50
     S is workspace to store a copy of the full matrix (nxn reals)
51
*/
52
PetscErrorCode ArrowTridFlip(PetscInt n,PetscInt l,PetscReal *d,PetscReal *e,PetscReal *Q,PetscReal *S)
53
{
1548 slepc 54
#if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR) || defined(SLEPC_MISSING_LAPACK_STEQR)
55
  PetscFunctionBegin;
56
  SETERRQ(PETSC_ERR_SUP,"SYTRD/ORGTR/STEQR - Lapack routine is unavailable.");
57
#else
1542 slepc 58
  PetscInt       i,j;
59
  PetscBLASInt   n1,n2,lwork,info;
60
 
61
  PetscFunctionBegin;
62
 
63
  n1 = l+1;    /* size of leading block, including residuals */
64
  n2 = n-l-1;  /* size of trailing block */
65
 
66
 /* Clean matrix S */
67
  for (i=0;i<n;i++)
68
    for (j=0;j<n;j++)
69
      S[i+j*n] = 0.0;
70
 
71
  /* Flip matrix S, copying the values saved in Q */
72
  for (i=0;i<n;i++)
73
    S[(n-1-i)+(n-1-i)*n] = d[i];
74
  for (i=0;i<l;i++)
75
    S[(n-1-i)+(n-1-l)*n] = e[i];
76
  for (i=l;i<n;i++)
77
    S[(n-1-i)+(n-1-i-1)*n] = e[i];
78
 
79
  /* Reduce (2,2)-block of flipped S to tridiagonal form */
80
  lwork = n*n-n;
81
  LAPACKsytrd_("L",&n1,S+n2*(n+1),&n,d,e,Q,Q+n,&lwork,&info);
82
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);
83
 
84
  /* Flip back diag and subdiag, put them in d and e */
85
  for (i=0;i<n-1;i++) {
86
    d[n-i-1] = S[i+i*n];
87
    e[n-i-2] = S[i+1+i*n];
88
  }
89
  d[0] = S[n-1+(n-1)*n];
90
 
91
  /* Compute the orthogonal matrix used for tridiagonalization */
92
  LAPACKorgtr_("L",&n1,S+n2*(n+1),&n,Q,Q+n,&lwork,&info);
93
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);
94
 
95
  /* Create full-size Q, flipped back to original order */
96
  for (i=0;i<n;i++)
97
    for (j=0;j<n;j++)
98
      Q[i+j*n] = 0.0;
99
  for (i=n1;i<n;i++)
100
    Q[i+i*n] = 1.0;
101
  for (i=0;i<n1;i++)
102
    for (j=0;j<n1;j++)
103
      Q[i+j*n] = S[n-i-1+(n-j-1)*n];
104
 
105
  /* Solve the tridiagonal eigenproblem */
106
  LAPACKsteqr_("V",&n,d,e,Q,&n,S,&info);
107
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
108
 
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:
120
     l is the number of vectors kept in previous restart (0 means first restart)
1477 slepc 121
     S is the projected matrix (order n, leading dimension is lds)
1474 slepc 122
 
123
   On output:
1475 slepc 124
     S is diagonal with diagonal elements (eigenvalues) sorted appropriately
1477 slepc 125
     eig is the sorted list of eigenvalues
126
     Q is the eigenvector matrix (order n)
1474 slepc 127
 
128
   Workspace:
1477 slepc 129
     work is workspace to store 2n reals and 2n integers
1474 slepc 130
*/
1509 slepc 131
PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *eig,PetscScalar *Q,PetscReal *work)
1474 slepc 132
{
133
  PetscErrorCode ierr;
1509 slepc 134
  PetscInt       i,j;
1474 slepc 135
  PetscReal      *ritz = work;
1542 slepc 136
  PetscReal      *e = work+n;
1509 slepc 137
  PetscInt       *perm = ((PetscInt*)(work+n))+n;
1542 slepc 138
  PetscReal      *Sreal = (PetscReal*)S, *Qreal = (PetscReal*)Q;
1474 slepc 139
 
140
  PetscFunctionBegin;
1477 slepc 141
 
142
  /* Compute eigendecomposition of S, S <- Q S Q' */
1474 slepc 143
  if (l==0) {
1477 slepc 144
    ierr = EPSDenseTridiagonal(n,S,lds,ritz,Q);CHKERRQ(ierr);
1474 slepc 145
  } else {
1542 slepc 146
    for (i=0;i<n;i++)
147
      ritz[i] = S[i+i*lds];
148
    for (i=0;i<l;i++)
149
      e[i] = S[l+i*lds];
150
    for (i=l;i<n;i++)
151
      e[i] = S[i+1+i*lds];
152
    ierr = ArrowTridFlip(n,l,ritz,e,Qreal,Sreal);CHKERRQ(ierr);
153
#if defined(PETSC_USE_COMPLEX)
154
    for (j=n-1;j>=0;j--)
155
      for (i=n-1;i>=0;i--)
156
        Q[i+j*n] = Qreal[i+j*n];
157
#endif
1474 slepc 158
  }
1477 slepc 159
 
160
  /* Sort eigendecomposition according to eps->which */
1542 slepc 161
  ierr = EPSSortEigenvaluesReal(n,ritz,eps->which,n,perm,e);CHKERRQ(ierr);
1477 slepc 162
  for (i=0;i<n;i++)
163
    eig[i] = ritz[perm[i]];
164
  for (j=0;j<n;j++)
1474 slepc 165
    for (i=0;i<n;i++)
1477 slepc 166
      S[i+j*lds] = Q[i+j*n];
167
  for (j=0;j<n;j++)
1474 slepc 168
    for (i=0;i<n;i++)
1477 slepc 169
      Q[i+j*n] = S[i+perm[j]*lds];
170
 
171
  /* Rebuild S from eig */
172
  for (i=0;i<n;i++) {
173
    S[i+i*lds] = eig[i];
174
    for (j=i+1;j<n;j++)
175
      S[j+i*lds] = 0.0;
1474 slepc 176
  }
177
  PetscFunctionReturn(0);
178
}
179
 
180
#undef __FUNCT__  
181
#define __FUNCT__ "EPSSolve_KRYLOVSCHUR_SYMM"
182
PetscErrorCode EPSSolve_KRYLOVSCHUR_SYMM(EPS eps)
183
{
184
  PetscErrorCode ierr;
1509 slepc 185
  PetscInt       i,k,l,n,lwork;
1474 slepc 186
  Vec            u=eps->work[1];
187
  PetscScalar    *S=eps->T,*Q;
188
  PetscReal      beta,*work;
189
  PetscTruth     breakdown;
190
 
191
  PetscFunctionBegin;
192
  ierr = PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));CHKERRQ(ierr);
193
  ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
1503 slepc 194
  lwork = 2*eps->ncv*sizeof(PetscReal) + 2*eps->ncv*sizeof(PetscInt);
1474 slepc 195
  ierr = PetscMalloc(lwork,&work);CHKERRQ(ierr);
196
 
197
  /* Get the starting Arnoldi vector */
198
  ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
199
  l = 0;
200
 
201
  /* Restart loop */
202
  while (eps->reason == EPS_CONVERGED_ITERATING) {
203
    eps->its++;
204
 
205
    /* Compute an nv-step Arnoldi factorization */
1576 slepc 206
    eps->nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
1574 slepc 207
    ierr = EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&eps->nv,u,&beta,&breakdown);CHKERRQ(ierr);
1474 slepc 208
 
209
    /* Solve projected problem and compute residual norm estimates */
210
    n = eps->nv-eps->nconv;
1477 slepc 211
    ierr = EPSProjectedKSSym(eps,n,l,S+eps->nconv*(eps->ncv+1),eps->ncv,eps->eigr+eps->nconv,Q,work);CHKERRQ(ierr);
1474 slepc 212
    for (i=eps->nconv;i<eps->nv;i++)
1477 slepc 213
      eps->errest[i] = beta*PetscAbsScalar(Q[(i-eps->nconv+1)*n-1]) / PetscAbsScalar(eps->eigr[i]);
1474 slepc 214
 
215
    /* Check convergence */
216
    k = eps->nconv;
217
    while (k<eps->nv && eps->errest[k]<eps->tol) k++;    
218
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
219
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
220
 
221
    /* Update l */
222
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
1475 slepc 223
    else l = (eps->nv-k)/2;
1574 slepc 224
 
1474 slepc 225
    if (eps->reason == EPS_CONVERGED_ITERATING) {
226
      if (breakdown) {
227
        /* Start a new Arnoldi factorization */
228
        PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
229
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
230
        if (breakdown) {
231
          eps->reason = EPS_DIVERGED_BREAKDOWN;
232
          PetscInfo(eps,"Unable to generate more start vectors\n");
233
        }
234
      } else {
235
        /* Prepare the Rayleigh quotient for restart */
236
        for (i=k;i<k+l;i++) {
1477 slepc 237
          S[i*eps->ncv+k+l] = Q[(i-eps->nconv+1)*n-1]*beta;
1474 slepc 238
        }
239
      }
240
    }
241
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1582 slepc 242
    ierr = EPSUpdateVectors(n,eps->V+eps->nconv,0,k+l-eps->nconv,Q,eps->AV);CHKERRQ(ierr);
1475 slepc 243
    /* Normalize u and append it to V */
1474 slepc 244
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1475 slepc 245
      ierr = VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);CHKERRQ(ierr);
1474 slepc 246
    }
247
    eps->nconv = k;
248
 
249
    EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nv);
250
 
251
  }
252
 
253
  ierr = PetscFree(Q);CHKERRQ(ierr);
254
  ierr = PetscFree(work);CHKERRQ(ierr);
255
  PetscFunctionReturn(0);
256
}
257