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-2010, Universidad 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
 25:  #include private/slepcimpl.h
 26:  #include slepcblaslapack.h

 30: PetscErrorCode QEPDestroy_Default(QEP qep)
 31: {

 36:   PetscFree(qep->data);

 38:   /* free work vectors */
 39:   QEPDefaultFreeWork(qep);
 40:   return(0);
 41: }

 45: /*
 46:   QEPDefaultGetWork - Gets a number of work vectors.
 47:  */
 48: PetscErrorCode QEPDefaultGetWork(QEP qep, PetscInt nw)
 49: {
 51:   PetscInt       i;


 55:   if (qep->nwork != nw) {
 56:     if (qep->nwork > 0) {
 57:       VecDestroyVecs(qep->work,qep->nwork);
 58:     }
 59:     qep->nwork = nw;
 60:     PetscMalloc(nw*sizeof(Vec),&qep->work);
 61:     for (i=0;i<nw;i++) {
 62:       MatGetVecs(qep->M,PETSC_NULL,qep->work+i);
 63:     }
 64:     PetscLogObjectParents(qep,nw,qep->work);
 65:   }
 66: 
 67:   return(0);
 68: }

 72: /*
 73:   QEPDefaultFreeWork - Free work vectors.
 74:  */
 75: PetscErrorCode QEPDefaultFreeWork(QEP qep)
 76: {

 81:   if (qep->work)  {
 82:     VecDestroyVecs(qep->work,qep->nwork);
 83:   }
 84:   return(0);
 85: }

 89: /*
 90:   QEPDefaultConverged - Checks convergence relative to the eigenvalue.
 91: */
 92: PetscErrorCode QEPDefaultConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 93: {
 94:   PetscReal w;
 96:   w = SlepcAbsEigenvalue(eigr,eigi);
 97:   *errest = res;
 98:   if (w > res) *errest = res / w;
 99:   return(0);
100: }

104: /*
105:   QEPAbsoluteConverged - Checks convergence absolutely.
106: */
107: PetscErrorCode QEPAbsoluteConverged(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
108: {
110:   *errest = res;
111:   return(0);
112: }

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
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;
131: 
133:   ncv = PetscBLASIntCast(qep->ncv);
134:   nconv = PetscBLASIntCast(qep->nconv);
135: 
136:   PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
137:   PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
138: #if defined(PETSC_USE_COMPLEX)
139:   PetscMalloc(nconv*sizeof(PetscReal),&rwork);
140: #endif

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);

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:     }
167:   }
168: 
169:   /* AV = V * Z */
170:   SlepcUpdateVectors(qep->nconv,qep->V,0,qep->nconv,Z,qep->nconv,PETSC_FALSE);
171: 
172:   PetscFree(Z);
173:   PetscFree(work);
174: #if defined(PETSC_USE_COMPLEX)
175:   PetscFree(rwork);
176: #endif
177:    return(0);
178: #endif 
179: }

183: /*
184:    QEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
185:    for quadratic Krylov methods.

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: {
196:   PetscInt       k,marker;
197:   PetscScalar    re,im,*Z,*work2;
198:   PetscReal      resnorm;
199:   PetscTruth     iscomplex;

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:     DenseSelectedEvec(S,lds,Q,Z,k,iscomplex,nv,work2);
212:     if (iscomplex) resnorm = beta*SlepcAbsEigenvalue(Z[nv-1],Z[2*nv-1]);
213:     else resnorm = beta*PetscAbsScalar(Z[nv-1]);
214:     /* error estimate */
215:     (*qep->conv_func)(qep,re,im,resnorm,&qep->errest[k],qep->conv_ctx);
216:     if (marker==-1 && qep->errest[k] >= qep->tol) marker = k;
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;

223:   return(0);
224: }