Actual source code: power.c

  1: /*                       

  3:    SLEPc eigensolver: "power"

  5:    Method: Power Iteration

  7:    Algorithm:

  9:        This solver implements the power iteration for finding dominant
 10:        eigenpairs. It also includes the following well-known methods:
 11:        - Inverse Iteration: when used in combination with shift-and-invert
 12:          spectral transformation.
 13:        - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
 14:          a variable shift.

 16:    References:

 18:        [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report STR-2, 
 19:            available at http://www.grycap.upv.es/slepc.

 21:    Last update: June 2005

 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24:       SLEPc - Scalable Library for Eigenvalue Problem Computations
 25:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

 27:       This file is part of SLEPc. See the README file for conditions of use
 28:       and additional information.
 29:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 30: */

 32:  #include src/eps/epsimpl.h
 33:  #include slepcblaslapack.h

 35: typedef struct {
 36:   EPSPowerShiftType shift_type;
 37: } EPS_POWER;

 41: PetscErrorCode EPSSetUp_POWER(EPS eps)
 42: {
 44:   EPS_POWER      *power = (EPS_POWER *)eps->data;
 45:   PetscInt       N;
 46:   PetscTruth     flg;
 47:   STMatMode      mode;

 50:   VecGetSize(eps->vec_initial,&N);
 51:   if (eps->ncv) {
 52:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 53:   }
 54:   else eps->ncv = eps->nev;
 55:   if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
 56:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 57:     SETERRQ(1,"Wrong value of eps->which");
 58:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
 59:     PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
 60:     if (!flg)
 61:       SETERRQ(PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert ST");
 62:     STGetMatMode(eps->OP,&mode);
 63:     if (mode == STMATMODE_INPLACE)
 64:       SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 65:   }
 66:   EPSAllocateSolution(eps);
 67:   if (eps->solverclass==EPS_TWO_SIDE) {
 68:     EPSDefaultGetWork(eps,1);
 69:   } else {
 70:     EPSDefaultGetWork(eps,2);
 71:   }
 72:   return(0);
 73: }

 77: PetscErrorCode EPSSolve_POWER(EPS eps)
 78: {
 80:   EPS_POWER      *power = (EPS_POWER *)eps->data;
 81:   int            i, nsv;
 82:   Vec            v, y, e, *SV;
 83:   Mat            A;
 84:   PetscReal      relerr, norm, rt1, rt2, cs1, anorm;
 85:   PetscScalar    theta, rho, delta, sigma, alpha2, beta1, sn1;
 86:   PetscTruth     breakdown;

 89:   v = eps->V[0];
 90:   y = eps->AV[0];
 91:   e = eps->work[0];

 93:   /* prepare for selective orthogonalization of converged vectors */
 94:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
 95:     PetscMalloc(eps->nev*sizeof(Vec),&SV);
 96:     for (i=0;i<eps->nds;i++) SV[i]=eps->DS[i];
 97:     if (eps->nev>1) {
 98:       STGetOperators(eps->OP,&A,PETSC_NULL);
 99:       MatNorm(A,NORM_INFINITY,&anorm);
100:     }
101:   }

103:   EPSGetStartVector(eps,0,v,PETSC_NULL);
104:   STGetShift(eps->OP,&sigma);    /* original shift */
105:   rho = sigma;

