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-2012, 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 <slepc-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:   QEPConvergedDefault - Checks convergence relative to the eigenvalue.
 64: */
 65: PetscErrorCode QEPConvergedDefault(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 66: {
 67:   PetscReal w;

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

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

 89: PetscErrorCode QEPComputeVectors_Schur(QEP qep)
 90: {
 92:   PetscInt       n,ld;
 93:   PetscScalar    *Z;
 94: 
 96:   DSGetLeadingDimension(qep->ds,&ld);
 97:   DSGetDimensions(qep->ds,&n,PETSC_NULL,PETSC_NULL,PETSC_NULL);

 99:   /* right eigenvectors */
100:   DSVectors(qep->ds,DS_MAT_X,PETSC_NULL,PETSC_NULL);

102:   /* AV = V * Z */
103:   DSGetArray(qep->ds,DS_MAT_X,&Z);
104:   SlepcUpdateVectors(n,qep->V,0,n,Z,ld,PETSC_FALSE);
105:   DSRestoreArray(qep->ds,DS_MAT_X,&Z);
106:   return(0);
107: }

111: /*
112:    QEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
113:    for quadratic Krylov methods.

115:    Differences:
116:    - Always non-symmetric
117:    - Does not check for STSHIFT
118:    - No correction factor
119:    - No support for true residual
120: */
121: PetscErrorCode QEPKrylovConvergence(QEP qep,PetscBool getall,PetscInt kini,PetscInt nits,PetscInt nv,PetscReal beta,PetscInt *kout)
122: {
124:   PetscInt       k,newk,marker,ld;
125:   PetscScalar    re,im;
126:   PetscReal      resnorm;

129:   DSGetLeadingDimension(qep->ds,&ld);
130:   marker = -1;
131:   if (qep->trackall) getall = PETSC_TRUE;
132:   for (k=kini;k<kini+nits;k++) {
133:     /* eigenvalue */
134:     re = qep->eigr[k];
135:     im = qep->eigi[k];
136:     newk = k;
137:     DSVectors(qep->ds,DS_MAT_X,&newk,&resnorm);
138:     resnorm *= beta;
139:     /* error estimate */
140:     (*qep->conv_func)(qep,re,im,resnorm,&qep->errest[k],qep->conv_ctx);
141:     if (marker==-1 && qep->errest[k] >= qep->tol) marker = k;
142:     if (newk==k+1) { qep->errest[k+1] = qep->errest[k]; k++; }
143:     if (marker!=-1 && !getall) break;
144:   }
145:   if (marker!=-1) k = marker;
146:   *kout = k;
147:   return(0);
148: }