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-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/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: {

 51:   STBackTransform(eps->OP,eps->nconv,eps->eigr,eps->eigi);
 52:   return(0);
 53: }

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

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

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

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

134:   PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
135:   PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
136: #if defined(PETSC_USE_COMPLEX)
137:   PetscMalloc(nconv*sizeof(PetscReal),&rwork);
138: #endif

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

148:   /* normalize eigenvectors (when not using purification nor balancing)*/
149:   if (!(eps->ispositive || (eps->balance!=EPS_BALANCE_NONE && eps->D))) {
150:     for (i=0;i<eps->nconv;i++) {
151: #if !defined(PETSC_USE_COMPLEX)
152:       if (eps->eigi[i] != 0.0) {
153:         norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
154:         normi = BLASnrm2_(&nconv,Z+(i+1)*nconv,&one);
155:         tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
156:         BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
157:         BLASscal_(&nconv,&tmp,Z+(i+1)*nconv,&one);
158:         i++;
159:       } else
160: #endif
161:       {
162:         norm = BLASnrm2_(&nconv,Z+i*nconv,&one);
163:         tmp = 1.0 / norm;
164:         BLASscal_(&nconv,&tmp,Z+i*nconv,&one);
165:       }
166:     }
167:   }
168: 
169:   /* AV = V * Z */
170:   SlepcUpdateVectors(eps->nconv,eps->V,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);

172:   /* Purify eigenvectors */
173:   if (eps->ispositive) {
174:     VecDuplicate(eps->V[0],&w);
175:     for (i=0;i<eps->nconv;i++) {
176:       VecCopy(eps->V[i],w);
177:       STApply(eps->OP,w,eps->V[i]);
178:     }
179:     VecDestroy(w);
180:   }

182:   /* Fix eigenvectors if balancing was used */
183:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
184:     for (i=0;i<eps->nconv;i++) {
185:       VecPointwiseDivide(eps->V[i],eps->V[i],eps->D);
186:     }
187:   }

189:   /* normalize eigenvectors (when using purification or balancing) */
190:   if (eps->ispositive || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
191:     for (i=0;i<eps->nconv;i++) {
192: #if !defined(PETSC_USE_COMPLEX)
193:       if (eps->eigi[i] != 0.0) {
194:         VecNorm(eps->V[i],NORM_2,&norm);
195:         VecNorm(eps->V[i+1],NORM_2,&normi);
196:         tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
197:         VecScale(eps->V[i],tmp);
198:         VecScale(eps->V[i+1],tmp);
199:         i++;
200:       } else
201: #endif
202:       {
203:         VecNormalize(eps->V[i],PETSC_NULL);
204:       }
205:     }
206:   }
207: 
208:   /* left eigenvectors */
209:   if (eps->leftvecs) {
210: #if !defined(PETSC_USE_COMPLEX)
211:     LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
212: #else
213:     LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->Tl,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
214: #endif
215:     if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

217:     /* AW = W * Z */
218:     SlepcUpdateVectors(eps->nconv,eps->W,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
219:   }
220: 
221:   PetscFree(Z);
222:   PetscFree(work);
223: #if defined(PETSC_USE_COMPLEX)
224:   PetscFree(rwork);
225: #endif
226:   eps->evecsavailable = PETSC_TRUE;
227:   return(0);
228: #endif 
229: }

233: /*
234:   EPSDefaultGetWork - Gets a number of work vectors.
235:  */
236: PetscErrorCode EPSDefaultGetWork(EPS eps, PetscInt nw)
237: {


242:   if (eps->nwork != nw) {
243:     if (eps->nwork > 0) {
244:       VecDestroyVecs(eps->work,eps->nwork);
245:     }
246:     eps->nwork = nw;
247:     VecDuplicateVecs(eps->V[0],nw,&eps->work);
248:     PetscLogObjectParents(eps,nw,eps->work);
249:   }
250: 
251:   return(0);
252: }

256: /*
257:   EPSDefaultFreeWork - Free work vectors.
258:  */
259: PetscErrorCode EPSDefaultFreeWork(EPS eps)
260: {

265:   if (eps->work)  {
266:     VecDestroyVecs(eps->work,eps->nwork);
267:   }
268:   return(0);
269: }

273: /*
274:   EPSDefaultConverged - Checks convergence relative to the eigenvalue.
275: */
276: PetscErrorCode EPSEigRelativeConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
277: {
278:   PetscReal w;
280:   w = SlepcAbsEigenvalue(eigr,eigi);
281:   *errest = res;
282:   if (w > res) *errest = res / w;
283:   return(0);
284: }

288: /*
289:   EPSAbsoluteConverged - Checks convergence absolutely.
290: */
291: PetscErrorCode EPSAbsoluteConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
292: {
294:   *errest = res;
295:   return(0);
296: }

300: /*
301:   EPSNormRelativeConverged - Checks convergence relative to the eigenvalue and 
302:   the matrix norms.
303: */
304: PetscErrorCode EPSNormRelativeConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
305: {
306:   PetscReal w;
308:   w = SlepcAbsEigenvalue(eigr,eigi);
309:   *errest = res / (eps->nrma + w*eps->nrmb);
310:   return(0);
311: }

