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-2012, Universitat 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 <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 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:   PetscBool      flg;
 60:   STMatMode      mode;

 63:   if (eps->ncv) {
 64:     if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,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) { EPSDefaultSetWhich(eps); }
 70:   if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
 71:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
 72:     PetscObjectTypeCompareAny((PetscObject)eps->OP,&flg,STSINVERT,STCAYLEY,"");
 73:     if (!flg) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
 74:     STGetMatMode(eps->OP,&mode);
 75:     if (mode == ST_MATMODE_INPLACE) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 76:   }
 77:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 78:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Balancing not supported in this solver");
 79:   if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 80:   EPSAllocateSolution(eps);
 81:   if (eps->leftvecs) {
 82:     EPSDefaultGetWork(eps,3);
 83:   } else {
 84:     EPSDefaultGetWork(eps,2);
 85:   }

 87:   /* dispatch solve method */
 88:   if (eps->leftvecs) eps->ops->solve = EPSSolve_TS_Power;
 89:   else eps->ops->solve = EPSSolve_Power;
 90:   return(0);
 91: }

 95: PetscErrorCode EPSSolve_Power(EPS eps)
 96: {
 97: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
 99:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
100: #else 
102:   EPS_POWER      *power = (EPS_POWER *)eps->data;
103:   PetscInt       i;
104:   Vec            v,y,e;
105:   Mat            A;
106:   PetscReal      relerr,norm,rt1,rt2,cs1,anorm;
107:   PetscScalar    theta,rho,delta,sigma,alpha2,beta1,sn1;
108:   PetscBool      breakdown,*select = PETSC_NULL,hasnorm;

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

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

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

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

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

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

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

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

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

149:     } else {  /* RQI */

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

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

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

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

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

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

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

215:     /* v = y/||y||_B */
216:     VecCopy(y,v);
217:     VecScale(v,1.0/norm);

219:     /* if relerr<tol, accept eigenpair */
220:     if (relerr<eps->tol) {
221:       eps->nconv = eps->nconv + 1;
222:       if (eps->nconv==eps->nev) eps->reason = EPS_CONVERGED_TOL;
223:       else {
224:         v = eps->V[eps->nconv];
225:         EPSGetStartVector(eps,eps->nconv,v,&breakdown);
226:         if (breakdown) {
227:           eps->reason = EPS_DIVERGED_BREAKDOWN;
228:           PetscInfo(eps,"Unable to generate more start vectors\n");
229:         }
230:       }
231:     }
232:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
233:   }
234:   PetscFree(select);
235:   return(0);
236: #endif 
237: }

241: PetscErrorCode EPSSolve_TS_Power(EPS eps)
242: {
243: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
245:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
246: #else 
248:   EPS_POWER      *power = (EPS_POWER *)eps->data;
249:   Vec            v,w,y,z,e;
250:   Mat            A;
251:   PetscReal      relerr,norm,rt1,rt2,cs1;
252:   PetscScalar    theta,alpha,beta,rho,delta,sigma,alpha2,beta1,sn1;

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

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

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

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

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

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

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

293:     } else {  /* RQI */

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

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

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

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

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

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

352:     /* purge previously converged eigenvectors */
353:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
354:     IPBiOrthogonalize(eps->ip,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);

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

367:     /* if relerr<tol (both right and left), accept eigenpair */
368:     if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
369:       eps->nconv = eps->nconv + 1;
370:       if (eps->nconv==eps->nev) break;
371:       v = eps->V[eps->nconv];
372:       EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
373:       w = eps->W[eps->nconv];
374:       EPSGetStartVectorLeft(eps,eps->nconv,w,PETSC_NULL);
375:     }
376:   }
377:   if (eps->nconv == eps->nev) eps->reason = EPS_CONVERGED_TOL;
378:   else eps->reason = EPS_DIVERGED_ITS;
379:   return(0);
380: #endif 
381: }

