Actual source code: qepdefault.c

  1: /*
  2:      This file contains some simple default routines for common QEP operations.  

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

  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.

 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.

 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: */

 24: #include <private/qepimpl.h>     /*I "slepcqep.h" I*/
 25: #include <private/slepcimpl.h>   /*I "slepcsys.h" I*/
 26: #include <slepcblaslapack.h>

 30: /*
 31:   QEPDefaultGetWork - Gets a number of work vectors.
 32:  */
 33: PetscErrorCode QEPDefaultGetWork(QEP qep,PetscInt nw)
 34: {

 38:   if (qep->nwork != nw) {
 39:     VecDestroyVecs(qep->nwork,&qep->work);
 40:     qep->nwork = nw;
 41:     VecDuplicateVecs(qep->t,nw,&qep->work);
 42:     PetscLogObjectParents(qep,nw,qep->work);
 43:   }
 44:   return(0);
 45: }

 49: /*
 50:   QEPDefaultFreeWork - Free work vectors.
 51:  */
 52: PetscErrorCode QEPDefaultFreeWork(QEP qep)
 53: {

 57:   VecDestroyVecs(qep->nwork,&qep->work);
 58:   return(0);
 59: }

 63: /*
 64:   QEPDefaultConverged - Checks convergence relative to the eigenvalue.
 65: */
 66: PetscErrorCode QEPDefaultConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 67: {
 68:   PetscReal w;

 71:   w = SlepcAbsEigenvalue(eigr,eigi);
 72:   *errest = res;
 73:   if (w > res) *errest = res / w;
 74:   return(0);
 75: }

 79: /*
 80:   QEPAbsoluteConverged - Checks convergence absolutely.
 81: */
 82: PetscErrorCode QEPAbsoluteConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 83: {
 85:   *errest = res;
 86:   return(0);
 87: }

 91: PetscErrorCode QEPComputeVectors_Schur(QEP qep)
 92: {
 93: #if defined(SLEPC_MISSING_LAPACK_TREVC)
 94:   SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
 95: #else
 97:   PetscInt       i;
 98:   PetscBLASInt   ncv,nconv,mout,info,one = 1;
 99:   PetscScalar    *Z,*work,tmp;
100: #if defined(PETSC_USE_COMPLEX)
101:   PetscReal      *rwork;
102: #else 
103:   PetscReal      normi;
104: #endif
105:   PetscReal      norm;
106: 
108:   ncv = PetscBLASIntCast(qep->ncv);
109:   nconv = PetscBLASIntCast(qep->nconv);
110:   PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
111:   PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
112: #if defined(PETSC_USE_COMPLEX)
113:   PetscMalloc(nconv*sizeof(PetscReal),&rwork);
114: #endif

116:   /* right eigenvectors */
117: #if !defined(PETSC_USE_COMPLEX)
118:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,qep->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
119: #else
120:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,qep->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
121: #endif
122:   if (info) SETERRQ1(((PetscObject)qep)->comm,PETSC_ERR_LIB,"Error in Lapack xTREVC %d",info);

124:   /* normalize eigenvectors */
125:   for (i=0;i<qep->nconv;i++) {
126: #if !defined(PETSC_USE_COMPLEX)
127:     if (qep->eigi[i] != 0.0) {
128:       norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
129:       normi = BLASnrm2_(&nconv,Z+(i+1)*nconv,&one);
130:       tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
131:       BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
132:       BLASscal_(&nconv,&tmp,Z+(i+1)*nconv,&one);
133:       i++;
134:     } else
135: #endif
136:     {
137:       norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
138:       tmp = 1.0 / norm;
139:       BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
140:     }
141:   }
142: 
143:   /* AV = V * Z */
144:   SlepcUpdateVectors(qep->nconv,qep->V,0,qep->nconv,Z,qep->nconv,PETSC_FALSE);
145: 
146:   PetscFree(Z);
147:   PetscFree(work);
148: #if defined(PETSC_USE_COMPLEX)
149:   PetscFree(rwork);
150: #endif
151:    return(0);
152: #endif 
153: }

157: /*
158:    QEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
159:    for quadratic Krylov methods.

161:    Differences:
162:    - Always non-symmetric
163:    - Does not check for STSHIFT
164:    - No correction factor
165:    - No support for true residual
166: */
167: PetscErrorCode QEPKrylovConvergence(QEP qep,PetscInt kini,PetscInt nits,PetscScalar *S,PetscInt lds,PetscScalar *Q,PetscInt nv,PetscReal beta,PetscInt *kout,PetscScalar *work)
168: {
170:   PetscInt       k,marker;
171:   PetscScalar    re,im,*Z,*work2;
172:   PetscReal      resnorm;
173:   PetscBool      iscomplex;

176:   Z = work; work2 = work+2*nv;
177:   marker = -1;
178:   for (k=kini;k<kini+nits;k++) {
179:     /* eigenvalue */
180:     re = qep->eigr[k];
181:     im = qep->eigi[k];
182:     iscomplex = PETSC_FALSE;
183:     if (k<nv-1 && S[k+1+k*lds] != 0.0) iscomplex = PETSC_TRUE;
184:     /* residual norm */
185:     DenseSelectedEvec(S,lds,Q,Z,k,iscomplex,nv,work2);
186:     if (iscomplex) resnorm = beta*SlepcAbsEigenvalue(Z[nv-1],Z[2*nv-1]);
187:     else resnorm = beta*PetscAbsScalar(Z[nv-1]);
188:     /* error estimate */
189:     (*qep->conv_func)(qep,re,im,resnorm,&qep->errest[k],qep->conv_ctx);
190:     if (marker==-1 && qep->errest[k] >= qep->tol) marker = k;
191:     if (iscomplex) { qep->errest[k+1] = qep->errest[k]; k++; }
192:     if (marker!=-1 && !qep->trackall) break;
193:   }
194:   if (marker!=-1) k = marker;
195:   *kout = k;
196:   return(0);
197: }