315: /*
316:   EPSComputeTrueResidual - Computes the true residual norm of a given Ritz pair:
317:     ||r|| = ||A*x - lambda*B*x||
318:   where lambda is the Ritz value and x is the Ritz vector.

320:   Real lambda:
321:     lambda = eigr
322:     x = V*Z  (V is an array of nv vectors, Z has length nv)

324:   Complex lambda:
325:     lambda = eigr+i*eigi
326:     x = V*Z[0*nv]+i*V*Z[1*nv]  (Z has length 2*nv)
327: */
328: PetscErrorCode EPSComputeTrueResidual(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscScalar *Z,Vec *V,PetscInt nv,PetscReal *resnorm)
329: {
331:   Vec            x,y,z;
332:   PetscReal      norm;
333: 
335: 
336:   /* allocate workspace */
337:   VecDuplicate(V[0],&x);
338:   VecDuplicate(V[0],&y);
339:   if (!eps->ishermitian && eps->ispositive) { VecDuplicate(V[0],&z); }

341:   /* compute eigenvector */
342:   SlepcVecMAXPBY(x,0.0,1.0,nv,Z,V);

344:   /* purify eigenvector in positive generalized problems */
345:   if (eps->ispositive) {
346:     STApply(eps->OP,x,y);
347:     if (eps->ishermitian) {
348:       IPNorm(eps->ip,y,&norm);
349:     } else {
350:       VecNorm(y,NORM_2,&norm);
351:     }
352:     VecScale(y,1.0/norm);
353:     VecCopy(y,x);
354:   }
355:   /* fix eigenvector if balancing is used */
356:   if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
357:     VecPointwiseDivide(x,x,eps->D);
358:     VecNormalize(x,&norm);
359:   }
360: #ifndef PETSC_USE_COMPLEX      
361:   /* compute imaginary part of eigenvector */
362:   if (!eps->ishermitian && eigi != 0.0) {
363:     SlepcVecMAXPBY(y,0.0,1.0,nv,Z+nv,V);
364:     if (eps->ispositive) {
365:       STApply(eps->OP,y,z);
366:       VecNorm(z,NORM_2,&norm);
367:       VecScale(z,1.0/norm);
368:       VecCopy(z,y);
369:     }
370:     if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
371:       VecPointwiseDivide(y,y,eps->D);
372:       VecNormalize(y,&norm);
373:     }
374:   }
375: #endif
376:   /* compute relative error and update convergence flag */
377:   EPSComputeResidualNorm_Private(eps,eigr,eigi,x,y,resnorm);

379:   /* free workspace */
380:   VecDestroy(x);
381:   VecDestroy(y);
382:   if (!eps->ishermitian && eps->ispositive) { VecDestroy(z); }
383:   return(0);
384: }

388: /*
389:   EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
390:   diagonal matrix to be applied for balancing in non-Hermitian problems.
391: */
392: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
393: {
394:   Vec            z, p, r;
395:   PetscInt       i, j;
396:   PetscReal      norma;
397:   PetscScalar    *pz, *pr, *pp, *pD;

401:   VecDuplicate(eps->V[0],&r);
402:   VecDuplicate(eps->V[0],&p);
403:   VecDuplicate(eps->V[0],&z);
404:   VecSet(eps->D,1.0);

406:   for (j=0;j<eps->balance_its;j++) {

408:     /* Build a random vector of +-1's */
409:     SlepcVecSetRandom(z,eps->rand);
410:     VecGetArray(z,&pz);
411:     for (i=0;i<eps->nloc;i++) {
412:       if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
413:       else pz[i]=1.0;
414:     }
415:     VecRestoreArray(z,&pz);

417:     /* Compute p=DA(D\z) */
418:     VecPointwiseDivide(r,z,eps->D);
419:     STApply(eps->OP,r,p);
420:     VecPointwiseMult(p,p,eps->D);
421:     if (j==0) {
422:       /* Estimate the matrix inf-norm */
423:       VecAbs(p);
424:       VecMax(p,PETSC_NULL,&norma);
425:     }
426:     if (eps->balance == EPS_BALANCE_TWOSIDE) {
427:       /* Compute r=D\(A'Dz) */
428:       VecPointwiseMult(z,z,eps->D);
429:       STApplyTranspose(eps->OP,z,r);
430:       VecPointwiseDivide(r,r,eps->D);
431:     }
432: 
433:     /* Adjust values of D */
434:     VecGetArray(r,&pr);
435:     VecGetArray(p,&pp);
436:     VecGetArray(eps->D,&pD);
437:     for (i=0;i<eps->nloc;i++) {
438:       if (eps->balance == EPS_BALANCE_TWOSIDE) {
439:         if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
440:           pD[i] *= sqrt(PetscAbsScalar(pr[i]/pp[i]));
441:       } else {
442:         if (pp[i]!=0.0) pD[i] *= 1.0/PetscAbsScalar(pp[i]);
443:       }
444:     }
445:     VecRestoreArray(r,&pr);
446:     VecRestoreArray(p,&pp);
447:     VecRestoreArray(eps->D,&pD);
448:   }

450:   VecDestroy(r);
451:   VecDestroy(p);
452:   VecDestroy(z);
453:   return(0);
454: }