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: Feb 2009

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

 27:    This file is part of SLEPc.
 28:       
 29:    SLEPc is free software: you can redistribute it and/or modify it under  the
 30:    terms of version 3 of the GNU Lesser General Public License as published by
 31:    the Free Software Foundation.

 33:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 34:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 35:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 36:    more details.

 38:    You  should have received a copy of the GNU Lesser General  Public  License
 39:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 40:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: */

 43:  #include private/epsimpl.h
 44:  #include slepcblaslapack.h

 46: PetscErrorCode EPSSolve_POWER(EPS);
 47: PetscErrorCode EPSSolve_TS_POWER(EPS);

 49: typedef struct {
 50:   EPSPowerShiftType shift_type;
 51: } EPS_POWER;

 55: PetscErrorCode EPSSetUp_POWER(EPS eps)
 56: {
 58:   EPS_POWER      *power = (EPS_POWER *)eps->data;
 59:   PetscTruth     flg;
 60:   STMatMode      mode;

 63:   if (eps->ncv) {
 64:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 65:   }
 66:   else eps->ncv = eps->nev;
 67:   if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
 68:   if (!eps->max_it) eps->max_it = PetscMax(2000,100*eps->n);
 69:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
 70:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 71:     SETERRQ(1,"Wrong value of eps->which");
 72:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
 73:     PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&flg);
 74:     if (!flg)
 75:       SETERRQ(PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert ST");
 76:     STGetMatMode(eps->OP,&mode);
 77:     if (mode == ST_MATMODE_INPLACE)
 78:       SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 79:   }
 80:   if (eps->extraction) {
 81:      PetscInfo(eps,"Warning: extraction type ignored\n");
 82:   }
 83:   if (eps->balance!=EPS_BALANCE_NONE)
 84:     SETERRQ(PETSC_ERR_SUP,"Balancing not supported in this solver");
 85:   EPSAllocateSolution(eps);
 86:   if (eps->leftvecs) {
 87:     EPSDefaultGetWork(eps,3);
 88:   } else {
 89:     EPSDefaultGetWork(eps,2);
 90:   }

 92:   /* dispatch solve method */
 93:   if (eps->leftvecs) eps->ops->solve = EPSSolve_TS_POWER;
 94:   else eps->ops->solve = EPSSolve_POWER;
 95:   return(0);
 96: }

100: PetscErrorCode EPSSolve_POWER(EPS eps)
101: {
103:   EPS_POWER      *power = (EPS_POWER *)eps->data;
104:   PetscInt       i;
105:   Vec            v, y, e;
106:   Mat            A;
107:   PetscReal      relerr, norm, rt1, rt2, cs1, anorm;
108:   PetscScalar    theta, rho, delta, sigma, alpha2, beta1, sn1;
109:   PetscTruth     breakdown,*select = PETSC_NULL,hasnorm;

112:   v = eps->V[0];
113:   y = eps->work[1];
114:   e = eps->work[0];

116:   /* prepare for selective orthogonalization of converged vectors */
117:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT && eps->nev>1) {
118:     STGetOperators(eps->OP,&A,PETSC_NULL);
119:     MatHasOperation(A,MATOP_NORM,&hasnorm);
120:     if (hasnorm) {
121:       MatNorm(A,NORM_INFINITY,&anorm);
122:       PetscMalloc(eps->nev*sizeof(PetscTruth),&select);
123:     }
124:   }

126:   EPSGetStartVector(eps,0,v,PETSC_NULL);
127:   STGetShift(eps->OP,&sigma);    /* original shift */
128:   rho = sigma;

