Subversion Repositories slepc-dev

Rev

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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 23
   SLEPc - Scalable Library for Eigenvalue Problem Computations
24
   Copyright (c) 2002-2009, 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;
74
 
1643 slepc 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 */
1617 slepc 78
  ierr = PetscMemzero(S,n*n*sizeof(PetscReal));CHKERRQ(ierr);
1542 slepc 79
 
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];
1625 slepc 85
  for (i=l;i<n-1;i++)
1542 slepc 86
    S[(n-1-i)+(n-1-i-1)*n] = e[i];
87
 
88
  /* Reduce (2,2)-block of flipped S to tridiagonal form */
1643 slepc 89
  lwork = PetscBLASIntCast(n_*n_-n_);
1542 slepc 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);
92
 
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];
99
 
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);
103
 
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];
113
 
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);
117
 
118
  PetscFunctionReturn(0);
1548 slepc 119
#endif
1542 slepc 120
}
121
 
122
#undef __FUNCT__  
1475 slepc 123
#define __FUNCT__ "EPSProjectedKSSym"
1474 slepc 124
/*
125
   EPSProjectedKSSym - Solves the projected eigenproblem in the Krylov-Schur
126
   method (symmetric case).
127
 
128
   On input:
1655 slepc 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)
1474 slepc 133
 
134
   On output:
1477 slepc 135
     eig is the sorted list of eigenvalues
1655 slepc 136
     Q is the eigenvector matrix (order n, leading dimension n)
1474 slepc 137
 
138
   Workspace:
1655 slepc 139
     work is workspace to store a real square matrix of order n
140
     perm is workspace to store 2n integers
1474 slepc 141
*/
1655 slepc 142
PetscErrorCode EPSProjectedKSSym(EPS eps,PetscInt n,PetscInt l,PetscReal *a,PetscReal *b,PetscScalar *eig,PetscScalar *Q,PetscReal *work,PetscInt *perm)
1474 slepc 143
{
144
  PetscErrorCode ierr;
1655 slepc 145
  PetscInt       i,j,k,p;
146
  PetscReal      rtmp,*Qreal = (PetscReal*)Q;
1474 slepc 147
 
148
  PetscFunctionBegin;
1655 slepc 149
  /* Compute eigendecomposition of projected matrix */
150
  ierr = ArrowTridFlip(n,l,a,b,Qreal,work);CHKERRQ(ierr);
1477 slepc 151
 
1655 slepc 152
  /* Sort eigendecomposition according to eps->which */
1658 slepc 153
  ierr = EPSSortEigenvaluesReal(n,a,eps->which,n,perm,work);CHKERRQ(ierr);
1616 slepc 154
  for (i=0;i<n;i++)
1655 slepc 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++) {
1669 slepc 164
        rtmp = Qreal[k+p*n]; Qreal[k+p*n] = Qreal[k+i*n]; Qreal[k+i*n] = rtmp;
1655 slepc 165
      }
166
    }
167
  }
168
 
1542 slepc 169
#if defined(PETSC_USE_COMPLEX)
1616 slepc 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];
1542 slepc 173
#endif
1474 slepc 174
  PetscFunctionReturn(0);
175
}
176
 
