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 <slepcblaslapack.h>

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

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

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

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

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

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

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

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

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

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

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

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

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