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
1217 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
 
1376 slepc 20
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 21
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 22
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1376 slepc 23
 
1672 slepc 24
   This file is part of SLEPc.
25
 
26
   SLEPc is free software: you can redistribute it and/or modify it under  the
27
   terms of version 3 of the GNU Lesser General Public License as published by
28
   the Free Software Foundation.
29
 
30
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
31
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
32
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
33
   more details.
34
 
35
   You  should have received a copy of the GNU Lesser General  Public  License
36
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 37
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1217 slepc 38
*/
1376 slepc 39
 
2283 jroman 40
#include <private/epsimpl.h>                /*I "slepceps.h" I*/
41
#include <slepcblaslapack.h>
1171 slepc 42
 
2324 jroman 43
PetscErrorCode EPSSolve_KrylovSchur_Default(EPS);
44
extern PetscErrorCode EPSSolve_KrylovSchur_Harmonic(EPS);
45
extern PetscErrorCode EPSSolve_KrylovSchur_Symm(EPS);
2404 jroman 46
extern PetscErrorCode EPSSolve_KrylovSchur_Slice(EPS);
1474 slepc 47
 
1171 slepc 48
#undef __FUNCT__  
2324 jroman 49
#define __FUNCT__ "EPSSetUp_KrylovSchur"
50
PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
1171 slepc 51
{
52
  PetscErrorCode ierr;
2582 jroman 53
  PetscBool      issinv;
1171 slepc 54
 
55
  PetscFunctionBegin;
2404 jroman 56
  /* spectrum slicing requires special treatment of default values */
57
  if (eps->which==EPS_ALL) {
58
    if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(((PetscObject)eps)->comm,1,"Must define a computational interval when using EPS_ALL");
59
    if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,1,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
60
    if (!((PetscObject)(eps->OP))->type_name) { /* default to shift-and-invert */
61
      ierr = STSetType(eps->OP,STSINVERT);CHKERRQ(ierr);
62
    }
2582 jroman 63
    ierr = PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&issinv);CHKERRQ(ierr);
64
    if (!issinv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Shift-and-invert ST is needed for spectrum slicing");
2619 jroman 65
#if defined(PETSC_USE_REAL_DOUBLE)
66
    if (eps->tol==PETSC_DEFAULT) eps->tol = 1e-10;  /* use tighter tolerance */
67
#endif
2582 jroman 68
    if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
69
      if (eps->inta <= PETSC_MIN_REAL) SETERRQ(((PetscObject)eps)->comm,1,"The defined computational interval should have at least one of their sides bounded");
2459 carcamgo 70
      ierr = STSetDefaultShift(eps->OP,eps->inta);CHKERRQ(ierr);
2582 jroman 71
    }
72
    else { ierr = STSetDefaultShift(eps->OP,eps->intb);CHKERRQ(ierr); }
2459 carcamgo 73
 
2595 carcamgo 74
    if (eps->nev==1) eps->nev = 40;  /* nev not set, use default value */
2404 jroman 75
    if (eps->nev<10) SETERRQ(((PetscObject)eps)->comm,1,"nev cannot be less than 10 in spectrum slicing runs");
2459 carcamgo 76
    eps->ops->backtransform = PETSC_NULL;
2404 jroman 77
  }
78
 
79
  /* proceed with the general case */
1576 slepc 80
  if (eps->ncv) { /* ncv set */
2214 jroman 81
    if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
2404 jroman 82
  } else if (eps->mpd) { /* mpd set */
1928 jroman 83
    eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
2404 jroman 84
  } else { /* neither set: defaults depend on nev being small or large */
1928 jroman 85
    if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
86
    else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
1576 slepc 87
  }
88
  if (!eps->mpd) eps->mpd = eps->ncv;
2214 jroman 89
  if (eps->ncv>eps->nev+eps->mpd) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must not be larger than nev+mpd");
1928 jroman 90
  if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
2613 jroman 91
  if (!eps->which) { ierr = EPSDefaultSetWhich(eps);CHKERRQ(ierr); }