177
#undef __FUNCT__  
178
#define __FUNCT__ "EPSSolve_KRYLOVSCHUR_SYMM"
179
PetscErrorCode EPSSolve_KRYLOVSCHUR_SYMM(EPS eps)
180
{
181
  PetscErrorCode ierr;
1655 slepc 182
  PetscInt       i,k,l,lds,nv,m;
1474 slepc 183
  Vec            u=eps->work[1];
1655 slepc 184
  PetscScalar    *Q;
185
  PetscReal      *a,*b,*work,beta;
186
  PetscInt       *iwork;
1474 slepc 187
  PetscTruth     breakdown;
188
 
189
  PetscFunctionBegin;
1612 slepc 190
  lds = PetscMin(eps->nev+eps->mpd,eps->ncv);
1655 slepc 191
  ierr = PetscMalloc(lds*lds*sizeof(PetscReal),&work);CHKERRQ(ierr);
1612 slepc 192
  ierr = PetscMalloc(lds*lds*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
1622 slepc 193
  ierr = PetscMalloc(lds*sizeof(PetscReal),&a);CHKERRQ(ierr);  
194
  ierr = PetscMalloc(lds*sizeof(PetscReal),&b);CHKERRQ(ierr);  
1655 slepc 195
  ierr = PetscMalloc(2*lds*sizeof(PetscInt),&iwork);CHKERRQ(ierr);
1474 slepc 196
 
1610 slepc 197
  /* Get the starting Lanczos vector */
1474 slepc 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
 
1610 slepc 205
    /* Compute an nv-step Lanczos factorization */
1622 slepc 206
    m = PetscMin(eps->nconv+eps->mpd,eps->ncv);
1655 slepc 207
    ierr = EPSFullLanczos(eps,a+l,b+l,eps->V,eps->nconv+l,&m,u,&breakdown);CHKERRQ(ierr);
1622 slepc 208
    nv = m - eps->nconv;
1655 slepc 209
    beta = b[nv-1];
1474 slepc 210
 
211
    /* Solve projected problem and compute residual norm estimates */
1655 slepc 212
    ierr = EPSProjectedKSSym(eps,nv,l,a,b,eps->eigr+eps->nconv,Q,work,iwork);CHKERRQ(ierr);
1610 slepc 213
    for (i=0;i<nv;i++)
214
      eps->errest[i+eps->nconv] = beta*PetscAbsScalar(Q[(i+1)*nv-1]) / PetscAbsScalar(eps->eigr[i+eps->nconv]);
1474 slepc 215
 
216
    /* Check convergence */
217
    k = eps->nconv;
1610 slepc 218
    while (k<eps->nconv+nv && eps->errest[k]<eps->tol) k++;    
1474 slepc 219
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
220
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
221
 
222
    /* Update l */
223
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
1610 slepc 224
    else l = (eps->nconv+nv-k)/2;
1574 slepc 225
 
1474 slepc 226
    if (eps->reason == EPS_CONVERGED_ITERATING) {
227
      if (breakdown) {
1610 slepc 228
        /* Start a new Lanczos factorization */
1474 slepc 229
        PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
230
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
231
        if (breakdown) {
232
          eps->reason = EPS_DIVERGED_BREAKDOWN;
233
          PetscInfo(eps,"Unable to generate more start vectors\n");
234
        }
235
      } else {
236
        /* Prepare the Rayleigh quotient for restart */
1612 slepc 237
        for (i=0;i<l;i++) {
1655 slepc 238
          a[i] = eps->eigr[i+k];
239
          b[i] = Q[nv-1+(i+k-eps->nconv)*nv]*beta;
1474 slepc 240
        }
241
      }
242
    }
243
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1610 slepc 244
    ierr = SlepcUpdateVectors(nv,eps->V+eps->nconv,0,k+l-eps->nconv,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
1475 slepc 245
    /* Normalize u and append it to V */
1474 slepc 246
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
1475 slepc 247
      ierr = VecAXPBY(eps->V[k+l],1.0/beta,0.0,u);CHKERRQ(ierr);
1474 slepc 248
    }
1610 slepc 249
 
250
    EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv+eps->nconv);
1474 slepc 251
    eps->nconv = k;
252
 
253
  }
254
 
255
  ierr = PetscFree(Q);CHKERRQ(ierr);
1622 slepc 256
  ierr = PetscFree(a);CHKERRQ(ierr);
257
  ierr = PetscFree(b);CHKERRQ(ierr);
1474 slepc 258
  ierr = PetscFree(work);CHKERRQ(ierr);
1655 slepc 259
  ierr = PetscFree(iwork);CHKERRQ(ierr);
1474 slepc 260
  PetscFunctionReturn(0);
261
}
262