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:  #include src/eps/epsimpl.h
 25:  #include slepcblaslapack.h

 27: typedef struct {
 28:   EPSPowerShiftType shift_type;
 29: } EPS_POWER;

 33: PetscErrorCode EPSSetUp_POWER(EPS eps)
 34: {
 36:   EPS_POWER      *power = (EPS_POWER *)eps->data;
 37:   PetscInt       N;
 38:   PetscTruth     flg;
 39:   STMatMode      mode;

 42:   VecGetSize(eps->vec_initial,&N);
 43:   if (eps->ncv) {
 44:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 45:   }
 46:   else eps->ncv = eps->nev;
 47:   if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
 48:   if (eps->which!=EPS_LARGEST_MAGNITUDE)
 49:     SETERRQ(1,"Wrong value of eps->which");
 50:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
 51:     PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
 52:     if (!flg)
 53:       SETERRQ(PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert ST");
 54:     STGetMatMode(eps->OP,&mode);
 55:     if (mode == STMATMODE_INPLACE)
 56:       SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 57:   }
 58:   EPSAllocateSolution(eps);
 59:   EPSDefaultGetWork(eps,1);
 60:   return(0);
 61: }

 65: PetscErrorCode EPSSolve_POWER(EPS eps)
 66: {
 68:   EPS_POWER      *power = (EPS_POWER *)eps->data;
 69:   int            i, nsv;
 70:   Vec            v, y, e, *SV;
 71:   Mat            A;
 72:   PetscReal      relerr, norm, rt1, rt2, cs1, anorm;
 73:   PetscScalar    theta, rho, delta, sigma, alpha2, beta1, sn1;
 74:   PetscTruth     breakdown;

 77:   v = eps->V[0];
 78:   y = eps->AV[0];
 79:   e = eps->work[0];

 81:   /* prepare for selective orthogonalization of converged vectors */
 82:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
 83:     PetscMalloc(eps->nev*sizeof(Vec),&SV);
 84:     for (i=0;i<eps->nds;i++) SV[i]=eps->DS[i];
 85:     if (eps->nev>1) {
 86:       STGetOperators(eps->OP,&A,PETSC_NULL);
 87:       MatNorm(A,NORM_INFINITY,&anorm);
 88:     }
 89:   }

 91:   EPSGetStartVector(eps,0,v,PETSC_NULL);
 92:   STGetShift(eps->OP,&sigma);    /* original shift */
 93:   rho = sigma;

 95:   while (eps->reason == EPS_CONVERGED_ITERATING) {
 96:     eps->its = eps->its + 1;

 98:     /* y = OP v */
 99:     STApply(eps->OP,v,y);

101:     /* theta = (v,y)_B */
102:     STInnerProduct(eps->OP,v,y,&theta);

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

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

109:       /* compute relative error as ||y-theta v||_2/|theta| */
110:       VecCopy(y,e);
111:       VecAXPY(e,-theta,v);
112:       VecNorm(e,NORM_2,&norm);
113:       relerr = norm / PetscAbsScalar(theta);

115:     } else {  /* RQI */

117:       /* delta = ||y||_B */
118:       STNorm(eps->OP,y,&norm);
119:       delta = norm;

121:       /* compute relative error */
122:       if (rho == 0.0) relerr = PETSC_MAX;
123:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

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

128:       /* compute new shift */
129:       if (relerr<eps->tol) {
130:         rho = sigma; /* if converged, restore original shift */
131:         STSetShift(eps->OP,rho);
132:       } else {
133:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v) */
134:         if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
135: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
136:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
137: #else 
138:           /* beta1 is the norm of the residual associated to R(v) */
139:           VecAXPY(v,-theta/(delta*delta),y);
140:           VecScale(v,1.0/delta);
141:           STNorm(eps->OP,v,&norm);
142:           beta1 = norm;
143: 
144:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
145:           STGetOperators(eps->OP,&A,PETSC_NULL);
146:           MatMult(A,v,e);
147:           VecDot(v,e,&alpha2);
148:           alpha2 = alpha2 / (beta1 * beta1);

150:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
151:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
152:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
153:           else rho = rt2;
154: #endif 
155:         }
156:         /* update operator according to new shift */
157:         PetscPushErrorHandler(SlepcQuietErrorHandler,PETSC_NULL);
158:         STSetShift(eps->OP,rho);
159:         PetscPopErrorHandler();
160:         if (ierr) {
161:           eps->eigr[eps->nconv] = rho;
162:           relerr = PETSC_MACHINE_EPSILON;
163:           rho = sigma;
164:           STSetShift(eps->OP,rho);
165:         }
166:       }
167:     }

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

172:     /* purge previously converged eigenvectors */
173:     if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
174:       nsv = eps->nds;
175:       for (i=0;i<eps->nconv;i++) {
176:         if(PetscAbsScalar(rho-eps->eigr[i])>eps->its*anorm/1000) SV[nsv++]=eps->V[i];
177:       }
178:       EPSOrthogonalize(eps,nsv,PETSC_NULL,SV,y,PETSC_NULL,&norm,PETSC_NULL);
179:     } else {
180:       EPSOrthogonalize(eps,eps->nds+eps->nconv,PETSC_NULL,eps->DSV,y,PETSC_NULL,&norm,PETSC_NULL);
181:     }