130:   while (eps->reason == EPS_CONVERGED_ITERATING) {
131:     eps->its = eps->its + 1;

133:     /* y = OP v */
134:     STApply(eps->OP,v,y);

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

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

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

144:       /* compute relative error as ||y-theta v||_2/|theta| */
145:       VecCopy(y,e);
146:       VecAXPY(e,-theta,v);
147:       VecNorm(e,NORM_2,&norm);
148:       relerr = norm / PetscAbsScalar(theta);

150:     } else {  /* RQI */

152:       /* delta = ||y||_B */
153:       IPNorm(eps->ip,y,&norm);
154:       delta = norm;

156:       /* compute relative error */
157:       if (rho == 0.0) relerr = PETSC_MAX;
158:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

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

163:       /* compute new shift */
164:       if (relerr<eps->tol) {
165:         rho = sigma; /* if converged, restore original shift */
166:         STSetShift(eps->OP,rho);
167:       } else {
168:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v) */
169:         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
170: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
171:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
172: #else 
173:           /* beta1 is the norm of the residual associated to R(v) */
174:           VecAXPY(v,-theta/(delta*delta),y);
175:           VecScale(v,1.0/delta);
176:           IPNorm(eps->ip,v,&norm);
177:           beta1 = norm;
178: 
179:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
180:           STGetOperators(eps->OP,&A,PETSC_NULL);
181:           MatMult(A,v,e);
182:           VecDot(v,e,&alpha2);
183:           alpha2 = alpha2 / (beta1 * beta1);

185:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
186:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
187:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
188:           else rho = rt2;
189: #endif 
190:         }
191:         /* update operator according to new shift */
192:         PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
193:         STSetShift(eps->OP,rho);
194:         PetscPopErrorHandler();
195:         if (ierr) {
196:           eps->eigr[eps->nconv] = rho;
197:           relerr = PETSC_MACHINE_EPSILON;
198:           rho = sigma;
199:           STSetShift(eps->OP,rho);
200:         }
201:       }
202:     }

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

207:     /* purge previously converged eigenvectors */
208:     if (select) {
209:       for (i=0;i<eps->nconv;i++) {
210:         if(PetscAbsScalar(rho-eps->eigr[i])>eps->its*anorm/1000) select[i] = PETSC_TRUE;
211:         else select[i] = PETSC_FALSE;
212:       }
213:       IPOrthogonalize(eps->ip,eps->nds,eps->DS,eps->nconv,select,eps->V,y,PETSC_NULL,&norm,PETSC_NULL);
214:     } else {
215:       IPOrthogonalize(eps->ip,eps->nds,eps->DS,eps->nconv,PETSC_NULL,eps->V,y,PETSC_NULL,&norm,PETSC_NULL);
216:     }

218:     /* v = y/||y||_B */
219:     VecCopy(y,v);
220:     VecScale(v,1.0/norm);

222:     /* if relerr<tol, accept eigenpair */
223:     if (relerr<eps->tol) {
224:       eps->nconv = eps->nconv + 1;
225:       if (eps->nconv==eps->nev) eps->reason = EPS_CONVERGED_TOL;
226:       else {
227:         v = eps->V[eps->nconv];
228:         EPSGetStartVector(eps,eps->nconv,v,&breakdown);
229:         if (breakdown) {
230:           eps->reason = EPS_DIVERGED_BREAKDOWN;
231:           PetscInfo(eps,"Unable to generate more start vectors\n");
232:         }
233:       }
234:     }

236:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
237:   }

239:   PetscFree(select);

241:   return(0);
242: }

