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
 
1673 slepc 20
   Last update: Feb 2009
1217 slepc 21
 
1376 slepc 22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 23
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 24
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1376 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/>.
1376 slepc 39
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1217 slepc 40
*/
1376 slepc 41
 
1521 slepc 42
#include "private/epsimpl.h"                /*I "slepceps.h" I*/
1171 slepc 43
#include "slepcblaslapack.h"
44
 
1947 jroman 45
PetscErrorCode EPSSolve_KRYLOVSCHUR_DEFAULT(EPS);
46
extern PetscErrorCode EPSSolve_KRYLOVSCHUR_HARMONIC(EPS);
47
extern PetscErrorCode EPSSolve_KRYLOVSCHUR_SYMM(EPS);
1474 slepc 48
 
1171 slepc 49
#undef __FUNCT__  
50
#define __FUNCT__ "EPSSetUp_KRYLOVSCHUR"
51
PetscErrorCode EPSSetUp_KRYLOVSCHUR(EPS eps)
52
{
53
  PetscErrorCode ierr;
54
 
55
  PetscFunctionBegin;
1576 slepc 56
  if (eps->ncv) { /* ncv set */
2214 jroman 57
    if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
1171 slepc 58
  }
1576 slepc 59
  else if (eps->mpd) { /* mpd set */
1928 jroman 60
    eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
1576 slepc 61
  }
62
  else { /* neither set: defaults depend on nev being small or large */
1928 jroman 63
    if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
64
    else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
1576 slepc 65
  }
66
  if (!eps->mpd) eps->mpd = eps->ncv;
2214 jroman 67
  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 68
  if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
1942 jroman 69
  if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
1177 slepc 70
  if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
2214 jroman 71
    SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
1220 slepc 72
 
1560 slepc 73
  if (!eps->extraction) {
74
    ierr = EPSSetExtraction(eps,EPS_RITZ);CHKERRQ(ierr);
75
  } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC) {
2214 jroman 76
    SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1426 slepc 77
  }
78
 
1171 slepc 79
  ierr = EPSAllocateSolution(eps);CHKERRQ(ierr);
80
  ierr = PetscFree(eps->T);CHKERRQ(ierr);
1661 slepc 81
  if (!eps->ishermitian || eps->extraction==EPS_HARMONIC) {
82
    ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&eps->T);CHKERRQ(ierr);
1830 antodo 83
  }
1755 antodo 84
  ierr = EPSDefaultGetWork(eps,1);CHKERRQ(ierr);
1947 jroman 85
 
86
  /* dispatch solve method */
2214 jroman 87
  if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
1947 jroman 88
  if (eps->ishermitian) {
89
    switch (eps->extraction) {
90
      case EPS_RITZ:     eps->ops->solve = EPSSolve_KRYLOVSCHUR_SYMM; break;
91
      case EPS_HARMONIC: eps->ops->solve = EPSSolve_KRYLOVSCHUR_HARMONIC; break;
2214 jroman 92
      default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1947 jroman 93
    }
94
  } else {
95
    switch (eps->extraction) {
96
      case EPS_RITZ: eps->ops->solve = EPSSolve_KRYLOVSCHUR_DEFAULT; break;
97
      case EPS_HARMONIC: eps->ops->solve = EPSSolve_KRYLOVSCHUR_HARMONIC; break;
2214 jroman 98
      default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
1947 jroman 99
    }
100
  }
1171 slepc 101
  PetscFunctionReturn(0);
102
}
103
 
104
#undef __FUNCT__  
1476 slepc 105
#define __FUNCT__ "EPSProjectedKSNonsym"
1424 slepc 106
/*
107
   EPSProjectedKSNonsym - Solves the projected eigenproblem in the Krylov-Schur
108
   method (non-symmetric case).
109
 
110
   On input:
111
     l is the number of vectors kept in previous restart (0 means first restart)
112
     S is the projected matrix (leading dimension is lds)
113
 
114
   On output:
115
     S has (real) Schur form with diagonal blocks sorted appropriately
1494 slepc 116
     Q contains the corresponding Schur vectors (order n, leading dimension n)
1424 slepc 117
*/
1509 slepc 118
PetscErrorCode EPSProjectedKSNonsym(EPS eps,PetscInt l,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt n)
1424 slepc 119
{
120
  PetscErrorCode ierr;
1509 slepc 121
  PetscInt       i;
1424 slepc 122
 
123
  PetscFunctionBegin;
1498 slepc 124
  if (l==0) {
1424 slepc 125
    ierr = PetscMemzero(Q,n*n*sizeof(PetscScalar));CHKERRQ(ierr);
126
    for (i=0;i<n;i++)
127
      Q[i*(n+1)] = 1.0;
128
  } else {
129
    /* Reduce S to Hessenberg form, S <- Q S Q' */
130
    ierr = EPSDenseHessenberg(n,eps->nconv,S,lds,Q);CHKERRQ(ierr);
131
  }
132
  /* Reduce S to (quasi-)triangular form, S <- Q S Q' */
133
  ierr = EPSDenseSchur(n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);
134
  /* Sort the remaining columns of the Schur form */
1823 antodo 135
  ierr = EPSSortDenseSchur(eps,n,eps->nconv,S,lds,Q,eps->eigr,eps->eigi);CHKERRQ(ierr);    
1424 slepc 136
  PetscFunctionReturn(0);
137
}
138
 