183:     /* v = y/||y||_B */
184:     VecCopy(y,v);
185:     VecScale(v,1.0/norm);

187:     /* if relerr<tol, accept eigenpair */
188:     if (relerr<eps->tol) {
189:       eps->nconv = eps->nconv + 1;
190:       if (eps->nconv==eps->nev) eps->reason = EPS_CONVERGED_TOL;
191:       else {
192:         v = eps->V[eps->nconv];
193:         EPSGetStartVector(eps,eps->nconv,v,&breakdown);
194:         if (breakdown) {
195:           eps->reason = EPS_DIVERGED_BREAKDOWN;
196:           PetscInfo(eps,"Unable to generate more start vectors\n");
197:         }
198:       }
199:     }

201:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
202:   }

204:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
205:     PetscFree(SV);
206:   }

208:   return(0);
209: }

213: PetscErrorCode EPSSolve_TS_POWER(EPS eps)
214: {
216:   EPS_POWER      *power = (EPS_POWER *)eps->data;
217:   Vec            v, w, y, z, e;
218:   Mat            A;
219:   PetscReal      relerr, norm, rt1, rt2, cs1;
220:   PetscScalar    theta, alpha, beta, rho, delta, sigma, alpha2, beta1, sn1;

223:   v = eps->V[0];
224:   y = eps->AV[0];
225:   e = eps->work[0];
226:   w = eps->W[0];
227:   z = eps->AW[0];

229:   EPSGetStartVector(eps,0,v,PETSC_NULL);
230:   EPSGetLeftStartVector(eps,0,w);
231:   STGetShift(eps->OP,&sigma);    /* original shift */
232:   rho = sigma;

234:   while (eps->its<eps->max_it) {
235:     eps->its++;
236: 
237:     /* y = OP v, z = OP' w */
238:     STApply(eps->OP,v,y);
239:     STApplyTranspose(eps->OP,w,z);

241:     /* theta = (v,z)_B */
242:     STInnerProduct(eps->OP,v,z,&theta);

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

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

249:       /* compute relative errors (right and left) */
250:       VecCopy(y,e);
251:       VecAXPY(e,-theta,v);
252:       VecNorm(e,NORM_2,&norm);
253:       relerr = norm / PetscAbsScalar(theta);
254:       eps->errest[eps->nconv] = relerr;
255:       VecCopy(z,e);
256:       VecAXPY(e,-theta,w);
257:       VecNorm(e,NORM_2,&norm);
258:       relerr = norm / PetscAbsScalar(theta);
259:       eps->errest_left[eps->nconv] = relerr;

261:     } else {  /* RQI */

263:       /* delta = sqrt(y,z)_B */
264:       STInnerProduct(eps->OP,y,z,&alpha);
265:       if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
266:       delta = PetscSqrtScalar(alpha);

268:       /* compute relative error */
269:       if (rho == 0.0) relerr = PETSC_MAX;
270:       else relerr = 1.0 / (PetscAbsScalar(delta*rho));
271:       eps->errest[eps->nconv] = relerr;
272:       eps->errest_left[eps->nconv] = relerr;

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

277:       /* compute new shift */
278:       if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
279:         rho = sigma; /* if converged, restore original shift */
280:         STSetShift(eps->OP,rho);
281:       } else {
282:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v,w) */
283:         if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
284: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
285:           SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
286: #else 
287:           /* beta1 is the norm of the residual associated to R(v,w) */
288:           VecAXPY(v,-theta/(delta*delta),y);
289:           VecScale(v,1.0/delta);
290:           STNorm(eps->OP,v,&norm);
291:           beta1 = norm;
292: 
293:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
294:           STGetOperators(eps->OP,&A,PETSC_NULL);
295:           MatMult(A,v,e);
296:           VecDot(v,e,&alpha2);
297:           alpha2 = alpha2 / (beta1 * beta1);

299:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
300:           LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
301:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
302:           else rho = rt2;
303: #endif 
304:         }
305:         /* update operator according to new shift */
306:         PetscPushErrorHandler(SlepcQuietErrorHandler,PETSC_NULL);
307:         STSetShift(eps->OP,rho);
308:         PetscPopErrorHandler();
309:         if (ierr) {
310:           eps->eigr[eps->nconv] = rho;
311:           eps->errest[eps->nconv] = PETSC_MACHINE_EPSILON;
312:           eps->errest_left[eps->nconv] = PETSC_MACHINE_EPSILON;
313:           rho = sigma;
314:           STSetShift(eps->OP,rho);
315:         }
316:       }
317:     }

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

