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
1887 jroman 1
/*
2
     This file contains some simple default routines for common QEP operations.  
3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
6
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
7
 
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/>.
21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22
*/
23
 
24
#include "private/qepimpl.h"   /*I "slepcqep.h" I*/
2066 jroman 25
#include "private/slepcimpl.h"
2044 antodo 26
#include "slepcblaslapack.h"
1887 jroman 27
 
28
#undef __FUNCT__  
29
#define __FUNCT__ "QEPDestroy_Default"
30
PetscErrorCode QEPDestroy_Default(QEP qep)
31
{
32
  PetscErrorCode ierr;
33
 
34
  PetscFunctionBegin;
35
  PetscValidHeaderSpecific(qep,QEP_COOKIE,1);
36
  ierr = PetscFree(qep->data);CHKERRQ(ierr);
37
 
38
  /* free work vectors */
39
  ierr = QEPDefaultFreeWork(qep);CHKERRQ(ierr);
40
  PetscFunctionReturn(0);
41
}
42
 
43
#undef __FUNCT__  
44
#define __FUNCT__ "QEPDefaultGetWork"
45
/*
46
  QEPDefaultGetWork - Gets a number of work vectors.
47
 */
48
PetscErrorCode QEPDefaultGetWork(QEP qep, PetscInt nw)
49
{
50
  PetscErrorCode ierr;
2044 antodo 51
  PetscInt       i;
1887 jroman 52
 
53
  PetscFunctionBegin;
54
 
55
  if (qep->nwork != nw) {
56
    if (qep->nwork > 0) {
57
      ierr = VecDestroyVecs(qep->work,qep->nwork); CHKERRQ(ierr);
58
    }
59
    qep->nwork = nw;
2044 antodo 60
    ierr = PetscMalloc(nw*sizeof(Vec),&qep->work);CHKERRQ(ierr);
61
    for (i=0;i<nw;i++) {
62
      ierr = MatGetVecs(qep->M,PETSC_NULL,qep->work+i); CHKERRQ(ierr);
63
    }
1887 jroman 64
    ierr = PetscLogObjectParents(qep,nw,qep->work);
65
  }
66
 
67
  PetscFunctionReturn(0);
68
}
69
 
70
#undef __FUNCT__  
71
#define __FUNCT__ "QEPDefaultFreeWork"
72
/*
73
  QEPDefaultFreeWork - Free work vectors.
74
 */
75
PetscErrorCode QEPDefaultFreeWork(QEP qep)
76
{
77
  PetscErrorCode ierr;
78
 
79
  PetscFunctionBegin;
80
  PetscValidHeaderSpecific(qep,QEP_COOKIE,1);
81
  if (qep->work)  {
82
    ierr = VecDestroyVecs(qep->work,qep->nwork);CHKERRQ(ierr);
83
  }
84
  PetscFunctionReturn(0);
85
}
86
 
87
#undef __FUNCT__  
88
#define __FUNCT__ "QEPDefaultConverged"
89
/*
2070 jroman 90
  QEPDefaultConverged - Checks convergence relative to the eigenvalue.
1887 jroman 91
*/
2070 jroman 92
PetscErrorCode QEPDefaultConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
1887 jroman 93
{
94
  PetscReal w;
95
  PetscFunctionBegin;
2044 antodo 96
  w = SlepcAbsEigenvalue(eigr,eigi);
2070 jroman 97
  *errest = res;
98
  if (w > res) *errest = res / w;
1887 jroman 99
  PetscFunctionReturn(0);
100
}
101
 
102
#undef __FUNCT__  
103
#define __FUNCT__ "QEPAbsoluteConverged"
104
/*
2070 jroman 105
  QEPAbsoluteConverged - Checks convergence absolutely.
1887 jroman 106
*/
2070 jroman 107
PetscErrorCode QEPAbsoluteConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
1887 jroman 108
{
2044 antodo 109
  PetscFunctionBegin;
2070 jroman 110
  *errest = res;
2044 antodo 111
  PetscFunctionReturn(0);
112
}
113
 
114
#undef __FUNCT__  
115
#define __FUNCT__ "QEPComputeVectors_Schur"
116
PetscErrorCode QEPComputeVectors_Schur(QEP qep)
117
{
118
#if defined(SLEPC_MISSING_LAPACK_TREVC)
119
  SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
120
#else
121
  PetscErrorCode ierr;
122
  PetscInt       i;
123
  PetscBLASInt   ncv,nconv,mout,info,one = 1;
124
  PetscScalar    *Z,*work,tmp;
125
#if defined(PETSC_USE_COMPLEX)
126
  PetscReal      *rwork;
127
#else
128
  PetscReal      normi;
129
#endif
130
  PetscReal      norm;
1887 jroman 131
 
132
  PetscFunctionBegin;
2044 antodo 133
  ncv = PetscBLASIntCast(qep->ncv);
134
  nconv = PetscBLASIntCast(qep->nconv);
135
 
136
  ierr = PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);CHKERRQ(ierr);