1177 slepc 92
  if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
2214 jroman 93
    SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
1220 slepc 94
 
1560 slepc 95
  if (!eps->extraction) {
96
    ierr = EPSSetExtraction(eps,EPS_RITZ);CHKERRQ(ierr);
97
  } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
2214 jroman 98
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1426 slepc 99
  }
100
 
1171 slepc 101
  ierr = EPSAllocateSolution(eps);CHKERRQ(ierr);
102
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
1661 slepc 103
  if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
104
    ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);CHKERRQ(ierr);
1830 antodo 105
  }
1755 antodo 106
  ierr = EPSDefaultGetWork(eps,1);CHKERRQ(ierr);
1947 jroman 107
 
108
  /* dispatch solve method */
2214 jroman 109
  if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
1947 jroman 110
  if (eps->ishermitian) {
2404 jroman 111
    if (eps->which==EPS_ALL) eps->ops->solve = EPSSolve_KrylovSchur_Slice;
112
    else {
113
      switch (eps->extraction) {
114
        case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Symm; break;
115
        case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
116
        default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
117
      }
1947 jroman 118
    }
119
  } else {
120
    switch (eps->extraction) {
2404 jroman 121
      case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Default; break;
2324 jroman 122
      case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
2214 jroman 123
      default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1947 jroman 124
    }
125
  }
1171 slepc 126
  PetscFunctionReturn(0);
127
}
128
 
129
#undef __FUNCT__  
1476 slepc 130
#define __FUNCT__ "EPSProjectedKSNonsym"
1424 slepc 131
/*
132
   EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
133
   method (non-symmetric case).
134
 
135
   On input:
136
     l is the number of vectors kept in previous restart (0 means first restart)
137
     S is the projected matrix (leading dimension is lds)
138
 
139
   On output:
140
     S has (real) Schur form with diagonal blocks sorted appropriately
1494 slepc 141
     Q contains the corresponding Schur vectors (order n, leading dimension n)
1424 slepc 142
*/
1509 slepc 143
PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
1424 slepc 144
{
145
  PetscErrorCode ierr;
1509 slepc 146
  PetscInt       i;
1424 slepc 147
 
148
  PetscFunctionBegin;
1498 slepc 149
  if (l==0) {
1424 slepc 150
    ierr = PetscMemzero(Q,n*n*sizeof(PetscScalar));CHKERRQ(ierr);
151
    for (i=0;i<n;i++)
152
      Q[i*(n+1)] = 1.0;
153
  } else {
154
    /* Reduce S to Hessenberg form, S <- Q S Q' */
155
    ierr = EPSDenseHessenberg(n,eps->nconv,S,lds,Q);CHKERRQ(ierr);
156
  }
157
  /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
158
  ierr = EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);
159
  /* Sort the remaining columns of the Schur form */
1823 antodo 160
  ierr = EPSSortDenseSchur(eps,n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);    
1424 slepc 161
  PetscFunctionReturn(0);
162
}
163
 
164
#undef __FUNCT__  
2324 jroman 165
#define __FUNCT__ "EPSSolve_KrylovSchur_Default"
166
PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
1474 slepc 167
{
168
  PetscErrorCode ierr;
2059 jroman 169
  PetscInt       i,k,l,lwork,nv;
1755 antodo 170
  Vec            u=eps->work[0];
2059 jroman 171
  PetscScalar    *S=eps->T,*Q,*work;
172
  PetscReal      beta;
2216 jroman 173
  PetscBool      breakdown;
1171 slepc 174
 
175
  PetscFunctionBegin;
176
  ierr = PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));CHKERRQ(ierr);
1830 antodo 177
  ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
2059 jroman 178
  lwork = 7*eps->ncv;
1474 slepc 179
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
1171 slepc 180
 
181
  /* Get the starting Arnoldi vector */
182
  ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
1172 slepc 183
  l = 0;
1171 slepc 184
 
185
  /* Restart loop */