246: PetscErrorCode EPSSolve_TS_POWER(EPS eps)
247: {
249:   EPS_POWER      *power = (EPS_POWER *)eps->data;
250:   Vec            v, w, y, z, e;
251:   Mat            A;
252:   PetscReal      relerr, norm, rt1, rt2, cs1;
253:   PetscScalar    theta, alpha, beta, rho, delta, sigma, alpha2, beta1, sn1;

256:   v = eps->V[0];
257:   y = eps->work[1];
258:   e = eps->work[0];
259:   w = eps->W[0];
260:   z = eps->work[2];

262:   EPSGetStartVector(eps,0,v,PETSC_NULL);
263:   EPSGetStartVectorLeft(eps,0,w,PETSC_NULL);
264:   STGetShift(eps->OP,&sigma);    /* original shift */
265:   rho = sigma;

267:   while (eps->its<eps->max_it) {
268:     eps->its++;
269: 
270:     /* y = OP v, z = OP' w */
271:     STApply(eps->OP,v,y);
272:     STApplyTranspose(eps->OP,w,z);

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

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

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

282:       /* compute relative errors (right and left) */
283:       VecCopy(y,e);
284:       VecAXPY(e,-theta,v);
285:       VecNorm(e,NORM_2,&norm);
286:       relerr = norm / PetscAbsScalar(theta);
287:       eps->errest[eps->nconv] = relerr;
288:       VecCopy(z,e);
289:       VecAXPY(e,-theta,w);
290:       VecNorm(e,NORM_2,&norm);
291:       relerr = norm / PetscAbsScalar(theta);
292:       eps->errest_left[eps->nconv] = relerr;

294:     } else {  /* RQI */

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

301:       /* compute relative error */
302:       if (rho == 0.0) relerr = PETSC_MAX;
303:       else relerr = 1.0 / (PetscAbsScalar(delta*rho));
304:       eps->errest[eps->nconv] = relerr;
305:       eps->errest_left[eps->nconv] = relerr;

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

310:       /* compute new shift */
311:       if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
312:         rho = sigma; /* if converged, restore original shift */
313:         STSetShift(eps->OP,rho);
314:       } else {
315:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v,w) */
316:         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
317: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
318:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
319: #else 
320:           /* beta1 is the norm of the residual associated to R(v,w) */
321:           VecAXPY(v,-theta/(delta*delta),y);
322:           VecScale(v,1.0/delta);
323:           IPNorm(eps->ip,v,&norm);
324:           beta1 = norm;
325: 
326:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
327:           STGetOperators(eps->OP,&A,PETSC_NULL);
328:           MatMult(A,v,e);
329:           VecDot(v,e,&alpha2);
330:           alpha2 = alpha2 / (beta1 * beta1);

332:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
333:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
334:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
335:           else rho = rt2;
336: #endif 
337:         }
338:         /* update operator according to new shift */
339:         PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
340:         STSetShift(eps->OP,rho);
341:         PetscPopErrorHandler();
342:         if (ierr) {
343:           eps->eigr[eps->nconv] = rho;
344:           eps->errest[eps->nconv] = PETSC_MACHINE_EPSILON;
345:           eps->errest_left[eps->nconv] = PETSC_MACHINE_EPSILON;
346:           rho = sigma;
347:           STSetShift(eps->OP,rho);
348:         }
349:       }
350:     }

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

355:     /* purge previously converged eigenvectors */
356:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
357:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);

359:     /* normalize so that (y,z)_B=1  */
360:     VecCopy(y,v);
361:     VecCopy(z,w);
362:     IPInnerProduct(eps->ip,y,z,&alpha);
363:     if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
364:     delta = PetscSqrtScalar(PetscAbsScalar(alpha));
365:     beta = 1.0/PetscConj(alpha/delta);
366:     delta = 1.0/delta;
367:     VecScale(w,beta);
368:     VecScale(v,delta);

370:     /* if relerr<tol (both right and left), accept eigenpair */
371:     if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
372:       eps->nconv = eps->nconv + 1;
373:       if (eps->nconv==eps->nev) break;
374:       v = eps->V[eps->nconv];
375:       EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
376:       w = eps->W[eps->nconv];
377:       EPSGetStartVectorLeft(eps,eps->nconv,w,PETSC_NULL);
378:     }
379:   }

381:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
382:   else eps->reason = EPS_DIVERGED_ITS;

384:   return(0);
385: }

389: PetscErrorCode EPSBackTransform_POWER(EPS eps)
390: {
392:   EPS_POWER *power = (EPS_POWER *)eps->data;

395:   if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
396:     EPSBackTransform_Default(eps);
397:   }
398:   return(0);
399: }

