Actual source code: power.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  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
 19:            STR-2, available at http://slepc.upv.es.

 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 23:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 25:    This file is part of SLEPc.

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

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

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

 41: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 42: #include <slepcblaslapack.h>

 44: typedef struct {
 45:   EPSPowerShiftType shift_type;
 46: } EPS_POWER;

 50: PetscErrorCode EPSSetUp_Power(EPS eps)
 51: {
 53:   EPS_POWER      *power = (EPS_POWER*)eps->data;
 54:   PetscBool      flg,istrivial;
 55:   STMatMode      mode;

 58:   if (eps->ncv) {
 59:     if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
 60:   } else eps->ncv = eps->nev;
 61:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 62:   if (!eps->max_it) eps->max_it = PetscMax(2000,100*eps->n);
 63:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 64:   if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 65:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
 66:     PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
 67:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
 68:     STGetMatMode(eps->st,&mode);
 69:     if (mode == ST_MATMODE_INPLACE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 70:   }
 71:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 72:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
 73:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 74:   RGIsTrivial(eps->rg,&istrivial);
 75:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
 76:   EPSAllocateSolution(eps,0);
 77:   EPS_SetInnerProduct(eps);
 78:   EPSSetWorkVecs(eps,2);
 79:   return(0);
 80: }

 84: PetscErrorCode EPSSolve_Power(EPS eps)
 85: {
 86: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
 88:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
 89: #else
 91:   EPS_POWER      *power = (EPS_POWER*)eps->data;
 92:   PetscInt       k;
 93:   Vec            v,y,e;
 94:   Mat            A;
 95:   PetscReal      relerr,norm,rt1,rt2,cs1;
 96:   PetscScalar    theta,rho,delta,sigma,alpha2,beta1,sn1;
 97:   PetscBool      breakdown;

100:   y = eps->work[1];
101:   e = eps->work[0];

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

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

111:     /* y = OP v */
112:     BVGetColumn(eps->V,k,&v);
113:     STApply(eps->st,v,y);
114:     BVRestoreColumn(eps->V,k,&v);

116:     /* theta = (v,y)_B */
117:     BVSetActiveColumns(eps->V,k,k+1);
118:     BVDotVec(eps->V,y,&theta);

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

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

125:       /* compute relative error as ||y-theta v||_2/|theta| */
126:       VecCopy(y,e);
127:       BVGetColumn(eps->V,k,&v);
128:       VecAXPY(e,-theta,v);
129:       BVRestoreColumn(eps->V,k,&v);
130:       VecNorm(e,NORM_2,&norm);
131:       relerr = norm / PetscAbsScalar(theta);

133:     } else {  /* RQI */

135:       /* delta = ||y||_B */
136:       BVNormVec(eps->V,y,NORM_2,&norm);
137:       delta = norm;

139:       /* compute relative error */
140:       if (rho == 0.0) relerr = PETSC_MAX_REAL;
141:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

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

146:       /* compute new shift */
147:       if (relerr<eps->tol) {
148:         rho = sigma; /* if converged, restore original shift */
149:         STSetShift(eps->st,rho);
150:       } else {
151:         rho = rho + theta/(delta*delta);  /* Rayleigh quotient R(v) */
152:         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
153:           /* beta1 is the norm of the residual associated to R(v) */
154:           BVGetColumn(eps->V,k,&v);
155:           VecAXPY(v,-theta/(delta*delta),y);
156:           BVRestoreColumn(eps->V,k,&v);
157:           BVScaleColumn(eps->V,k,1.0/delta);
158:           BVNormColumn(eps->V,k,NORM_2,&norm);
159:           beta1 = norm;

161:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
162:           STGetOperators(eps->st,0,&A);
163:           BVGetColumn(eps->V,k,&v);
164:           MatMult(A,v,e);
165:           VecDot(v,e,&alpha2);
166:           BVRestoreColumn(eps->V,k,&v);
167:           alpha2 = alpha2 / (beta1 * beta1);

169:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
170:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
171:           PetscStackCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
172:           PetscFPTrapPop();
173:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
174:           else rho = rt2;
175:         }
176:         /* update operator according to new shift */
177:         PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
178:         STSetShift(eps->st,rho);
179:         PetscPopErrorHandler();
180:         if (ierr) {
181:           eps->eigr[eps->nconv] = rho;
182:           relerr = PETSC_MACHINE_EPSILON;
183:           rho = sigma;
184:           STSetShift(eps->st,rho);
185:         }
186:       }
187:     }
188:     eps->errest[eps->nconv] = relerr;