186
  while (eps->reason == EPS_CONVERGED_ITERATING) {
1220 slepc 187
    eps->its++;
1171 slepc 188
 
189
    /* Compute an nv-step Arnoldi factorization */
1830 antodo 190
    nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
191
    ierr = EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);CHKERRQ(ierr);
1171 slepc 192
    ierr = VecScale(u,1.0/beta);CHKERRQ(ierr);
193
 
1958 jroman 194
    /* Solve projected problem */
1830 antodo 195
    ierr = EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);CHKERRQ(ierr);
1171 slepc 196
 
2059 jroman 197
    /* Check convergence */
2583 jroman 198
    ierr = EPSKrylovConvergence(eps,PETSC_FALSE,eps->trackall,eps->nconv,nv-eps->nconv,S,eps->ncv,Q,eps->V,nv,beta,1.0,&k,work);CHKERRQ(ierr);
1172 slepc 199
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
200
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
201
 
1175 slepc 202
    /* Update l */
1172 slepc 203
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
204
    else {
1830 antodo 205
      l = (nv-k)/2;
1172 slepc 206
#if !defined(PETSC_USE_COMPLEX)
1185 slepc 207
      if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
1830 antodo 208
        if (k+l<nv-1) l = l+1;
1468 slepc 209
        else l = l-1;
1172 slepc 210
      }
211
#endif
212
    }
1958 jroman 213
 
1172 slepc 214
    if (eps->reason == EPS_CONVERGED_ITERATING) {
1171 slepc 215
      if (breakdown) {
1474 slepc 216
        /* Start a new Arnoldi factorization */
2499 jroman 217
        ierr = PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);CHKERRQ(ierr);
1468 slepc 218
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
219
        if (breakdown) {
1172 slepc 220
          eps->reason = EPS_DIVERGED_BREAKDOWN;
2499 jroman 221
          ierr = PetscInfo(eps,"Unable to generate more start vectors\n");CHKERRQ(ierr);
1468 slepc 222
        }
1172 slepc 223
      } else {
1474 slepc 224
        /* Prepare the Rayleigh quotient for restart */
1468 slepc 225
        for (i=k;i<k+l;i++) {
1830 antodo 226
          S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
1468 slepc 227
        }
1171 slepc 228
      }
1175 slepc 229
    }
1467 slepc 230
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1830 antodo 231
    ierr = SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
1582 slepc 232
 
1467 slepc 233
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
234
      ierr = VecCopy(u,eps->V[k+l]);CHKERRQ(ierr);
235
    }
236
    eps->nconv = k;
237
 
2313 jroman 238
    ierr = EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);CHKERRQ(ierr);
1172 slepc 239
  }
1958 jroman 240
  ierr = PetscFree(Q);CHKERRQ(ierr);
1474 slepc 241
  ierr = PetscFree(work);CHKERRQ(ierr);
1171 slepc 242
  PetscFunctionReturn(0);
243
}
244
 
2352 jroman 245
#undef __FUNCT__  
246
#define __FUNCT__ "EPSReset_KrylovSchur"
247
PetscErrorCode EPSReset_KrylovSchur(EPS eps)
248
{
249
  PetscErrorCode ierr;
250
 
251
  PetscFunctionBegin;
252
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
253
  ierr = EPSReset_Default(eps);CHKERRQ(ierr);
254
  PetscFunctionReturn(0);
255
}
256
 
1171 slepc 257
EXTERN_C_BEGIN
258
#undef __FUNCT__  
2324 jroman 259
#define __FUNCT__ "EPSCreate_KrylovSchur"
260
PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1171 slepc 261
{
262
  PetscFunctionBegin;
2324 jroman 263
  eps->ops->setup          = EPSSetUp_KrylovSchur;
2352 jroman 264
  eps->ops->reset          = EPSReset_KrylovSchur;
2317 jroman 265
  eps->ops->backtransform  = EPSBackTransform_Default;
266
  eps->ops->computevectors = EPSComputeVectors_Schur;
1171 slepc 267
  PetscFunctionReturn(0);
268
}
269
EXTERN_C_END
270