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
538 dsic.upv.es!jroman 1
/*                      
1707 jroman 2
   Common subroutines for all Krylov-type solvers.
6 dsic.upv.es!jroman 3
 
1376 slepc 4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 6
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6 dsic.upv.es!jroman 22
*/
1376 slepc 23
 
1521 slepc 24
#include "private/epsimpl.h"                /*I "slepceps.h" I*/
2066 jroman 25
#include "private/slepcimpl.h"
6 dsic.upv.es!jroman 26
#include "slepcblaslapack.h"
27
 
28
#undef __FUNCT__  
431 dsic.upv.es!antodo 29
#define __FUNCT__ "EPSBasicArnoldi"
538 dsic.upv.es!jroman 30
/*
31
   EPSBasicArnoldi - Computes an m-step Arnoldi factorization. The first k
32
   columns are assumed to be locked and therefore they are not modified. On
33
   exit, the following relation is satisfied:
34
 
35
                    OP * V - V * H = f * e_m^T
36
 
37
   where the columns of V are the Arnoldi vectors (which are B-orthonormal),
591 dsic.upv.es!jroman 38
   H is an upper Hessenberg matrix, f is the residual vector and e_m is
39
   the m-th vector of the canonical basis. The vector f is B-orthogonal to
40
   the columns of V. On exit, beta contains the B-norm of f and the next
41
   Arnoldi vector can be computed as v_{m+1} = f / beta.
538 dsic.upv.es!jroman 42
*/
1574 slepc 43
PetscErrorCode EPSBasicArnoldi(EPS eps,PetscTruth trans,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscTruth *breakdown)
431 dsic.upv.es!antodo 44
{
476 dsic.upv.es!antodo 45
  PetscErrorCode ierr;
1755 antodo 46
  PetscInt       j,m = *M;
476 dsic.upv.es!antodo 47
  PetscReal      norm;
431 dsic.upv.es!antodo 48
 
49
  PetscFunctionBegin;
1538 slepc 50
 
574 dsic.upv.es!antodo 51
  for (j=k;j<m-1;j++) {
879 ono.com!jroman 52
    if (trans) { ierr = STApplyTranspose(eps->OP,V[j],V[j+1]);CHKERRQ(ierr); }
53
    else { ierr = STApply(eps->OP,V[j],V[j+1]);CHKERRQ(ierr); }
1755 antodo 54
    ierr = IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,PETSC_NULL,V,V[j+1],H+ldh*j,&norm,breakdown);CHKERRQ(ierr);
1574 slepc 55
    H[j+1+ldh*j] = norm;
1113 slepc 56
    if (*breakdown) {
1049 slepc 57
      *M = j+1;
58
      *beta = norm;
59
      PetscFunctionReturn(0);
828 dsic.upv.es!antodo 60
    } else {
61
      ierr = VecScale(V[j+1],1/norm);CHKERRQ(ierr);
549 dsic.upv.es!antodo 62
    }
431 dsic.upv.es!antodo 63
  }
1600 slepc 64
  if (trans) { ierr = STApplyTranspose(eps->OP,V[m-1],f);CHKERRQ(ierr); }
65
  else { ierr = STApply(eps->OP,V[m-1],f);CHKERRQ(ierr); }
1755 antodo 66
  ierr = IPOrthogonalize(eps->ip,eps->nds,eps->DS,m,PETSC_NULL,V,f,H+ldh*(m-1),beta,PETSC_NULL);CHKERRQ(ierr);
1615 slepc 67
 
431 dsic.upv.es!antodo 68
  PetscFunctionReturn(0);
69
}
70
 
71
#undef __FUNCT__  
2059 jroman 72
#define __FUNCT__ "EPSKrylovConvergence"
73
/*
74
   EPSKrylovConvergence - Implements the loop that checks for convergence
75
   in Krylov methods.
76
 
77
   Input Parameters:
78
     eps   - the eigensolver; some error estimates are updated in eps->errest
79
     issym - whether the projected problem is symmetric or not
80
     kini  - initial value of k (the loop variable)
81
     nits  - number of iterations of the loop
82
     S     - Schur form of projected matrix (not referenced if issym)
83
     lds   - leading dimension of S
84
     Q     - Schur vectors of projected matrix (eigenvectors if issym)
85
     V     - set of basis vectors (used only if trueresidual is activated)
86
     nv    - number of vectors to process (dimension of Q, columns of V)
87
     beta  - norm of f (the residual vector of the Arnoldi/Lanczos factorization)
88
     corrf - correction factor for residual estimates (only in harmonic KS)
89
 
90
   Output Parameters:
91
     kout  - the first index where the convergence test failed
92
 
93
   Workspace:
94
     work is workspace to store 5*nv scalars, nv booleans and nv reals (only if !issym)
95
*/
96
PetscErrorCode EPSKrylovConvergence(EPS eps,PetscTruth issym,PetscInt kini,PetscInt nits,PetscScalar *S,PetscInt lds,PetscScalar *Q,Vec *V,PetscInt nv,PetscReal beta,PetscReal corrf,PetscInt *kout,PetscScalar *work)
97
{
98
  PetscErrorCode ierr;
99
  PetscInt       k,marker;
2140 jroman 100
  PetscScalar    re,im,*Z = work,*work2 = work;
2059 jroman 101
  PetscReal      resnorm;
2070 jroman 102
  PetscTruth     iscomplex,isshift;
2059 jroman 103
 
104
  PetscFunctionBegin;
105
  if (!issym) { Z = work; work2 = work+2*nv; }
106
  ierr = PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isshift);CHKERRQ(ierr);
107
  marker = -1;
