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
 
1673 slepc 20
   Last update: Feb 2009
1474 slepc 21
 
22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 23
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 24
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1474 slepc 25
 
1672 slepc 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.
31
 
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.
36
 
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/>.
1474 slepc 39
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40
*/
41
 
1521 slepc 42
#include "private/epsimpl.h"                /*I "slepceps.h" I*/
1474 slepc 43
#include "slepcblaslapack.h"
44
 
45
#undef __FUNCT__  
1542 slepc 46
#define __FUNCT__ "ArrowTridFlip"
47
/*
48
   ArrowTridFlip - Solves the arrowhead-tridiagonal eigenproblem by flipping
49
   the matrix and tridiagonalizing the bottom part.
50
 
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)
55
 
56
   On output:
57
     d contains the eigenvalues in ascending order
58
     Q is the eigenvector matrix (order n)
59
 
60
   Workspace:
61
     S is workspace to store a copy of the full matrix (nxn reals)
62
*/
1643 slepc 63
PetscErrorCode ArrowTridFlip(PetscInt n_,PetscInt l,PetscReal *d,PetscReal *e,PetscReal *Q,PetscReal *S)
1542 slepc 64
{
1548 slepc 65
#if defined(SLEPC_MISSING_LAPACK_SYTRD) || defined(SLEPC_MISSING_LAPACK_ORGTR) || defined(SLEPC_MISSING_LAPACK_STEQR)
66
  PetscFunctionBegin;
67
  SETERRQ(PETSC_ERR_SUP,"SYTRD/ORGTR/STEQR - Lapack routine is unavailable.");
68
#else
1617 slepc 69
  PetscErrorCode ierr;
1542 slepc 70
  PetscInt       i,j;
1643 slepc 71
  PetscBLASInt   n,n1,n2,lwork,info;
1542 slepc 72
 
73
  PetscFunctionBegin;
1840 antodo 74
  ierr = PetscLogEventBegin(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1542 slepc 75
 
1643 slepc 76
  n = PetscBLASIntCast(n_);
1718 antodo 77
  /* quick return */
78
  if (n == 1) {
79
    Q[0] = 1.0;
80
    PetscFunctionReturn(0);    
81
  }
1643 slepc 82
  n1 = PetscBLASIntCast(l+1);    /* size of leading block, including residuals */
83
  n2 = PetscBLASIntCast(n-l-1);  /* size of trailing block */
1617 slepc 84
  ierr = PetscMemzero(S,n*n*sizeof(PetscReal));CHKERRQ(ierr);
1542 slepc 85
 
86
  /* Flip matrix S, copying the values saved in Q */
87
  for (i=0;i<n;i++)
88
    S[(n-1-i)+(n-1-i)*n] = d[i];
89
  for (i=0;i<l;i++)
90
    S[(n-1-i)+(n-1-l)*n] = e[i];
1625 slepc 91
  for (i=l;i<n-1;i++)
1542 slepc 92
    S[(n-1-i)+(n-1-i-1)*n] = e[i];
93
 
94
  /* Reduce (2,2)-block of flipped S to tridiagonal form */
1643 slepc 95
  lwork = PetscBLASIntCast(n_*n_-n_);
1542 slepc 96
  LAPACKsytrd_("L",&n1,S+n2*(n+1),&n,d,e,Q,Q+n,&lwork,&info);
97
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSYTRD %d",info);
98
 
99
  /* Flip back diag and subdiag, put them in d and e */
100
  for (i=0;i<n-1;i++) {
101
    d[n-i-1] = S[i+i*n];
102
    e[n-i-2] = S[i+1+i*n];
103
  }
104
  d[0] = S[n-1+(n-1)*n];
105
 
106
  /* Compute the orthogonal matrix used for tridiagonalization */
107
  LAPACKorgtr_("L",&n1,S+n2*(n+1),&n,Q,Q+n,&lwork,&info);
108
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGTR %d",info);
109
 
110
  /* Create full-size Q, flipped back to original order */
111
  for (i=0;i<n;i++)
112
    for (j=0;j<n;j++)
113
      Q[i+j*n] = 0.0;
114
  for (i=n1;i<n;i++)
115
    Q[i+i*n] = 1.0;
116
  for (i=0;i<n1;i++)
117
    for (j=0;j<n1;j++)
118
      Q[i+j*n] = S[n-i-1+(n-j-1)*n];
119
 
120
  /* Solve the tridiagonal eigenproblem */
121
  LAPACKsteqr_("V",&n,d,e,Q,&n,S,&info);
122
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xSTEQR %d",info);
123
 
1840 antodo 124
  ierr = PetscLogEventEnd(EPS_Dense,0,0,0,0);CHKERRQ(ierr);
1542 slepc 125
  PetscFunctionReturn(0);
1548 slepc 126
#endif
1542 slepc 127
}
128
 
129
#undef __FUNCT__  
1475 slepc 130
#define __FUNCT__ "EPSProjectedKSSym"
1474 slepc 131
/*
132
   EPSProjectedKSSym - Solves the projected eigenproblem in the Krylov-Schur
133
   method (symmetric case).
134
 
135
   On input:
1655 slepc 136
     n is the matrix dimension
137
     l is the number of vectors kept in previous restart
138
     a contains diagonal elements (length n)
139
     b contains offdiagonal elements (length n-1)
1474 slepc 140
 
141
   On output:
1477 slepc 142
     eig is the sorted list of eigenvalues
1655 slepc 143
     Q is the eigenvector matrix (order n, leading dimension n)
1474 slepc 144
 
145
   Workspace:
1655 slepc 146
     work is workspace to store a real square matrix of order n
147
     perm is workspace to store 2n integers
1474 slepc 148
*/
1655 slepc 149
PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscReal *a,PetscReal *b,PetscScalar *eig,PetscScalar *Q,PetscReal *work,PetscInt *perm)
1474 slepc 150
{
151
  PetscErrorCode ierr;
1655 slepc 152
  PetscInt       i,j,k,p;
153
  PetscReal      rtmp,*Qreal = (PetscReal*)Q;
1474 slepc 154
 
155
  PetscFunctionBegin;
1655 slepc 156
  /* Compute eigendecomposition of projected matrix */
157
  ierr = ArrowTridFlip(n,l,a,b,Qreal,work);CHKERRQ(ierr);
1477 slepc 158
 
1655 slepc 159
  /* Sort eigendecomposition according to eps->which */
1782 antodo 160
  ierr = EPSSortEigenvaluesReal(eps,n,a,perm);CHKERRQ(ierr);
1616 slepc 161
  for (i=0;i<n;i++)
1655 slepc 162
    eig[i] = a[perm[i]];
163
  for (i=0;i<n;i++) {
164
    p = perm[i];
165
    if (p != i) {
166
      j = i + 1;
167
      while (perm[j] != i) j++;
168
      perm[j] = p; perm[i] = i;
169
      /* swap eigenvectors i and j */
170
      for (k=0;k<n;k++) {
1669 slepc 171
        rtmp = Qreal[k+p*n]; Qreal[k+p*n] = Qreal[k+i*n]; Qreal[k+i*n] = rtmp;
1655 slepc 172
      }
173
    }
174
  }
175
 
1542 slepc 176
#if defined(PETSC_USE_COMPLEX)
1616 slepc 177
  for (j=n-1;j>=0;j--)
178
    for (i=n-1;i>=0;i--)
179
      Q[i+j*n] = Qreal[i+j*n];
1542 slepc 180
#endif
1474 slepc 181
  PetscFunctionReturn(0);
182
}
183
 