107:   while (eps->reason == EPS_CONVERGED_ITERATING) {
108:     eps->its = eps->its + 1;

110:     /* y = OP v */
111:     STApply(eps->OP,v,y);

113:     /* theta = (v,y)_B */
114:     IPInnerProduct(eps->ip,v,y,&theta);

116:     if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */

118:       /* approximate eigenvalue is the Rayleigh quotient */
119:       eps->eigr[eps->nconv] = theta;

121:       /* compute relative error as ||y-theta v||_2/|theta| */
122:       VecCopy(y,e);
123:       VecAXPY(e,-theta,v);
124:       VecNorm(e,NORM_2,&norm);
125:       relerr = norm / PetscAbsScalar(theta);

127:     } else {  /* RQI */

129:       /* delta = ||y||_B */
130:       IPNorm(eps->ip,y,&norm);
131:       delta = norm;

133:       /* compute relative error */
134:       if (rho == 0.0) relerr = PETSC_MAX;
135:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

137:       /* approximate eigenvalue is the shift */
138:       eps->eigr[eps->nconv] = rho;

140:       /* compute new shift */
141:       if (relerr<eps->tol) {
142:         rho = sigma; /* if converged, restore original shift */
143:         STSetShift(eps->OP,rho);
144:       } else {
145:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v) */
146:         if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
147: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
148:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
149: #else 
150:           /* beta1 is the norm of the residual associated to R(v) */
151:           VecAXPY(v,-theta/(delta*delta),y);
152:           VecScale(v,1.0/delta);
153:           IPNorm(eps->ip,v,&norm);
154:           beta1 = norm;
155: 
156:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
157:           STGetOperators(eps->OP,&A,PETSC_NULL);
158:           MatMult(A,v,e);
159:           VecDot(v,e,&alpha2);
160:           alpha2 = alpha2 / (beta1 * beta1);

162:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
163:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
164:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
165:           else rho = rt2;
166: #endif 
167:         }
168:         /* update operator according to new shift */
169:         PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
170:         STSetShift(eps->OP,rho);
171:         PetscPopErrorHandler();
172:         if (ierr) {
173:           eps->eigr[eps->nconv] = rho;
174:           relerr = PETSC_MACHINE_EPSILON;
175:           rho = sigma;
176:           STSetShift(eps->OP,rho);
177:         }
178:       }
179:     }

181:     eps->errest[eps->nconv] = relerr;
182:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);

184:     /* purge previously converged eigenvectors */
185:     if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
186:       nsv = eps->nds;
187:       for (i=0;i<eps->nconv;i++) {
188:         if(PetscAbsScalar(rho-eps->eigr[i])>eps->its*anorm/1000) SV[nsv++]=eps->V[i];
189:       }
190:       IPOrthogonalize(eps->ip,nsv,PETSC_NULL,SV,y,PETSC_NULL,&norm,PETSC_NULL,eps->work[1]);
191:     } else {
192:       IPOrthogonalize(eps->ip,eps->nds+eps->nconv,PETSC_NULL,eps->DSV,y,PETSC_NULL,&norm,PETSC_NULL,eps->work[1]);
193:     }

195:     /* v = y/||y||_B */
196:     VecCopy(y,v);
197:     VecScale(v,1.0/norm);

199:     /* if relerr<tol, accept eigenpair */
200:     if (relerr<eps->tol) {
201:       eps->nconv = eps->nconv + 1;
202:       if (eps->nconv==eps->nev) eps->reason = EPS_CONVERGED_TOL;
203:       else {
204:         v = eps->V[eps->nconv];
205:         EPSGetStartVector(eps,eps->nconv,v,&breakdown);
206:         if (breakdown) {
207:           eps->reason = EPS_DIVERGED_BREAKDOWN;
208:           PetscInfo(eps,"Unable to generate more start vectors\n");
209:         }
210:       }
211:     }

213:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
214:   }

216:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
217:     PetscFree(SV);
218:   }

220:   return(0);
221: }

