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-2011, 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 <private/epsimpl.h>   /*I "slepceps.h" I*/
 25: #include <slepcblaslapack.h>

 29: PetscErrorCode EPSReset_Default(EPS eps)
 30: {

 34:   EPSDefaultFreeWork(eps);
 35:   EPSFreeSolution(eps);
 36:   return(0);
 37: }

 41: PetscErrorCode EPSBackTransform_Default(EPS eps)
 42: {

 46:   STBackTransform(eps->OP,eps->nconv,eps->eigr,eps->eigi);
 47:   return(0);
 48: }

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

 66: /*
 67:   EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
 68:   using purification for generalized eigenproblems.
 69:  */
 70: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
 71: {
 73:   PetscInt       i;
 74:   PetscReal      norm;
 75:   Vec            w;

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

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

129:   PetscMalloc(nconv*nconv*sizeof(PetscScalar),&Z);
130:   PetscMalloc(3*nconv*sizeof(PetscScalar),&work);
131: #if defined(PETSC_USE_COMPLEX)
132:   PetscMalloc(nconv*sizeof(PetscReal),&rwork);
133: #endif

135:   /* right eigenvectors */
136: #if !defined(PETSC_USE_COMPLEX)
137:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,&info);
138: #else
139:   LAPACKtrevc_("R","A",PETSC_NULL,&nconv,eps->T,&ncv,PETSC_NULL,&nconv,Z,&nconv,&nconv,&mout,work,rwork,&info);
140: #endif
141:   if (info) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

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

167:   /* Purify eigenvectors */
168:   if (eps->ispositive) {
169:     VecDuplicate(eps->V[0],&w);
170:     for (i=0;i<eps->nconv;i++) {
171:       VecCopy(eps->V[i],w);
172:       STApply(eps->OP,w,eps->V[i]);
173:     }
174:     VecDestroy(&w);
175:   }

177:   /* Fix eigenvectors if balancing was used */
178:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
179:     for (i=0;i<eps->nconv;i++) {
180:       VecPointwiseDivide(eps->V[i],eps->V[i],eps->D);
181:     }
182:   }

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

212:     /* AW = W * Z */
213:     SlepcUpdateVectors(eps->nconv,eps->W,0,eps->nconv,Z,eps->nconv,PETSC_FALSE);
214:   }
215: 
216:   PetscFree(Z);
217:   PetscFree(work);
218: #if defined(PETSC_USE_COMPLEX)
219:   PetscFree(rwork);
220: #endif
221:   eps->evecsavailable = PETSC_TRUE;
222:   return(0);
223: #endif 
224: }

228: /*
229:   EPSDefaultGetWork - Gets a number of work vectors.
230:  */
231: PetscErrorCode EPSDefaultGetWork(EPS eps,PetscInt nw)
232: {

236:   if (eps->nwork != nw) {
237:     VecDestroyVecs(eps->nwork,&eps->work);
238:     eps->nwork = nw;
239:     VecDuplicateVecs(eps->t,nw,&eps->work);
240:     PetscLogObjectParents(eps,nw,eps->work);
241:   }
242:   return(0);
243: }

247: /*
248:   EPSDefaultFreeWork - Free work vectors.
249:  */
250: PetscErrorCode EPSDefaultFreeWork(EPS eps)
251: {

255:   VecDestroyVecs(eps->nwork,&eps->work);
256:   eps->nwork = 0;
257:   return(0);
258: }

262: /*
263:   EPSDefaultSetWhich - Sets the default value for which, depending on the ST.
264:  */
265: PetscErrorCode EPSDefaultSetWhich(EPS eps)
266: {
267:   PetscBool      target;

271:   PetscTypeCompareAny((PetscObject)eps->OP,&target,STSINVERT,STCAYLEY,STFOLD,"");
272:   if (target) eps->which = EPS_TARGET_MAGNITUDE;
273:   else eps->which = EPS_LARGEST_MAGNITUDE;
274:   return(0);
275: }

279: /*
280:   EPSDefaultConverged - Checks convergence relative to the eigenvalue.
281: */
282: PetscErrorCode EPSEigRelativeConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
283: {
284:   PetscReal w;

287:   w = SlepcAbsEigenvalue(eigr,eigi);
288:   *errest = res;
289:   if (w > res) *errest = res / w;
290:   return(0);
291: }

295: /*
296:   EPSAbsoluteConverged - Checks convergence absolutely.
297: */
298: PetscErrorCode EPSAbsoluteConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
299: {
301:   *errest = res;
302:   return(0);
303: }

307: /*
308:   EPSNormRelativeConverged - Checks convergence relative to the eigenvalue and 
309:   the matrix norms.
310: */
311: PetscErrorCode EPSNormRelativeConverged(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
312: {
313:   PetscReal w;

316:   w = SlepcAbsEigenvalue(eigr,eigi);
317:   *errest = res / (eps->nrma + w*eps->nrmb);
318:   return(0);
319: }

323: /*
324:   EPSComputeTrueResidual - Computes the true residual norm of a given Ritz pair:
325:     ||r|| = ||A*x - lambda*B*x||
326:   where lambda is the Ritz value and x is the Ritz vector.

328:   Real lambda:
329:     lambda = eigr
330:     x = V*Z  (V is an array of nv vectors, Z has length nv)

332:   Complex lambda:
333:     lambda = eigr+i*eigi
334:     x = V*Z[0*nv]+i*V*Z[1*nv]  (Z has length 2*nv)
335: */
336: PetscErrorCode EPSComputeTrueResidual(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscScalar *Z,Vec *V,PetscInt nv,PetscReal *resnorm)
337: {
339:   Vec            x,y,z=0;
340:   PetscReal      norm;
341: 
343:   /* allocate workspace */
344:   VecDuplicate(V[0],&x);
345:   VecDuplicate(V[0],&y);
346:   if (!eps->ishermitian && eps->ispositive) { VecDuplicate(V[0],&z); }

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

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

386:   /* free workspace */
387:   VecDestroy(&x);
388:   VecDestroy(&y);
389:   VecDestroy(&z);
390:   return(0);
391: }

395: /*
396:   EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
397:   diagonal matrix to be applied for balancing in non-Hermitian problems.
398: */
399: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
400: {
401:   Vec               z,p,r;
402:   PetscInt          i,j;
403:   PetscReal         norma;
404:   PetscScalar       *pz,*pD;
405:   const PetscScalar *pr,*pp;
406:   PetscErrorCode    ierr;

409:   VecDuplicate(eps->V[0],&r);
410:   VecDuplicate(eps->V[0],&p);
411:   VecDuplicate(eps->V[0],&z);
412:   VecSet(eps->D,1.0);

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

416:     /* Build a random vector of +-1's */
417:     SlepcVecSetRandom(z,eps->rand);
418:     VecGetArray(z,&pz);
419:     for (i=0;i<eps->nloc;i++) {
420:       if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
421:       else pz[i]=1.0;
422:     }
423:     VecRestoreArray(z,&pz);

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

458:   VecDestroy(&r);
459:   VecDestroy(&p);
460:   VecDestroy(&z);
461:   return(0);
462: }