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

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

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

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

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

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

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

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

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

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

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

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

148:     } else {  /* RQI */

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

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

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

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

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

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

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

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

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

241: PetscErrorCode EPSSolve_TS_Power(EPS eps)
242: {
244:   EPS_POWER      *power = (EPS_POWER *)eps->data;
245:   Vec            v,w,y,z,e;
246:   Mat            A;
247:   PetscReal      relerr,norm,rt1,rt2,cs1;
248:   PetscScalar    theta,alpha,beta,rho,delta,sigma,alpha2,beta1,sn1;

251:   v = eps->V[0];
252:   y = eps->work[1];
253:   e = eps->work[0];
254:   w = eps->W[0];
255:   z = eps->work[2];

257:   EPSGetStartVector(eps,0,v,PETSC_NULL);
258:   EPSGetStartVectorLeft(eps,0,w,PETSC_NULL);
259:   STGetShift(eps->OP,&sigma);    /* original shift */
260:   rho = sigma;

262:   while (eps->its<eps->max_it) {
263:     eps->its++;
264: 
265:     /* y = OP v, z = OP' w */
266:     STApply(eps->OP,v,y);
267:     STApplyTranspose(eps->OP,w,z);

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

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

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

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

289:     } else {  /* RQI */

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

296:       /* compute relative error */
297:       if (rho == 0.0) relerr = PETSC_MAX_REAL;
298:       else relerr = 1.0 / (PetscAbsScalar(delta*rho));
299:       eps->errest[eps->nconv] = relerr;
300:       eps->errest_left[eps->nconv] = relerr;

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

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

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

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

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

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

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

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

388:   if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
389:     EPSBackTransform_Default(eps);
390:   }
391:   return(0);
392: }

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

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

414: EXTERN_C_BEGIN
417: PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
418: {
419:   EPS_POWER *power = (EPS_POWER *)eps->data;

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

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

442:    Logically Collective on EPS

444:    Input Parameters:
445: +  eps - the eigenproblem solver context
446: -  shift - the type of shift

448:    Options Database Key:
449: .  -eps_power_shift_type - Sets the shift type (either 'constant' or 
450:                            'rayleigh' or 'wilkinson')

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

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

463: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
464: @*/
465: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
466: {

472:   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
473:   return(0);
474: }

476: EXTERN_C_BEGIN
479: PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
480: {
481:   EPS_POWER  *power = (EPS_POWER *)eps->data;

484:   *shift = power->shift_type;
485:   return(0);
486: }
487: EXTERN_C_END

491: /*@C
492:    EPSPowerGetShiftType - Gets the type of shifts used during the power
493:    iteration. 

495:    Not Collective

497:    Input Parameter:
498: .  eps - the eigenproblem solver context

500:    Input Parameter:
501: .  shift - the type of shift

503:    Level: advanced

505: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
506: @*/
507: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
508: {

514:   PetscTryMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
515:   return(0);
516: }

520: PetscErrorCode EPSDestroy_Power(EPS eps)
521: {

525:   PetscFree(eps->data);
526:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","",PETSC_NULL);
527:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","",PETSC_NULL);
528:   return(0);
529: }

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

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

548: EXTERN_C_BEGIN
551: PetscErrorCode EPSCreate_Power(EPS eps)
552: {

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