403: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
404: {
406:   EPS_POWER      *power = (EPS_POWER *)eps->data;
407:   PetscTruth     flg;
408:   PetscInt       i;
409:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

412:   PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"POWER Options","EPS");
413:   PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],&i,&flg);
414:   if (flg ) power->shift_type = (EPSPowerShiftType)i;
415:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
416:     STSetType(eps->OP,STSINVERT);
417:   }
418:   PetscOptionsEnd();
419:   return(0);
420: }

425: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
426: {
427:   EPS_POWER *power = (EPS_POWER *)eps->data;

430:   switch (shift) {
431:     case EPS_POWER_SHIFT_CONSTANT:
432:     case EPS_POWER_SHIFT_RAYLEIGH:
433:     case EPS_POWER_SHIFT_WILKINSON:
434:       power->shift_type = shift;
435:       break;
436:     default:
437:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
438:   }
439:   return(0);
440: }

445: /*@
446:    EPSPowerSetShiftType - Sets the type of shifts used during the power
447:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
448:    (RQI) method.

450:    Collective on EPS

452:    Input Parameters:
453: +  eps - the eigenproblem solver context
454: -  shift - the type of shift

456:    Options Database Key:
457: .  -eps_power_shift_type - Sets the shift type (either 'constant' or 
458:                            'rayleigh' or 'wilkinson')

460:    Notes:
461:    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
462:    is the simple power method (or inverse iteration if a shift-and-invert
463:    transformation is being used).

465:    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
466:    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
467:    a cubic converging method as RQI. See the users manual for details.
468:    
469:    Level: advanced

471: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
472: @*/
473: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
474: {
475:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);

479:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
480:   if (f) {
481:     (*f)(eps,shift);
482:   }
483:   return(0);
484: }

489: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
490: {
491:   EPS_POWER  *power = (EPS_POWER *)eps->data;
493:   *shift = power->shift_type;
494:   return(0);
495: }

500: /*@C
501:    EPSPowerGetShiftType - Gets the type of shifts used during the power
502:    iteration. 

504:    Collective on EPS

506:    Input Parameter:
507: .  eps - the eigenproblem solver context

509:    Input Parameter:
510: .  shift - the type of shift

512:    Level: advanced

514: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
515: @*/
516: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
517: {
518:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);

522:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
523:   if (f) {
524:     (*f)(eps,shift);
525:   }
526:   return(0);
527: }

531: PetscErrorCode EPSDestroy_POWER(EPS eps)
532: {

537:   EPSDestroy_Default(eps);
538:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","",PETSC_NULL);
539:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","",PETSC_NULL);
540:   return(0);
541: }

545: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
546: {
548:   EPS_POWER      *power = (EPS_POWER *)eps->data;
549:   PetscTruth     isascii;
550:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

553:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
554:   if (!isascii) {
555:     SETERRQ1(1,"Viewer type %s not supported for EPS_POWER",((PetscObject)viewer)->type_name);
556:   }
557:   PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
558:   return(0);
559: }

564: PetscErrorCode EPSCreate_POWER(EPS eps)
565: {
567:   EPS_POWER      *power;

570:   PetscNew(EPS_POWER,&power);
571:   PetscLogObjectMemory(eps,sizeof(EPS_POWER));
572:   eps->data                      = (void *) power;
573:   eps->ops->setup                = EPSSetUp_POWER;
574:   eps->ops->setfromoptions       = EPSSetFromOptions_POWER;
575:   eps->ops->destroy              = EPSDestroy_POWER;
576:   eps->ops->view                 = EPSView_POWER;
577:   eps->ops->backtransform        = EPSBackTransform_POWER;
578:   eps->ops->computevectors       = EPSComputeVectors_Default;
579:   power->shift_type              = EPS_POWER_SHIFT_CONSTANT;
580:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
581:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
582:   return(0);
583: }