139
#undef __FUNCT__  
1474 slepc 140
#define __FUNCT__ "EPSSolve_KRYLOVSCHUR_DEFAULT"
141
PetscErrorCode EPSSolve_KRYLOVSCHUR_DEFAULT(EPS eps)
142
{
143
  PetscErrorCode ierr;
2059 jroman 144
  PetscInt       i,k,l,lwork,nv;
1755 antodo 145
  Vec            u=eps->work[0];
2059 jroman 146
  PetscScalar    *S=eps->T,*Q,*work;
147
  PetscReal      beta;
2216 jroman 148
  PetscBool      breakdown;
1171 slepc 149
 
150
  PetscFunctionBegin;
151
  ierr = PetscMemzero(S,eps->ncv*eps->ncv*sizeof(PetscScalar));CHKERRQ(ierr);
1830 antodo 152
  ierr = PetscMalloc(eps->ncv*eps->ncv*sizeof(PetscScalar),&Q);CHKERRQ(ierr);
2059 jroman 153
  lwork = 7*eps->ncv;
1474 slepc 154
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
1171 slepc 155
 
156
  /* Get the starting Arnoldi vector */
157
  ierr = EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);CHKERRQ(ierr);
1172 slepc 158
  l = 0;
1171 slepc 159
 
160
  /* Restart loop */
161
  while (eps->reason == EPS_CONVERGED_ITERATING) {
1220 slepc 162
    eps->its++;
1171 slepc 163
 
164
    /* Compute an nv-step Arnoldi factorization */
1830 antodo 165
    nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
166
    ierr = EPSBasicArnoldi(eps,PETSC_FALSE,S,eps->ncv,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);CHKERRQ(ierr);
1171 slepc 167
    ierr = VecScale(u,1.0/beta);CHKERRQ(ierr);
168
 
1958 jroman 169
    /* Solve projected problem */
1830 antodo 170
    ierr = EPSProjectedKSNonsym(eps,l,S,eps->ncv,Q,nv);CHKERRQ(ierr);
1171 slepc 171
 
2059 jroman 172
    /* Check convergence */
173
    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 174
    if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
175
    if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
176
 
1175 slepc 177
    /* Update l */
1172 slepc 178
    if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
179
    else {
1830 antodo 180
      l = (nv-k)/2;
1172 slepc 181
#if !defined(PETSC_USE_COMPLEX)
1185 slepc 182
      if (S[(k+l-1)*(eps->ncv+1)+1] != 0.0) {
1830 antodo 183
        if (k+l<nv-1) l = l+1;
1468 slepc 184
        else l = l-1;
1172 slepc 185
      }
186
#endif
187
    }
1958 jroman 188
 
1172 slepc 189
    if (eps->reason == EPS_CONVERGED_ITERATING) {
1171 slepc 190
      if (breakdown) {
1474 slepc 191
        /* Start a new Arnoldi factorization */
1468 slepc 192
        PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%i norm=%g)\n",eps->its,beta);
193
        ierr = EPSGetStartVector(eps,k,eps->V[k],&breakdown);CHKERRQ(ierr);
194
        if (breakdown) {
1172 slepc 195
          eps->reason = EPS_DIVERGED_BREAKDOWN;
1468 slepc 196
          PetscInfo(eps,"Unable to generate more start vectors\n");
197
        }
1172 slepc 198
      } else {
1474 slepc 199
        /* Prepare the Rayleigh quotient for restart */
1468 slepc 200
        for (i=k;i<k+l;i++) {
1830 antodo 201
          S[i*eps->ncv+k+l] = Q[(i+1)*nv-1]*beta;
1468 slepc 202
        }
1171 slepc 203
      }
1175 slepc 204
    }
1467 slepc 205
    /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
1830 antodo 206
    ierr = SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,nv,PETSC_FALSE);CHKERRQ(ierr);
1582 slepc 207
 
1467 slepc 208
    if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
209
      ierr = VecCopy(u,eps->V[k+l]);CHKERRQ(ierr);
210
    }
211
    eps->nconv = k;
212
 
1830 antodo 213
    EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
1467 slepc 214
 
1172 slepc 215
  }
1171 slepc 216
 
1958 jroman 217
  ierr = PetscFree(Q);CHKERRQ(ierr);
1474 slepc 218
  ierr = PetscFree(work);CHKERRQ(ierr);
1171 slepc 219
  PetscFunctionReturn(0);
220
}
221
 
222
EXTERN_C_BEGIN
223
#undef __FUNCT__  
224
#define __FUNCT__ "EPSCreate_KRYLOVSCHUR"
225
PetscErrorCode EPSCreate_KRYLOVSCHUR(EPS eps)
226
{
227
  PetscFunctionBegin;
228
  eps->data                      = PETSC_NULL;
229
  eps->ops->setup                = EPSSetUp_KRYLOVSCHUR;
230
  eps->ops->setfromoptions       = PETSC_NULL;
231
  eps->ops->destroy              = EPSDestroy_Default;
232
  eps->ops->view                 = PETSC_NULL;
233
  eps->ops->backtransform        = EPSBackTransform_Default;
234
  eps->ops->computevectors       = EPSComputeVectors_Schur;
235
  PetscFunctionReturn(0);
236
}
237
EXTERN_C_END
238