137
  ierr = PetscMalloc(3*nconv*sizeof(PetscScalar),&work);CHKERRQ(ierr);
138
#if defined(PETSC_USE_COMPLEX)
139
  ierr = PetscMalloc(nconv*sizeof(PetscReal),&rwork);CHKERRQ(ierr);
140
#endif
141
 
142
  /* right eigenvectors */
143
#if !defined(PETSC_USE_COMPLEX)
144
  LAPACKtrevc_("R","A",PETSC_NULL,&nconv,qep->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
145
#else
146
  LAPACKtrevc_("R","A",PETSC_NULL,&nconv,qep->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
147
#endif
148
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);
149
 
150
  /* normalize eigenvectors */
151
  for (i=0;i<qep->nconv;i++) {
152
#if !defined(PETSC_USE_COMPLEX)
153
    if (qep->eigi[i] != 0.0) {
154
      norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
155
      normi = BLASnrm2_(&nconv,Z+(i+1)*nconv,&one);
156
      tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
157
      BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
158
      BLASscal_(&nconv,&tmp,Z+(i+1)*nconv,&one);
159
      i++;    
160
    } else
161
#endif
162
    {
163
      norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
164
      tmp = 1.0 / norm;
165
      BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
166
    }
1887 jroman 167
  }
2044 antodo 168
 
169
  /* AV = V * Z */
170
  ierr = SlepcUpdateVectors(qep->nconv,qep->V,0,qep->nconv,Z,qep->nconv,PETSC_FALSE);CHKERRQ(ierr);
171
 
172
  ierr = PetscFree(Z);CHKERRQ(ierr);
173
  ierr = PetscFree(work);CHKERRQ(ierr);
174
#if defined(PETSC_USE_COMPLEX)
175
  ierr = PetscFree(rwork);CHKERRQ(ierr);
176
#endif
177
   PetscFunctionReturn(0);
178
#endif
1887 jroman 179
}
180
 
2066 jroman 181
#undef __FUNCT__  
182
#define __FUNCT__ "QEPKrylovConvergence"
183
/*
184
   QEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
185
   for quadratic Krylov methods.
186
 
187
   Differences:
188
   - Always non-symmetric
189
   - Does not check for STSHIFT
190
   - No correction factor
191
   - No support for true residual
192
*/
193
PetscErrorCode QEPKrylovConvergence(QEP qep,PetscInt kini,PetscInt nits,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt nv,PetscReal beta,PetscInt *kout,PetscScalar *work)
194
{
195
  PetscErrorCode ierr;
196
  PetscInt       k,marker;
197
  PetscScalar    re,im,*Z,*work2;
198
  PetscReal      resnorm;
2070 jroman 199
  PetscTruth     iscomplex;
2066 jroman 200
 
201
  PetscFunctionBegin;
202
  Z = work; work2 = work+2*nv;
203
  marker = -1;
204
  for (k=kini;k<kini+nits;k++) {
205
    /* eigenvalue */
206
    re = qep->eigr[k];
207
    im = qep->eigi[k];
208
    iscomplex = PETSC_FALSE;
209
    if (k<nv-1 && S[k+1+k*lds] != 0.0) iscomplex = PETSC_TRUE;
210
    /* residual norm */
211
    ierr = DenseSelectedEvec(S,lds,Q,Z,k,iscomplex,nv,work2);CHKERRQ(ierr);
212
    if (iscomplex) resnorm = beta*SlepcAbsEigenvalue(Z[nv-1],Z[2*nv-1]);
213
    else resnorm = beta*PetscAbsScalar(Z[nv-1]);
214
    /* error estimate */
2070 jroman 215
    ierr = (*qep->conv_func)(qep,re,im,resnorm,&qep->errest[k],qep->conv_ctx);CHKERRQ(ierr);
216
    if (marker==-1 && qep->errest[k] >= qep->tol) marker = k;
2066 jroman 217
    if (iscomplex) { qep->errest[k+1] = qep->errest[k]; k++; }
218
    if (marker!=-1 && !qep->trackall) break;
219
  }
220
  if (marker!=-1) k = marker;
221
  *kout = k;
222
 
223
  PetscFunctionReturn(0);
224
}
225
 
226