225: PetscErrorCode EPSSolve_TS_POWER(EPS eps)
226: {
228:   EPS_POWER      *power = (EPS_POWER *)eps->data;
229:   Vec            v, w, y, z, e;
230:   Mat            A;
231:   PetscReal      relerr, norm, rt1, rt2, cs1;
232:   PetscScalar    theta, alpha, beta, rho, delta, sigma, alpha2, beta1, sn1;

235:   v = eps->V[0];
236:   y = eps->AV[0];
237:   e = eps->work[0];
238:   w = eps->W[0];
239:   z = eps->AW[0];

241:   EPSGetStartVector(eps,0,v,PETSC_NULL);
242:   EPSGetLeftStartVector(eps,0,w);
243:   STGetShift(eps->OP,&sigma);    /* original shift */
244:   rho = sigma;

246:   while (eps->its<eps->max_it) {
247:     eps->its++;
248: 
249:     /* y = OP v, z = OP' w */
250:     STApply(eps->OP,v,y);
251:     STApplyTranspose(eps->OP,w,z);

253:     /* theta = (v,z)_B */
254:     IPInnerProduct(eps->ip,v,z,&theta);

256:     if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */

258:       /* approximate eigenvalue is the Rayleigh quotient */
259:       eps->eigr[eps->nconv] = theta;

261:       /* compute relative errors (right and left) */
262:       VecCopy(y,e);
263:       VecAXPY(e,-theta,v);
264:       VecNorm(e,NORM_2,&norm);
265:       relerr = norm / PetscAbsScalar(theta);
266:       eps->errest[eps->nconv] = relerr;
267:       VecCopy(z,e);
268:       VecAXPY(e,-theta,w);
269:       VecNorm(e,NORM_2,&norm);
270:       relerr = norm / PetscAbsScalar(theta);
271:       eps->errest_left[eps->nconv] = relerr;

273:     } else {  /* RQI */

275:       /* delta = sqrt(y,z)_B */
276:       IPInnerProduct(eps->ip,y,z,&alpha);
277:       if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
278:       delta = PetscSqrtScalar(alpha);

280:       /* compute relative error */
281:       if (rho == 0.0) relerr = PETSC_MAX;
282:       else relerr = 1.0 / (PetscAbsScalar(delta*rho));
283:       eps->errest[eps->nconv] = relerr;
284:       eps->errest_left[eps->nconv] = relerr;

286:       /* approximate eigenvalue is the shift */
287:       eps->eigr[eps->nconv] = rho;

289:       /* compute new shift */
290:       if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
291:         rho = sigma; /* if converged, restore original shift */
292:         STSetShift(eps->OP,rho);
293:       } else {
294:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v,w) */
295:         if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
296: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
297:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
298: #else 
299:           /* beta1 is the norm of the residual associated to R(v,w) */
300:           VecAXPY(v,-theta/(delta*delta),y);
301:           VecScale(v,1.0/delta);
302:           IPNorm(eps->ip,v,&norm);
303:           beta1 = norm;
304: 
305:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
306:           STGetOperators(eps->OP,&A,PETSC_NULL);
307:           MatMult(A,v,e);
308:           VecDot(v,e,&alpha2);
309:           alpha2 = alpha2 / (beta1 * beta1);

311:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
312:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
313:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
314:           else rho = rt2;
315: #endif 
316:         }
317:         /* update operator according to new shift */
318:         PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
319:         STSetShift(eps->OP,rho);
320:         PetscPopErrorHandler();
321:         if (ierr) {
322:           eps->eigr[eps->nconv] = rho;
323:           eps->errest[eps->nconv] = PETSC_MACHINE_EPSILON;
324:           eps->errest_left[eps->nconv] = PETSC_MACHINE_EPSILON;
325:           rho = sigma;
326:           STSetShift(eps->OP,rho);
327:         }
328:       }
329:     }

331:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
332:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,eps->nconv+1);

334:     /* purge previously converged eigenvectors */
335:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
336:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);

338:     /* normalize so that (y,z)_B=1  */
339:     VecCopy(y,v);
340:     VecCopy(z,w);
341:     IPInnerProduct(eps->ip,y,z,&alpha);
342:     if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
343:     delta = PetscSqrtScalar(PetscAbsScalar(alpha));
344:     beta = 1.0/PetscConj(alpha/delta);
345:     delta = 1.0/delta;
346:     VecScale(w,beta);
347:     VecScale(v,delta);

349:     /* if relerr<tol (both right and left), accept eigenpair */
350:     if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
351:       eps->nconv = eps->nconv + 1;
352:       if (eps->nconv==eps->nev) break;
353:       v = eps->V[eps->nconv];
354:       EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
355:       w = eps->W[eps->nconv];
356:       EPSGetLeftStartVector(eps,eps->nconv,w);
357:     }
358:   }

360:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
361:   else eps->reason = EPS_DIVERGED_ITS;

363:   return(0);
364: }

