Actual source code: default.c

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

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2009, 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/epsimpl.h
 25:  #include slepcblaslapack.h

 29: PetscErrorCode EPSDestroy_Default(EPS eps)
 30: {

 35:   PetscFree(eps->data);

 37:   /* free work vectors */
 38:   EPSDefaultFreeWork(eps);
 39:   EPSFreeSolution(eps);
 40:   return(0);
 41: }

 45: PetscErrorCode EPSBackTransform_Default(EPS eps)
 46: {
 48:   PetscInt       i;

 52:   for (i=0;i<eps->nconv;i++) {
 53:     STBackTransform(eps->OP,&eps->eigr[i],&eps->eigi[i]);
 54:   }
 55:   return(0);
 56: }

 60: /*
 61:   EPSComputeVectors_Default - Compute eigenvectors from the vectors
 62:   provided by the eigensolver. This version just copies the vectors
 63:   and is intended for solvers such as power that provide the eigenvector.
 64:  */
 65: PetscErrorCode EPSComputeVectors_Default(EPS eps)
 66: {
 68:   eps->evecsavailable = PETSC_TRUE;
 69:   return(0);
 70: }

 74: /*
 75:   EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
 76:   using purification for generalized eigenproblems.
 77:  */
 78: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
 79: {
 81:   PetscInt       i;
 82:   PetscReal      norm;
 83:   Vec            w;

 86:   if (eps->isgeneralized) {
 87:     /* Purify eigenvectors */
 88:     VecDuplicate(eps->V[0],&w);
 89:     for (i=0;i<eps->nconv;i++) {
 90:       VecCopy(eps->V[i],w);
 91:       STApply(eps->OP,w,eps->V[i]);
 92:       IPNorm(eps->ip,eps->V[i],&norm);
 93:       VecScale(eps->V[i],1.0/norm);
 94:     }
 95:     VecDestroy(w);
 96:   }
 97:   eps->evecsavailable = PETSC_TRUE;
 98:   return(0);
 99: }

103: /*
104:   EPSComputeVectors_Schur - Compute eigenvectors from the vectors
105:   provided by the eigensolver. This version is intended for solvers 
106:   that provide Schur vectors. Given the partial Schur decomposition
107:   OP*V=V*T, the following steps are performed:
108:       1) compute eigenvectors of T: T*Z=Z*D
109:       2) compute eigenvectors of OP: X=V*Z
110:   If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
111:  */
112: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
113: {
114: #if defined(SLEPC_MISSING_LAPACK_TREVC)
115:   SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
116: #else
118:   PetscInt       i;
119:   PetscBLASInt   ncv,nconv,mout,info;
120:   PetscScalar    *Z,*work;
121: #if defined(PETSC_USE_COMPLEX)
122:   PetscReal      *rwork;
123: #endif
124:   PetscReal      norm;
125:   Vec            w;
126: 
128:   ncv = PetscBLASIntCast(eps->ncv);
129:   nconv = PetscBLASIntCast(eps->nconv);
130:   if (eps->ishermitian) {
131:     EPSComputeVectors_Hermitian(eps);
132:     return(0);
133:   }
134:   if (eps->ispositive) {
135:     VecDuplicate(eps->V[0],&w);
136:   }

138:   PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
139:   PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
140: #if defined(PETSC_USE_COMPLEX)
141:   PetscMalloc(nconv*sizeof(PetscReal),&rwork);
142: #endif

144:   /* right eigenvectors */
145: #if !defined(PETSC_USE_COMPLEX)
146:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
147: #else
148:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
149: #endif
150:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

152:   /* AV = V * Z */
153:   SlepcUpdateVectors(eps->nconv,eps->V,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
154:   if (eps->ispositive) {
155:     /* Purify eigenvectors */
156:     for (i=0;i<eps->nconv;i++) {
157:       VecCopy(eps->V[i],w);
158:       STApply(eps->OP,w,eps->V[i]);
159:       VecNormalize(eps->V[i],&norm);
160:     }
161:   }
162: 
163:   /* left eigenvectors */
164:   if (eps->solverclass == EPS_TWO_SIDE) {
165: #if !defined(PETSC_USE_COMPLEX)
166:     LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
167: #else
168:     LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
169: #endif
170:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

172:     /* AW = W * Z */
173:     SlepcUpdateVectors(eps->nconv,eps->W,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
174:   }
175: 
176:   PetscFree(Z);
177:   PetscFree(work);
178: #if defined(PETSC_USE_COMPLEX)
179:   PetscFree(rwork);
180: #endif
181:   if (eps->ispositive) {
182:     VecDestroy(w);
183:   }
184:   eps->evecsavailable = PETSC_TRUE;
185:   return(0);
186: #endif 
187: }

191: /*
192:   EPSDefaultGetWork - Gets a number of work vectors.

194:   Input Parameters:
195: + eps  - eigensolver context
196: - nw   - number of work vectors to allocate

198:   Notes:
199:   Call this only if no work vectors have been allocated.

201:  */
202: PetscErrorCode EPSDefaultGetWork(EPS eps, PetscInt nw)
203: {


208:   if (eps->nwork != nw) {
209:     if (eps->nwork > 0) {
210:       VecDestroyVecs(eps->work,eps->nwork);
211:     }
212:     eps->nwork = nw;
213:     VecDuplicateVecs(eps->vec_initial,nw,&eps->work);
214:     PetscLogObjectParents(eps,nw,eps->work);
215:   }
216: 
217:   return(0);
218: }

222: /*
223:   EPSDefaultFreeWork - Free work vectors.

225:   Input Parameters:
226: . eps  - eigensolver context

228:  */
229: PetscErrorCode EPSDefaultFreeWork(EPS eps)
230: {

235:   if (eps->work)  {
236:     VecDestroyVecs(eps->work,eps->nwork);
237:   }
238:   return(0);
239: }