184
#undef __FUNCT__  
185
#define __FUNCT__ "EPSSolve_KRYLOVSCHUR_SYMM"
186
PetscErrorCode EPSSolve_KRYLOVSCHUR_SYMM(EPS eps)
187
{
188
  PetscErrorCode ierr;
2059 jroman 189
  PetscInt       i,k,l,lds,lt,nv,m;
1755 antodo 190
  Vec            u=eps->work[0];
2059 jroman 191
  PetscScalar    *Q;
192
  PetscReal      *a,*b,*work,beta;
1655 slepc 193
  PetscInt       *iwork;
2059 jroman 194
  PetscTruth     breakdown;
1474 slepc 195
 
196
  PetscFunctionBegin;
1691 slepc 197
  lds = PetscMin(eps->mpd,eps->ncv);
1655 slepc 198
  ierr = PetscMalloc(lds*lds*sizeof(PetscReal),&work);CHKERRQ(ierr);
2062 jroman 199
  ierr = PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
1655 slepc 200
  ierr = PetscMalloc(2*lds*sizeof(PetscInt),&iwork);CHKERRQ(ierr);
1691 slepc 201
  lt = PetscMin(eps->nev+eps->mpd,eps->ncv);
202
  ierr = PetscMalloc(lt*sizeof(PetscReal),&a);CHKERRQ(ierr);  
203
  ierr = PetscMalloc(lt*sizeof(PetscReal),&b);CHKERRQ(ierr);  
1474 slepc 204
 
1610 slepc 205
  /* Get the starting Lanczos vector */
1474 slepc 206
  ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
207
  l = 0;
208
 
209
  /* Restart loop */
210
  while (eps->reason == EPS_CONVERGED_ITERATING) {
211
    eps->its++;
212
 
1610 slepc 213
    /* Compute an nv-step Lanczos factorization */
1622 slepc 214
    m = PetscMin(eps->nconv+eps->mpd,eps->ncv);
1655 slepc 215
    ierr = EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);CHKERRQ(ierr);
1830 antodo 216
    nv = m - eps->nconv;
217
    beta = b[nv-1];
1474 slepc 218
 
219
    /* Solve projected problem and compute residual norm estimates */
1830 antodo 220
    ierr = EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);CHKERRQ(ierr);
