Subversion Repositories slepc-dev

Rev

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;
53
 
54
  PetscFunctionBegin;
2404 jroman 55
  /* spectrum slicing requires special treatment of default values */
56
  if (eps->which==EPS_ALL) {
57
    if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(((PetscObject)eps)->comm,1,"Must define a computational interval when using EPS_ALL");
58
    if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,1,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
59
    if (!eps->isgeneralized) SETERRQ(((PetscObject)eps)->comm,1,"Spectrum slicing not implemented for standard eigenproblems yet");
60
    if (!((PetscObject)(eps->OP))->type_name) { /* default to shift-and-invert */
61
      ierr = STSetType(eps->OP,STSINVERT);CHKERRQ(ierr);
62
    }
2459 carcamgo 63
    if(eps->intb >= PETSC_MAX_REAL){/* right-open interval */
64
      if(eps->inta <= PETSC_MIN_REAL) SETERRQ(((PetscObject)eps)->comm,1,"The defined computational interval should have at least one of their sides bounded");
65
      ierr = STSetDefaultShift(eps->OP,eps->inta);CHKERRQ(ierr);
66
    }else ierr = STSetDefaultShift(eps->OP,eps->intb);CHKERRQ(ierr);
67
 
2404 jroman 68
    if (eps->nev==1) eps->nev = 20;  /* nev not set, use default value */
69
    if (eps->nev<10) SETERRQ(((PetscObject)eps)->comm,1,"nev cannot be less than 10 in spectrum slicing runs");
70
    eps->mpd = 2*eps->nev;
71
    eps->ncv = 2*eps->nev;
2459 carcamgo 72
    eps->ops->backtransform = PETSC_NULL;
2404 jroman 73
  }
74
 
75
  /* proceed with the general case */
1576 slepc 76
  if (eps->ncv) { /* ncv set */
2214 jroman 77
    if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
2404 jroman 78
  } else if (eps->mpd) { /* mpd set */
1928 jroman 79
    eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
2404 jroman 80
  } else { /* neither set: defaults depend on nev being small or large */
1928 jroman 81
    if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
82
    else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
1576 slepc 83
  }
84
  if (!eps->mpd) eps->mpd = eps->ncv;
2214 jroman 85
  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 86
  if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
1942 jroman 87
  if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
1177 slepc 88
  if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
2214 jroman 89
    SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
1220 slepc 90
 
1560 slepc 91
  if (!eps->extraction) {
92
    ierr = EPSSetExtraction(eps,EPS_RITZ);CHKERRQ(ierr);
93
  } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
2214 jroman 94
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1426 slepc 95
  }
96
 
1171 slepc 97
  ierr = EPSAllocateSolution(eps);CHKERRQ(ierr);
98
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
1661 slepc 99
  if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
100
    ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);CHKERRQ(ierr);
1830 antodo 101
  }
1755 antodo 102
  ierr = EPSDefaultGetWork(eps,1);CHKERRQ(ierr);
1947 jroman 103
 
104
  /* dispatch solve method */
2214 jroman 105
  if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
1947 jroman 106
  if (eps->ishermitian) {
2404 jroman 107
    if (eps->which==EPS_ALL) eps->ops->solve = EPSSolve_KrylovSchur_Slice;
108
    else {
109
      switch (eps->extraction) {
110
        case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Symm; break;
111
        case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
112
        default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
113
      }
1947 jroman 114
    }
115
  } else {
116
    switch (eps->extraction) {
2404 jroman 117
      case EPS_RITZ:     eps->ops->solve = EPSSolve_KrylovSchur_Default; break;
2324 jroman 118
      case EPS_HARMONIC: eps->ops->solve = EPSSolve_KrylovSchur_Harmonic; break;
2214 jroman 119
      default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1947 jroman 120
    }
121
  }
1171 slepc 122
  PetscFunctionReturn(0);
123
}
124
 
125
#undef __FUNCT__  
1476 slepc 126
#define __FUNCT__ "EPSProjectedKSNonsym"
1424 slepc 127
/*
128
   EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
129
   method (non-symmetric case).
130
 
131
   On input:
132
     l is the number of vectors kept in previous restart (0 means first restart)
133
     S is the projected matrix (leading dimension is lds)
134
 
135
   On output:
136
     S has (real) Schur form with diagonal blocks sorted appropriately
1494 slepc 137
     Q contains the corresponding Schur vectors (order n, leading dimension n)
1424 slepc 138
*/
1509 slepc 139
PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
1424 slepc 140
{
141
  PetscErrorCode ierr;
1509 slepc 142
  PetscInt       i;
1424 slepc 143
 
144
  PetscFunctionBegin;
1498 slepc 145
  if (l==0) {
1424 slepc 146
    ierr = PetscMemzero(Q,n*n*sizeof(PetscScalar));CHKERRQ(ierr);
147
    for (i=0;i<n;i++)
148
      Q[i*(n+1)] = 1.0;
149
  } else {
150
    /* Reduce S to Hessenberg form, S <- Q S Q' */
151
    ierr = EPSDenseHessenberg(n,eps->nconv,S,lds,Q);CHKERRQ(ierr);
152
  }
153
  /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
154
  ierr = EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);