385: PetscErrorCode EPSBackTransform_Power(EPS eps)
386: {
388:   EPS_POWER      *power = (EPS_POWER *)eps->data;

391:   if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
392:     EPSBackTransform_Default(eps);
393:   }
394:   return(0);
395: }

399: PetscErrorCode EPSSetFromOptions_Power(EPS eps)
400: {
401:   PetscErrorCode    ierr;
402:   EPS_POWER         *power = (EPS_POWER *)eps->data;
403:   PetscBool         flg;
404:   EPSPowerShiftType shift;

407:   PetscOptionsHead("EPS Power Options");
408:   PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
409:   if (flg) { EPSPowerSetShiftType(eps,shift); }
410:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
411:     STSetType(eps->OP,STSINVERT);
412:   }
413:   PetscOptionsTail();
414:   return(0);
415: }

417: EXTERN_C_BEGIN
420: PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
421: {
422:   EPS_POWER *power = (EPS_POWER *)eps->data;

425:   switch (shift) {
426:     case EPS_POWER_SHIFT_CONSTANT:
427:     case EPS_POWER_SHIFT_RAYLEIGH:
428:     case EPS_POWER_SHIFT_WILKINSON:
429:       power->shift_type = shift;
430:       break;
431:     default:
432:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
433:   }
434:   return(0);
435: }
436: EXTERN_C_END

440: /*@
441:    EPSPowerSetShiftType - Sets the type of shifts used during the power
442:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
443:    (RQI) method.

445:    Logically Collective on EPS

447:    Input Parameters:
448: +  eps - the eigenproblem solver context
449: -  shift - the type of shift

451:    Options Database Key:
452: .  -eps_power_shift_type - Sets the shift type (either 'constant' or 
453:                            'rayleigh' or 'wilkinson')

455:    Notes:
456:    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
457:    is the simple power method (or inverse iteration if a shift-and-invert
458:    transformation is being used).

460:    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
461:    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
462:    a cubic converging method as RQI. See the users manual for details.
463:    
464:    Level: advanced

466: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
467: @*/
468: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
469: {

475:   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
476:   return(0);
477: }

479: EXTERN_C_BEGIN
482: PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
483: {
484:   EPS_POWER  *power = (EPS_POWER *)eps->data;

487:   *shift = power->shift_type;
488:   return(0);
489: }
490: EXTERN_C_END

494: /*@C
495:    EPSPowerGetShiftType - Gets the type of shifts used during the power
496:    iteration. 

498:    Not Collective

500:    Input Parameter:
501: .  eps - the eigenproblem solver context

503:    Input Parameter:
504: .  shift - the type of shift

506:    Level: advanced

508: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
509: @*/
510: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
511: {

517:   PetscTryMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
518:   return(0);
519: }

523: PetscErrorCode EPSDestroy_Power(EPS eps)
524: {

528:   PetscFree(eps->data);
529:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","",PETSC_NULL);
530:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","",PETSC_NULL);
531:   return(0);
532: }

536: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
537: {
539:   EPS_POWER      *power = (EPS_POWER *)eps->data;
540:   PetscBool      isascii;

543:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
544:   if (!isascii) SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPS Power",((PetscObject)viewer)->type_name);
545:   PetscViewerASCIIPrintf(viewer,"  Power: %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
546:   return(0);
547: }

549: EXTERN_C_BEGIN
552: PetscErrorCode EPSCreate_Power(EPS eps)
553: {

557:   PetscNewLog(eps,EPS_POWER,&eps->data);
558:   eps->ops->setup                = EPSSetUp_Power;
559:   eps->ops->setfromoptions       = EPSSetFromOptions_Power;
560:   eps->ops->destroy              = EPSDestroy_Power;
561:   eps->ops->reset                = EPSReset_Default;
562:   eps->ops->view                 = EPSView_Power;
563:   eps->ops->backtransform        = EPSBackTransform_Power;
564:   eps->ops->computevectors       = EPSComputeVectors_Default;
565:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_Power",EPSPowerSetShiftType_Power);
566:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_Power",EPSPowerGetShiftType_Power);
567:   return(0);
568: }
569: EXTERN_C_END