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-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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: PetscErrorCode QEPReset_Default(QEP qep)
 30: {

 34:   VecDestroyVecs(qep->nwork,&qep->work);
 35:   qep->nwork = 0;
 36:   QEPFreeSolution(qep);
 37:   return(0);
 38: }

 42: /*@
 43:    QEPSetWorkVecs - Sets a number of work vectors into a QEP object

 45:    Collective on QEP

 47:    Input Parameters:
 48: +  qep - quadratic eigensolver context
 49: -  nw  - number of work vectors to allocate

 51:    Developers Note:
 52:    This is PETSC_EXTERN because it may be required by user plugin QEP
 53:    implementations.

 55:    Level: developer
 56: @*/
 57: PetscErrorCode QEPSetWorkVecs(QEP qep,PetscInt nw)
 58: {

 62:   if (qep->nwork != nw) {
 63:     VecDestroyVecs(qep->nwork,&qep->work);
 64:     qep->nwork = nw;
 65:     VecDuplicateVecs(qep->t,nw,&qep->work);
 66:     PetscLogObjectParents(qep,nw,qep->work);
 67:   }
 68:   return(0);
 69: }

 73: /*
 74:   QEPConvergedDefault - Checks convergence relative to the eigenvalue.
 75: */
 76: PetscErrorCode QEPConvergedDefault(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 77: {
 78:   PetscReal w;

 81:   w = SlepcAbsEigenvalue(eigr,eigi);
 82:   *errest = res/w;
 83:   return(0);
 84: }

 88: /*
 89:   QEPConvergedAbsolute - Checks convergence absolutely.
 90: */
 91: PetscErrorCode QEPConvergedAbsolute(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 92: {
 94:   *errest = res;
 95:   return(0);
 96: }

100: PetscErrorCode QEPComputeVectors_Schur(QEP qep)
101: {
103:   PetscInt       n,ld;
104:   PetscScalar    *Z;

107:   DSGetLeadingDimension(qep->ds,&ld);
108:   DSGetDimensions(qep->ds,&n,NULL,NULL,NULL,NULL);

110:   /* right eigenvectors */
111:   DSVectors(qep->ds,DS_MAT_X,NULL,NULL);

113:   /* AV = V * Z */
114:   DSGetArray(qep->ds,DS_MAT_X,&Z);
115:   SlepcUpdateVectors(n,qep->V,0,n,Z,ld,PETSC_FALSE);
116:   DSRestoreArray(qep->ds,DS_MAT_X,&Z);
117:   return(0);
118: }

122: /*
123:    QEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
124:    for quadratic Krylov methods.

126:    Differences:
127:    - Always non-symmetric
128:    - Does not check for STSHIFT
129:    - No correction factor
130:    - No support for true residual
131: */
132: PetscErrorCode QEPKrylovConvergence(QEP qep,PetscBool getall,PetscInt kini,PetscInt nits,PetscInt nv,PetscReal beta,PetscInt *kout)
133: {
135:   PetscInt       k,newk,marker,ld;
136:   PetscScalar    re,im;
137:   PetscReal      resnorm;

140:   DSGetLeadingDimension(qep->ds,&ld);
141:   marker = -1;
142:   if (qep->trackall) getall = PETSC_TRUE;
143:   for (k=kini;k<kini+nits;k++) {
144:     /* eigenvalue */
145:     re = qep->eigr[k];
146:     im = qep->eigi[k];
147:     newk = k;
148:     DSVectors(qep->ds,DS_MAT_X,&newk,&resnorm);
149:     resnorm *= beta;
150:     /* error estimate */
151:     (*qep->converged)(qep,re,im,resnorm,&qep->errest[k],qep->convergedctx);
152:     if (marker==-1 && qep->errest[k] >= qep->tol) marker = k;
153:     if (newk==k+1) {
154:       qep->errest[k+1] = qep->errest[k];
155:       k++;
156:     }
157:     if (marker!=-1 && !getall) break;
158:   }
159:   if (marker!=-1) k = marker;
160:   *kout = k;
161:   return(0);
162: }