155
  /* Sort the remaining columns of the Schur form */
1823 antodo 156
  ierr = EPSSortDenseSchur(eps,n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);    
1424 slepc 157
  PetscFunctionReturn(0);
158
}
159
 
160
#undef __FUNCT__  
2324 jroman 161
#define __FUNCT__ "EPSSolve_KrylovSchur_Default"
162
PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
1474 slepc 163
{
164
  PetscErrorCode ierr;
2059 jroman 165
  PetscInt       i,k,l,lwork,nv;
1755 antodo 166
  Vec            u=eps->work[0];
2059 jroman 167
  PetscScalar    *S=eps->T,*Q,*work;
168
  PetscReal      beta;
2216 jroman 169
  PetscBool      breakdown;
1171 slepc 170
 
171
  PetscFunctionBegin;
172
  ierr = PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));CHKERRQ(ierr);
1830 antodo 173
  ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
2059 jroman 174
  lwork = 7*eps->ncv;
1474 slepc 175
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
1171 slepc 176
 
177
  /* Get the starting Arnoldi vector */
178
  ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
1172 slepc 179
  l = 0;
1171 slepc 180
 
181
  /* Restart loop */
182
  while (eps->reason == EPS_CONVERGED_ITERATING) {
1220 slepc 183
    eps->its++;
1171 slepc 184
 
185
    /* Compute an nv-step Arnoldi factorization */
1830 antodo 186
    nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
187
    ierr = EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);CHKERRQ(ierr);
1171 slepc 188
    ierr = VecScale(u,1.0/beta);CHKERRQ(ierr);
189
 
1958 jroman 190
    /* Solve projected problem */
1830 antodo 191
    ierr = EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);CHKERRQ(ierr);
1171 slepc 192
 
2059 jroman 193
    /* Check convergence */
194
    ierr = EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,S,eps->ncv,Q,eps->V,nv,beta,1.0,&k,work);CHKERRQ(ierr);
1172 slepc 195
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
196
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
197
 
1175 slepc 198
    /* Update l */
1172 slepc 199
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
200
    else {
1830 antodo 201
      l = (nv-k)/2;
1172 slepc 202
#if !defined(PETSC_USE_COMPLEX)
1185 slepc 203
      if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
1830 antodo 204
        if (k+l<nv-1) l = l+1;
1468 slepc 205
        else l = l-1;
1172 slepc 206
      }
207
#endif
208
    }
1958 jroman 209
 
1172 slepc 210
    if (eps->reason == EPS_CONVERGED_ITERATING) {
1171 slepc 211
      if (breakdown) {
1474 slepc 212
        /* Start a new Arnoldi factorization */
2499 jroman 213
        ierr = PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);CHKERRQ(ierr);
1468 slepc 214
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
215
        if (breakdown) {
1172 slepc 216
          eps->reason = EPS_DIVERGED_BREAKDOWN;
2499 jroman 217
          ierr = PetscInfo(eps,"Unable to generate more start vectors\n");CHKERRQ(ierr);
1468 slepc 218
        }
1172 slepc 219
      } else {
1474 slepc 220
        /* Prepare the Rayleigh quotient for restart */
1468 slepc 221
        for (i=k;i<k+l;i++) {
1830 antodo 222
          S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
1468 slepc 223
        }
1171 slepc 224
      }
1175 slepc 225
    }
1467 slepc 226
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1830 antodo 227
    ierr = SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
1582 slepc 228
 
1467 slepc 229
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
230
      ierr = VecCopy(u,eps->V[k+l]);CHKERRQ(ierr);
231
    }
232
    eps->nconv = k;
233
 
2313 jroman 234
    ierr = EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);CHKERRQ(ierr);
1172 slepc 235
  }
1958 jroman 236
  ierr = PetscFree(Q);CHKERRQ(ierr);
1474 slepc 237
  ierr = PetscFree(work);CHKERRQ(ierr);
1171 slepc 238
  PetscFunctionReturn(0);
239
}
240
 
2352 jroman 241
#undef __FUNCT__  
242
#define __FUNCT__ "EPSReset_KrylovSchur"
243
PetscErrorCode EPSReset_KrylovSchur(EPS eps)
244
{
245
  PetscErrorCode ierr;
246
 
247
  PetscFunctionBegin;
248
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
249
  ierr = EPSReset_Default(eps);CHKERRQ(ierr);
250
  PetscFunctionReturn(0);
251
}
252
 
1171 slepc 253
EXTERN_C_BEGIN
254
#undef __FUNCT__  
2324 jroman 255
#define __FUNCT__ "EPSCreate_KrylovSchur"
256
PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
1171 slepc 257
{
258
  PetscFunctionBegin;
2324 jroman 259
  eps->ops->setup          = EPSSetUp_KrylovSchur;
2352 jroman 260
  eps->ops->reset          = EPSReset_KrylovSchur;
2317 jroman 261
  eps->ops->backtransform  = EPSBackTransform_Default;
262
  eps->ops->computevectors = EPSComputeVectors_Schur;
1171 slepc 263
  PetscFunctionReturn(0);
264
}
265
EXTERN_C_END
266