368: PetscErrorCode EPSBackTransform_POWER(EPS eps)
369: {
371:   EPS_POWER *power = (EPS_POWER *)eps->data;

374:   if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) {
375:     EPSBackTransform_Default(eps);
376:   }
377:   return(0);
378: }

382: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
383: {
385:   EPS_POWER      *power = (EPS_POWER *)eps->data;
386:   PetscTruth     flg;
387:   PetscInt       i;
388:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

391:   PetscOptionsHead("POWER options");
392:   PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],&i,&flg);
393:   if (flg ) power->shift_type = (EPSPowerShiftType)i;
394:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
395:     STSetType(eps->OP,STSINV);
396:   }
397:   PetscOptionsTail();
398:   return(0);
399: }

404: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
405: {
406:   EPS_POWER *power = (EPS_POWER *)eps->data;

409:   switch (shift) {
410:     case EPSPOWER_SHIFT_CONSTANT:
411:     case EPSPOWER_SHIFT_RAYLEIGH:
412:     case EPSPOWER_SHIFT_WILKINSON:
413:       power->shift_type = shift;
414:       break;
415:     default:
416:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
417:   }
418:   return(0);
419: }

424: /*@
425:    EPSPowerSetShiftType - Sets the type of shifts used during the power
426:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
427:    (RQI) method.

429:    Collective on EPS

431:    Input Parameters:
432: +  eps - the eigenproblem solver context
433: -  shift - the type of shift

435:    Options Database Key:
436: .  -eps_power_shift_type - Sets the shift type (either 'constant' or 
437:                            'rayleigh' or 'wilkinson')

439:    Notes:
440:    By default, shifts are constant (EPSPOWER_SHIFT_CONSTANT) and the iteration
441:    is the simple power method (or inverse iteration if a shift-and-invert
442:    transformation is being used).

444:    A variable shift can be specified (EPSPOWER_SHIFT_RAYLEIGH or
445:    EPSPOWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
446:    a cubic converging method as RQI. See the users manual for details.
447:    
448:    Level: advanced

450: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
451: @*/
452: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
453: {
454:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);

458:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
459:   if (f) {
460:     (*f)(eps,shift);
461:   }
462:   return(0);
463: }

468: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
469: {
470:   EPS_POWER  *power = (EPS_POWER *)eps->data;
472:   *shift = power->shift_type;
473:   return(0);
474: }

479: /*@C
480:    EPSPowerGetShiftType - Gets the type of shifts used during the power
481:    iteration. 

483:    Collective on EPS

485:    Input Parameter:
486: .  eps - the eigenproblem solver context

488:    Input Parameter:
489: .  shift - the type of shift

491:    Level: advanced

493: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
494: @*/
495: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
496: {
497:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);

501:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
502:   if (f) {
503:     (*f)(eps,shift);
504:   }
505:   return(0);
506: }

510: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
511: {
513:   EPS_POWER      *power = (EPS_POWER *)eps->data;
514:   PetscTruth     isascii;
515:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

518:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
519:   if (!isascii) {
520:     SETERRQ1(1,"Viewer type %s not supported for EPSPOWER",((PetscObject)viewer)->type_name);
521:   }
522:   PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
523:   return(0);
524: }

529: PetscErrorCode EPSCreate_POWER(EPS eps)
530: {
532:   EPS_POWER      *power;

535:   PetscNew(EPS_POWER,&power);
536:   PetscLogObjectMemory(eps,sizeof(EPS_POWER));
537:   eps->data                      = (void *) power;
538:   eps->ops->solve                = EPSSolve_POWER;
539:   eps->ops->solvets              = EPSSolve_TS_POWER;
540:   eps->ops->setup                = EPSSetUp_POWER;
541:   eps->ops->setfromoptions       = EPSSetFromOptions_POWER;
542:   eps->ops->destroy              = EPSDestroy_Default;
543:   eps->ops->view                 = EPSView_POWER;
544:   eps->ops->backtransform        = EPSBackTransform_POWER;
545:   eps->ops->computevectors       = EPSComputeVectors_Default;
546:   power->shift_type              = EPSPOWER_SHIFT_CONSTANT;
547:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
548:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
549:   return(0);
550: }