108
  for (k=kini;k<kini+nits;k++) {
109
    /* eigenvalue */
110
    re = eps->eigr[k];
111
    im = eps->eigi[k];
112
    if (eps->trueres || isshift) {
113
      ierr = STBackTransform(eps->OP,1,&re,&im);CHKERRQ(ierr);
114
    }
115
    iscomplex = PETSC_FALSE;
116
    if (!issym && k<nv-1 && S[k+1+k*lds] != 0.0) iscomplex = PETSC_TRUE;
117
    /* residual norm */
118
    if (issym) {
119
      resnorm = beta*PetscAbsScalar(Q[(k-kini+1)*nv-1]);
120
    } else {
2065 jroman 121
      ierr = DenseSelectedEvec(S,lds,Q,Z,k,iscomplex,nv,work2);CHKERRQ(ierr);
122
      if (iscomplex) resnorm = beta*SlepcAbsEigenvalue(Z[nv-1],Z[2*nv-1]);
123
      else resnorm = beta*PetscAbsScalar(Z[nv-1]);
2059 jroman 124
    }
125
    if (eps->trueres) {
126
      if (issym) Z = Q+(k-kini)*nv;
127
      ierr = EPSComputeTrueResidual(eps,re,im,Z,V,nv,&resnorm);CHKERRQ(ierr);
128
    }
129
    else resnorm *= corrf;
130
    /* error estimate */
2070 jroman 131
    ierr = (*eps->conv_func)(eps,re,im,resnorm,&eps->errest[k],eps->conv_ctx);CHKERRQ(ierr);
132
    if (marker==-1 && eps->errest[k] >= eps->tol) marker = k;
2059 jroman 133
    if (iscomplex) { eps->errest[k+1] = eps->errest[k]; k++; }
134
    if (marker!=-1 && !eps->trackall) break;
135
  }
136
  if (marker!=-1) k = marker;
137
  *kout = k;
138
 
139
  PetscFunctionReturn(0);
140
}
141
 
142
#undef __FUNCT__  
1707 jroman 143
#define __FUNCT__ "EPSFullLanczos"
1484 slepc 144
/*
1707 jroman 145
   EPSFullLanczos - Computes an m-step Lanczos factorization with full
146
   reorthogonalization.  At each Lanczos step, the corresponding Lanczos
147
   vector is orthogonalized with respect to all previous Lanczos vectors.
148
   This is equivalent to computing an m-step Arnoldi factorization and
149
   exploting symmetry of the operator.
1484 slepc 150
 
1707 jroman 151
   The first k columns are assumed to be locked and therefore they are
152
   not modified. On exit, the following relation is satisfied:
1484 slepc 153
 
1707 jroman 154
                    OP * V - V * T = f * e_m^T
1484 slepc 155
 
1707 jroman 156
   where the columns of V are the Lanczos vectors (which are B-orthonormal),
157
   T is a real symmetric tridiagonal matrix, f is the residual vector and e_m
158
   is the m-th vector of the canonical basis. The tridiagonal is stored as
159
   two arrays: alpha contains the diagonal elements, beta the off-diagonal.
160
   The vector f is B-orthogonal to the columns of V. On exit, the last element
161
   of beta contains the B-norm of f and the next Lanczos vector can be
162
   computed as v_{m+1} = f / beta(end).
1484 slepc 163
 
164
*/
1707 jroman 165
PetscErrorCode EPSFullLanczos(EPS eps,PetscReal *alpha,PetscReal *beta,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscTruth *breakdown)
1484 slepc 166
{
167
  PetscErrorCode ierr;
1707 jroman 168
  PetscInt       j,m = *M;
169
  PetscReal      norm;
1755 antodo 170
  PetscScalar    *hwork,lhwork[100];
1484 slepc 171
 
172
  PetscFunctionBegin;
1755 antodo 173
  if (m > 100) {
1707 jroman 174
    ierr = PetscMalloc((eps->nds+m)*sizeof(PetscScalar),&hwork);CHKERRQ(ierr);
175
  } else {
176
    hwork = lhwork;
1484 slepc 177
  }
178
 
1707 jroman 179
  for (j=k;j<m-1;j++) {
180
    ierr = STApply(eps->OP,V[j],V[j+1]);CHKERRQ(ierr);
1755 antodo 181
    ierr = IPOrthogonalize(eps->ip,eps->nds,eps->DS,j+1,PETSC_NULL,V,V[j+1],hwork,&norm,breakdown);CHKERRQ(ierr);
182
    alpha[j-k] = PetscRealPart(hwork[j]);
1707 jroman 183
    beta[j-k] = norm;
184
    if (*breakdown) {
185
      *M = j+1;
1755 antodo 186
      if (m > 100) {
1707 jroman 187
        ierr = PetscFree(hwork);CHKERRQ(ierr);
188
      }
189
      PetscFunctionReturn(0);
1083 slepc 190
    } else {
1707 jroman 191
      ierr = VecScale(V[j+1],1.0/norm);CHKERRQ(ierr);
1083 slepc 192
    }
431 dsic.upv.es!antodo 193
  }
1707 jroman 194
  ierr = STApply(eps->OP,V[m-1],f);CHKERRQ(ierr);
1755 antodo 195
  ierr = IPOrthogonalize(eps->ip,eps->nds,eps->DS,m,PETSC_NULL,V,f,hwork,&norm,PETSC_NULL);CHKERRQ(ierr);
196
  alpha[m-1-k] = PetscRealPart(hwork[m-1]);
1707 jroman 197
  beta[m-1-k] = norm;
431 dsic.upv.es!antodo 198
 
1755 antodo 199
  if (m > 100) {
1707 jroman 200
    ierr = PetscFree(hwork);CHKERRQ(ierr);
1484 slepc 201
  }
431 dsic.upv.es!antodo 202
  PetscFunctionReturn(0);
203
}
204
 
1083 slepc 205