190:     /* purge previously converged eigenvectors */
191:     BVInsertVec(eps->V,k,y);
192:     BVOrthogonalizeColumn(eps->V,k,NULL,&norm,NULL);
193:     BVScaleColumn(eps->V,k,1.0/norm);

195:     /* if relerr<tol, accept eigenpair */
196:     if (relerr<eps->tol) {
197:       eps->nconv = eps->nconv + 1;
198:       if (eps->nconv<eps->nev) {
199:         EPSGetStartVector(eps,eps->nconv,&breakdown);
200:         if (breakdown) {
201:           eps->reason = EPS_DIVERGED_BREAKDOWN;
202:           PetscInfo(eps,"Unable to generate more start vectors\n");
203:           break;
204:         }
205:       }
206:     }
207:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
208:     (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
209:   }
210:   return(0);
211: #endif
212: }

216: PetscErrorCode EPSBackTransform_Power(EPS eps)
217: {
219:   EPS_POWER      *power = (EPS_POWER*)eps->data;

222:   if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
223:     EPSBackTransform_Default(eps);
224:   }
225:   return(0);
226: }

230: PetscErrorCode EPSSetFromOptions_Power(PetscOptionItems *PetscOptionsObject,EPS eps)
231: {
232:   PetscErrorCode    ierr;
233:   EPS_POWER         *power = (EPS_POWER*)eps->data;
234:   PetscBool         flg;
235:   EPSPowerShiftType shift;

238:   PetscOptionsHead(PetscOptionsObject,"EPS Power Options");
239:   PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
240:   if (flg) {
241:     EPSPowerSetShiftType(eps,shift);
242:   }
243:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
244:     STSetType(eps->st,STSINVERT);
245:   }
246:   PetscOptionsTail();
247:   return(0);
248: }

252: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
253: {
254:   EPS_POWER *power = (EPS_POWER*)eps->data;

257:   switch (shift) {
258:     case EPS_POWER_SHIFT_CONSTANT:
259:     case EPS_POWER_SHIFT_RAYLEIGH:
260:     case EPS_POWER_SHIFT_WILKINSON:
261:       power->shift_type = shift;
262:       break;
263:     default:
264:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
265:   }
266:   return(0);
267: }

271: /*@
272:    EPSPowerSetShiftType - Sets the type of shifts used during the power
273:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
274:    (RQI) method.

276:    Logically Collective on EPS

278:    Input Parameters:
279: +  eps - the eigenproblem solver context
280: -  shift - the type of shift

282:    Options Database Key:
283: .  -eps_power_shift_type - Sets the shift type (either 'constant' or
284:                            'rayleigh' or 'wilkinson')

286:    Notes:
287:    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
288:    is the simple power method (or inverse iteration if a shift-and-invert
289:    transformation is being used).

291:    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
292:    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
293:    a cubic converging method as RQI. See the users manual for details.

295:    Level: advanced

297: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
298: @*/
299: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
300: {

306:   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
307:   return(0);
308: }

312: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
313: {
314:   EPS_POWER  *power = (EPS_POWER*)eps->data;

317:   *shift = power->shift_type;
318:   return(0);
319: }

323: /*@
324:    EPSPowerGetShiftType - Gets the type of shifts used during the power
325:    iteration.

327:    Not Collective

329:    Input Parameter:
330: .  eps - the eigenproblem solver context

332:    Input Parameter:
333: .  shift - the type of shift

335:    Level: advanced

337: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
338: @*/
339: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
340: {

346:   PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
347:   return(0);
348: }

352: PetscErrorCode EPSDestroy_Power(EPS eps)
353: {

357:   PetscFree(eps->data);
358:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
359:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
360:   return(0);
361: }

365: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
366: {
368:   EPS_POWER      *power = (EPS_POWER*)eps->data;
369:   PetscBool      isascii;

372:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
373:   if (isascii) {
374:     PetscViewerASCIIPrintf(viewer,"  Power: %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
375:   }
376:   return(0);
377: }

381: PETSC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
382: {
383:   EPS_POWER      *ctx;

387:   PetscNewLog(eps,&ctx);
388:   eps->data = (void*)ctx;

390:   eps->ops->setup                = EPSSetUp_Power;
391:   eps->ops->solve                = EPSSolve_Power;
392:   eps->ops->setfromoptions       = EPSSetFromOptions_Power;
393:   eps->ops->destroy              = EPSDestroy_Power;
394:   eps->ops->view                 = EPSView_Power;
395:   eps->ops->backtransform        = EPSBackTransform_Power;
396:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
397:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
398:   return(0);
399: }