1474 slepc 221
 
222
    /* Check convergence */
2059 jroman 223
    ierr = EPSKrylovConvergence(eps,PETSC_TRUE,eps->nconv,nv,PETSC_NULL,nv,Q,eps->V+eps->nconv,nv,beta,1.0,&k,PETSC_NULL);CHKERRQ(ierr);
1474 slepc 224
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
225
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
226
 
227
    /* Update l */
228
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
1830 antodo 229
    else l = (eps->nconv+nv-k)/2;
1574 slepc 230
 
1474 slepc 231
    if (eps->reason == EPS_CONVERGED_ITERATING) {
232
      if (breakdown) {
1610 slepc 233
        /* Start a new Lanczos factorization */
1474 slepc 234
        PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
235
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
236
        if (breakdown) {
237
          eps->reason = EPS_DIVERGED_BREAKDOWN;
238
          PetscInfo(eps,"Unable to generate more start vectors\n");
239
        }
240
      } else {
241
        /* Prepare the Rayleigh quotient for restart */
1612 slepc 242
        for (i=0;i<l;i++) {
1682 slepc 243
          a[i] = PetscRealPart(eps->eigr[i+k]);
1830 antodo 244
          b[i] = PetscRealPart(Q[nv-1+(i+k-eps->nconv)*nv]*beta);
1474 slepc 245
        }
246
      }
247
    }
248
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1830 antodo 249
    ierr = SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
1475 slepc 250
    /* Normalize u and append it to V */
1474 slepc 251
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1475 slepc 252
      ierr = VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);CHKERRQ(ierr);
1474 slepc 253
    }
1610 slepc 254
 
1830 antodo 255
    EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv+eps->nconv);
1474 slepc 256
    eps->nconv = k;
257
 
258
  }
259
 
2062 jroman 260
  ierr = PetscFree(Q);CHKERRQ(ierr);
1622 slepc 261
  ierr = PetscFree(a);CHKERRQ(ierr);
262
  ierr = PetscFree(b);CHKERRQ(ierr);
1474 slepc 263
  ierr = PetscFree(work);CHKERRQ(ierr);
1655 slepc 264
  ierr = PetscFree(iwork);CHKERRQ(ierr);
1474 slepc 265
  PetscFunctionReturn(0);
266
}
267