322:     /* purge previously converged eigenvectors */
323:     EPSBiOrthogonalize(eps,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
324:     EPSBiOrthogonalize(eps,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);

326:     /* normalize so that (y,z)_B=1  */
327:     VecCopy(y,v);
328:     VecCopy(z,w);
329:     STInnerProduct(eps->OP,y,z,&alpha);
330:     if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
331:     delta = PetscSqrtScalar(PetscAbsScalar(alpha));
332:     beta = 1.0/PetscConj(alpha/delta);
333:     delta = 1.0/delta;
334:     VecScale(w,beta);
335:     VecScale(v,delta);

337:     /* if relerr<tol (both right and left), accept eigenpair */
338:     if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
339:       eps->nconv = eps->nconv + 1;
340:       if (eps->nconv==eps->nev) break;
341:       v = eps->V[eps->nconv];
342:       EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
343:       w = eps->W[eps->nconv];
344:       EPSGetLeftStartVector(eps,eps->nconv,w);
345:     }
346:   }

348:   if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
349:   else eps->reason = EPS_DIVERGED_ITS;

351:   return(0);
352: }

356: PetscErrorCode EPSBackTransform_POWER(EPS eps)
357: {
359:   EPS_POWER *power = (EPS_POWER *)eps->data;

362:   if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) {
363:     EPSBackTransform_Default(eps);
364:   }
365:   return(0);
366: }

370: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
371: {
373:   EPS_POWER      *power = (EPS_POWER *)eps->data;
374:   PetscTruth     flg;
375:   PetscInt       i;
376:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

379:   PetscOptionsHead("POWER options");
380:   PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],&i,&flg);
381:   if (flg ) power->shift_type = (EPSPowerShiftType)i;
382:   if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
383:     STSetType(eps->OP,STSINV);
384:   }
385:   PetscOptionsTail();
386:   return(0);
387: }

392: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
393: {
394:   EPS_POWER *power = (EPS_POWER *)eps->data;

397:   switch (shift) {
398:     case EPSPOWER_SHIFT_CONSTANT:
399:     case EPSPOWER_SHIFT_RAYLEIGH:
400:     case EPSPOWER_SHIFT_WILKINSON:
401:       power->shift_type = shift;
402:       break;
403:     default:
404:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
405:   }
406:   return(0);
407: }

412: /*@
413:    EPSPowerSetShiftType - Sets the type of shifts used during the power
414:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
415:    (RQI) method.

417:    Collective on EPS

419:    Input Parameters:
420: +  eps - the eigenproblem solver context
421: -  shift - the type of shift

423:    Options Database Key:
424: .  -eps_power_shift_type - Sets the shift type (either 'constant' or 
425:                            'rayleigh' or 'wilkinson')

427:    Notes:
428:    By default, shifts are constant (EPSPOWER_SHIFT_CONSTANT) and the iteration
429:    is the simple power method (or inverse iteration if a shift-and-invert
430:    transformation is being used).

432:    A variable shift can be specified (EPSPOWER_SHIFT_RAYLEIGH or
433:    EPSPOWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
434:    a cubic converging method as RQI. See the users manual for details.
435:    
436:    Level: advanced

438: .seealso: EPSGetShiftType(), STSetShift()
439: @*/
440: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
441: {
442:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);

446:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
447:   if (f) {
448:     (*f)(eps,shift);
449:   }
450:   return(0);
451: }

456: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
457: {
458:   EPS_POWER  *power = (EPS_POWER *)eps->data;
460:   *shift = power->shift_type;
461:   return(0);
462: }

467: /*@C
468:    EPSPowerGetShiftType - Gets the type of shifts used during the power
469:    iteration. 

471:    Collective on EPS

473:    Input Parameter:
474: .  eps - the eigenproblem solver context

476:    Input Parameter:
477: .  shift - the type of shift

479:    Level: advanced

481: .seealso: EPSSetShiftType()
482: @*/
483: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
484: {
485:   PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);

489:   PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
490:   if (f) {
491:     (*f)(eps,shift);
492:   }
493:   return(0);
494: }

498: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
499: {
501:   EPS_POWER      *power = (EPS_POWER *)eps->data;
502:   PetscTruth     isascii;
503:   const char     *shift_list[3] = { "constant", "rayleigh", "wilkinson" };

506:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
507:   if (!isascii) {
508:     SETERRQ1(1,"Viewer type %s not supported for EPSPOWER",((PetscObject)viewer)->type_name);
509:   }
510:   PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
511:   return(0);
512: }

517: PetscErrorCode EPSCreate_POWER(EPS eps)
518: {
520:   EPS_POWER      *power;

523:   PetscNew(EPS_POWER,&power);
524:   PetscLogObjectMemory(eps,sizeof(EPS_POWER));
525:   eps->data                      = (void *) power;
526:   eps->ops->solve                = EPSSolve_POWER;
527:   eps->ops->solvets              = EPSSolve_TS_POWER;
528:   eps->ops->setup                = EPSSetUp_POWER;
529:   eps->ops->setfromoptions       = EPSSetFromOptions_POWER;
530:   eps->ops->destroy              = EPSDestroy_Default;
531:   eps->ops->view                 = EPSView_POWER;
532:   eps->ops->backtransform        = EPSBackTransform_POWER;
533:   eps->ops->computevectors       = EPSComputeVectors_Default;
534:   power->shift_type              = EPSPOWER_SHIFT_CONSTANT;
535:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
536:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
537:   return(0);
538: }