Actual source code: power.c

slepc-3.12.0 2019-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc eigensolver: "power"

 13:    Method: Power Iteration

 15:    Algorithm:

 17:        This solver implements the power iteration for finding dominant
 18:        eigenpairs. It also includes the following well-known methods:
 19:        - Inverse Iteration: when used in combination with shift-and-invert
 20:          spectral transformation.
 21:        - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
 22:          a variable shift.

 24:        It can also be used for nonlinear inverse iteration on the problem
 25:        A(x)*x=lambda*B(x)*x, where A and B are not constant but depend on x.

 27:    References:

 29:        [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report
 30:            STR-2, available at http://slepc.upv.es.
 31: */

 33: #include <slepc/private/epsimpl.h>                /*I "slepceps.h" I*/
 34: #include <slepcblaslapack.h>
 35: /* petsc headers */
 36: #include <petscdm.h>
 37: #include <petscsnes.h>

 39: static PetscErrorCode EPSPowerFormFunction_Update(SNES,Vec,Vec,void*);
 40: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx);
 41: PetscErrorCode EPSSolve_Power(EPS);
 42: PetscErrorCode EPSSolve_TS_Power(EPS);

 44: typedef struct {
 45:   EPSPowerShiftType shift_type;
 46:   PetscBool         nonlinear;
 47:   PetscBool         update;
 48:   SNES              snes;
 49:   PetscErrorCode    (*formFunctionB)(SNES,Vec,Vec,void*);
 50:   void              *formFunctionBctx;
 51:   PetscErrorCode    (*formFunctionA)(SNES,Vec,Vec,void*);
 52:   void              *formFunctionActx;
 53:   PetscErrorCode    (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
 54:   PetscInt          idx;  /* index of the first nonzero entry in the iteration vector */
 55:   PetscInt          p;    /* process id of the owner of idx */
 56: } EPS_POWER;

 58: PetscErrorCode EPSSetUp_Power(EPS eps)
 59: {
 61:   EPS_POWER      *power = (EPS_POWER*)eps->data;
 62:   PetscBool      flg,istrivial;
 63:   STMatMode      mode;
 64:   Mat            A,B;
 65:   Vec            res;
 66:   PetscContainer container;
 67:   PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
 68:   PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
 69:   void           *ctx;
 70:   SNESLineSearch linesearch;

 73:   if (eps->ncv) {
 74:     if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
 75:   } else eps->ncv = eps->nev;
 76:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 77:   if (!eps->max_it) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
 78:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 79:   if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 80:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
 81:     if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
 82:     PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
 83:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
 84:     STGetMatMode(eps->st,&mode);
 85:     if (mode == ST_MATMODE_INPLACE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
 86:   }
 87:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 88:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
 89:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 90:   RGIsTrivial(eps->rg,&istrivial);
 91:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
 92:   EPSAllocateSolution(eps,0);
 93:   EPS_SetInnerProduct(eps);

 95:   if (power->nonlinear) {
 96:     if (eps->nev>1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
 97:     EPSSetWorkVecs(eps,4);

 99:     /* set up SNES solver */
100:     if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
101:     EPSGetOperators(eps,&A,&B);

103:     PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA);
104:     if (!formFunctionA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
105:     PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container);
106:     if (container) {
107:       PetscContainerGetPointer(container,&ctx);
108:     } else ctx = NULL;
109:     MatCreateVecs(A,&res,NULL);
110:     power->formFunctionA = formFunctionA;
111:     power->formFunctionActx = ctx;
112:     if (power->update) {
113:       SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx);
114:       PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB);
115:     }
116:     else { SNESSetFunction(power->snes,res,formFunctionA,ctx); }
117:     VecDestroy(&res);

119:     PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA);
120:     if (!formJacobianA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
121:     PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container);
122:     if (container) {
123:       PetscContainerGetPointer(container,&ctx);
124:     } else ctx = NULL;
125:     SNESSetJacobian(power->snes,A,A,formJacobianA,ctx);
126:     SNESSetFromOptions(power->snes);
127:     SNESGetLineSearch(power->snes,&linesearch);
128:     if (power->update) {
129:       SNESLineSearchSetPostCheck(linesearch,SNESLineSearchPostheckFunction,ctx);
130:     }
131:     SNESSetUp(power->snes);
132:     if (B) {
133:       PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB);
134:       PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container);
135:       if (power->formFunctionB && container) {
136:         PetscContainerGetPointer(container,&power->formFunctionBctx);
137:       } else power->formFunctionBctx = NULL;
138:     }
139:   } else {
140:     if (eps->twosided) {
141:       EPSSetWorkVecs(eps,3);
142:     } else {
143:       EPSSetWorkVecs(eps,2);
144:     }
145:     DSSetType(eps->ds,DSNHEP);
146:     DSAllocate(eps->ds,eps->nev);
147:   }
148:   /* dispatch solve method */
149:   if (eps->twosided) {
150:     if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration does not have two-sided variant");
151:     if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided variant does not support Wilkinson shifts");
152:     eps->ops->solve = EPSSolve_TS_Power;
153:   } else eps->ops->solve = EPSSolve_Power;
154:   return(0);
155: }

157: /*
158:    Find the index of the first nonzero entry of x, and the process that owns it.
159: */
160: static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscInt *p)
161: {
162:   PetscErrorCode    ierr;
163:   PetscInt          i,first,last,N;
164:   PetscLayout       map;
165:   const PetscScalar *xx;

168:   VecGetSize(x,&N);
169:   VecGetOwnershipRange(x,&first,&last);
170:   VecGetArrayRead(x,&xx);
171:   for (i=first;i<last;i++) {
172:     if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
173:   }
174:   VecRestoreArrayRead(x,&xx);
175:   if (i==last) i=N;
176:   MPI_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x));
177:   if (*idx==N) SETERRQ(PetscObjectComm((PetscObject)x),1,"Zero vector found");
178:   VecGetLayout(x,&map);
179:   PetscLayoutFindOwner(map,*idx,p);
180:   return(0);
181: }

183: /*
184:    Normalize a vector x with respect to a given norm as well as the
185:    sign of the first nonzero entry (assumed to be idx in process p).
186: */
187: static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscInt p,PetscScalar *sign)
188: {
189:   PetscErrorCode    ierr;
190:   PetscScalar       alpha=1.0;
191:   PetscInt          first,last;
192:   PetscReal         absv;
193:   const PetscScalar *xx;

196:   VecGetOwnershipRange(x,&first,&last);
197:   if (idx>=first && idx<last) {
198:     VecGetArrayRead(x,&xx);
199:     absv = PetscAbsScalar(xx[idx-first]);
200:     if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
201:     VecRestoreArrayRead(x,&xx);
202:   }
203:   MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x));
204:   if (sign) *sign = alpha;
205:   alpha *= norm;
206:   VecScale(x,1.0/alpha);
207:   return(0);
208: }

210: static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
211: {
213:   EPS_POWER      *power = (EPS_POWER*)eps->data;
214:   Mat            B;

217:   STResetMatrixState(eps->st);
218:   EPSGetOperators(eps,NULL,&B);
219:   if (B) {
220:     if (power->formFunctionB) {
221:       (*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx);
222:     } else {
223:       MatMult(B,x,Bx);
224:     }
225:   } else {
226:     VecCopy(x,Bx);
227:   }
228:   return(0);
229: }

231: static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
232: {
234:   EPS_POWER      *power = (EPS_POWER*)eps->data;
235:   Mat            A;

238:   STResetMatrixState(eps->st);
239:   EPSGetOperators(eps,&A,NULL);
240:   if (A) {
241:     if (power->formFunctionA) {
242:       (*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx);
243:     } else {
244:       MatMult(A,x,Ax);
245:     }
246:   } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem");
247:   return(0);
248: }

250: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
251: {
253:   SNES           snes;
254:   EPS            eps;
255:   Vec            oldx;

258:   SNESLineSearchGetSNES(linesearch,&snes);
259:   PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
260:   if (!eps) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
261:   oldx = eps->work[3];
262:   VecCopy(x,oldx);
263:   return(0);
264: }

266: static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
267: {
269:   EPS            eps;
270:   PetscReal      bx;
271:   Vec            Bx;
272:   EPS_POWER      *power;

275:   PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
276:   if (!eps) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
277:   power = (EPS_POWER*)eps->data;
278:   Bx = eps->work[2];
279:   if (power->formFunctionAB) {
280:     (*power->formFunctionAB)(snes,x,y,Bx,ctx);
281:   } else {
282:     EPSPowerUpdateFunctionA(eps,x,y);
283:     EPSPowerUpdateFunctionB(eps,x,Bx);
284:   }
285:   VecNorm(Bx,NORM_2,&bx);
286:   Normalize(Bx,bx,power->idx,power->p,NULL);
287:   VecAXPY(y,-1.0,Bx);
288:   return(0);
289: }

291: /*
292:    Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
293: */
294: static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
295: {
297:   EPS_POWER      *power = (EPS_POWER*)eps->data;
298:   Vec            Bx;

301:   VecCopy(x,y);
302:   if (power->update) {
303:     SNESSolve(power->snes,NULL,y);
304:   } else {
305:     Bx = eps->work[2];
306:     SNESSolve(power->snes,Bx,y);
307:   }
308:   return(0);
309: }

311: /*
312:    Use nonlinear inverse power to compute an initial guess
313: */
314: static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
315: {
316:   EPS            powereps;
317:   Mat            A,B;
318:   Vec            v1,v2;
319:   SNES           snes;
320:   DM             dm,newdm;

324:   EPSCreate(PetscObjectComm((PetscObject)eps),&powereps);
325:   EPSGetOperators(eps,&A,&B);
326:   EPSSetType(powereps,EPSPOWER);
327:   EPSSetOperators(powereps,A,B);
328:   EPSSetTolerances(powereps,1e-6,4);
329:   EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix);
330:   EPSAppendOptionsPrefix(powereps,"init_");
331:   EPSSetProblemType(powereps,EPS_GNHEP);
332:   EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE);
333:   EPSPowerSetNonlinear(powereps,PETSC_TRUE);
334:   STSetType(powereps->st,STSINVERT);
335:   /* attach dm to initial solve */
336:   EPSPowerGetSNES(eps,&snes);
337:   SNESGetDM(snes,&dm);
338:   /* use  dmshell to temporarily store snes context */
339:   DMCreate(PetscObjectComm((PetscObject)eps),&newdm);
340:   DMSetType(newdm,DMSHELL);
341:   DMSetUp(newdm);
342:   DMCopyDMSNES(dm,newdm);
343:   EPSPowerGetSNES(powereps,&snes);
344:   SNESSetDM(snes,dm);
345:   EPSSetFromOptions(powereps);
346:   EPSSolve(powereps);
347:   BVGetColumn(eps->V,0,&v2);
348:   BVGetColumn(powereps->V,0,&v1);
349:   VecCopy(v1,v2);
350:   BVRestoreColumn(powereps->V,0,&v1);
351:   BVRestoreColumn(eps->V,0,&v2);
352:   EPSDestroy(&powereps);
353:   /* restore context back to the old nonlinear solver */
354:   DMCopyDMSNES(newdm,dm);
355:   DMDestroy(&newdm);
356:   return(0);
357: }

359: PetscErrorCode EPSSolve_Power(EPS eps)
360: {
361: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
363:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
364: #else
365:   PetscErrorCode      ierr;
366:   EPS_POWER           *power = (EPS_POWER*)eps->data;
367:   PetscInt            k,ld;
368:   Vec                 v,y,e,Bx;
369:   Mat                 A;
370:   KSP                 ksp;
371:   PetscReal           relerr,norm,norm1,rt1,rt2,cs1;
372:   PetscScalar         theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
373:   PetscBool           breakdown;
374:   KSPConvergedReason  reason;
375:   SNESConvergedReason snesreason;

378:   e = eps->work[0];
379:   y = eps->work[1];
380:   if (power->nonlinear) Bx = eps->work[2];
381:   else Bx = NULL;

383:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
384:   if (power->nonlinear) {
385:     PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps);
386:     if (power->update) {
387:       EPSPowerComputeInitialGuess_Update(eps);
388:     }
389:   } else {
390:     DSGetLeadingDimension(eps->ds,&ld);
391:   }
392:   if (!power->update) {
393:     EPSGetStartVector(eps,0,NULL);
394:   }
395:   if (power->nonlinear) {
396:     BVGetColumn(eps->V,0,&v);
397:     EPSPowerUpdateFunctionB(eps,v,Bx);
398:     VecNorm(Bx,NORM_2,&norm);
399:     FirstNonzeroIdx(Bx,&power->idx,&power->p);
400:     Normalize(Bx,norm,power->idx,power->p,NULL);
401:     BVRestoreColumn(eps->V,0,&v);
402:   }

404:   STGetShift(eps->st,&sigma);    /* original shift */
405:   rho = sigma;

407:   while (eps->reason == EPS_CONVERGED_ITERATING) {
408:     eps->its++;
409:     k = eps->nconv;

411:     /* y = OP v */
412:     BVGetColumn(eps->V,k,&v);
413:     if (power->nonlinear) {
414:       VecCopy(v,eps->work[3]);
415:       EPSPowerApply_SNES(eps,v,y);
416:       VecCopy(eps->work[3],v);
417:     } else {
418:       STApply(eps->st,v,y);
419:     }
420:     BVRestoreColumn(eps->V,k,&v);

422:     /* purge previously converged eigenvectors */
423:     if (power->nonlinear) {
424:       EPSPowerUpdateFunctionB(eps,y,Bx);
425:       VecNorm(Bx,NORM_2,&norm);
426:       Normalize(Bx,norm,power->idx,power->p,&sign);
427:     } else {
428:       DSGetArray(eps->ds,DS_MAT_A,&T);
429:       BVSetActiveColumns(eps->V,0,k);
430:       BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL);
431:     }

433:     /* theta = (v,y)_B */
434:     BVSetActiveColumns(eps->V,k,k+1);
435:     BVDotVec(eps->V,y,&theta);
436:     if (!power->nonlinear) {
437:       T[k+k*ld] = theta;
438:       DSRestoreArray(eps->ds,DS_MAT_A,&T);
439:     }

441:     if (power->nonlinear) theta = norm*sign;

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

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

448:       /* compute relative error as ||y-theta v||_2/|theta| */
449:       VecCopy(y,e);
450:       BVGetColumn(eps->V,k,&v);
451:       VecAXPY(e,power->nonlinear?-1.0:-theta,v);
452:       BVRestoreColumn(eps->V,k,&v);
453:       VecNorm(e,NORM_2,&relerr);
454:       relerr /= PetscAbsScalar(theta);

456:     } else {  /* RQI */

458:       /* delta = ||y||_B */
459:       delta = norm;

461:       /* compute relative error */
462:       if (rho == 0.0) relerr = PETSC_MAX_REAL;
463:       else relerr = 1.0 / (norm*PetscAbsScalar(rho));

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

468:       /* compute new shift */
469:       if (relerr<eps->tol) {
470:         rho = sigma;  /* if converged, restore original shift */
471:         STSetShift(eps->st,rho);
472:       } else {
473:         rho = rho + PetscConj(theta)/(delta*delta);  /* Rayleigh quotient R(v) */
474:         if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
475:           /* beta1 is the norm of the residual associated with R(v) */
476:           BVGetColumn(eps->V,k,&v);
477:           VecAXPY(v,-PetscConj(theta)/(delta*delta),y);
478:           BVRestoreColumn(eps->V,k,&v);
479:           BVScaleColumn(eps->V,k,1.0/delta);
480:           BVNormColumn(eps->V,k,NORM_2,&norm1);
481:           beta1 = norm1;

483:           /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
484:           STGetMatrix(eps->st,0,&A);
485:           BVGetColumn(eps->V,k,&v);
486:           MatMult(A,v,e);
487:           VecDot(v,e,&alpha2);
488:           BVRestoreColumn(eps->V,k,&v);
489:           alpha2 = alpha2 / (beta1 * beta1);

491:           /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
492:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
493:           PetscStackCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
494:           PetscFPTrapPop();
495:           if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
496:           else rho = rt2;
497:         }
498:         /* update operator according to new shift */
499:         KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
500:         STSetShift(eps->st,rho);
501:         KSPGetConvergedReason(ksp,&reason);
502:         if (reason) {
503:           PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
504:           rho *= 1+10*PETSC_MACHINE_EPSILON;
505:           STSetShift(eps->st,rho);
506:           KSPGetConvergedReason(ksp,&reason);
507:           if (reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
508:         }
509:         KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
510:       }
511:     }
512:     eps->errest[eps->nconv] = relerr;

514:     /* normalize */
515:     if (!power->nonlinear) { Normalize(y,norm,power->idx,power->p,NULL); }
516:     BVInsertVec(eps->V,k,y);

518:     if (power->update) {
519:       SNESGetConvergedReason(power->snes,&snesreason);
520:     }
521:     /* if relerr<tol, accept eigenpair */
522:     if (relerr<eps->tol || (power->update && snesreason > 0)) {
523:       eps->nconv = eps->nconv + 1;
524:       if (eps->nconv<eps->nev) {
525:         EPSGetStartVector(eps,eps->nconv,&breakdown);
526:         if (breakdown) {
527:           eps->reason = EPS_DIVERGED_BREAKDOWN;
528:           PetscInfo(eps,"Unable to generate more start vectors\n");
529:           break;
530:         }
531:       }
532:     }
533:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
534:     (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
535:   }

537:   if (power->nonlinear) {
538:     PetscObjectCompose((PetscObject)power->snes,"eps",NULL);
539:   } else {
540:     DSSetDimensions(eps->ds,eps->nconv,0,0,0);
541:     DSSetState(eps->ds,DS_STATE_RAW);
542:   }
543:   return(0);
544: #endif
545: }

547: PetscErrorCode EPSSolve_TS_Power(EPS eps)
548: {
549:   PetscErrorCode     ierr;
550:   EPS_POWER          *power = (EPS_POWER*)eps->data;
551:   PetscInt           k,ld;
552:   Vec                v,w,y,e,z;
553:   KSP                ksp;
554:   PetscReal          relerr,relerrl,delta;
555:   PetscScalar        theta,rho,alpha,sigma;
556:   PetscBool          breakdown;
557:   KSPConvergedReason reason;

560:   e = eps->work[0];
561:   v = eps->work[1];
562:   w = eps->work[2];

564:   if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
565:   DSGetLeadingDimension(eps->ds,&ld);
566:   EPSGetStartVector(eps,0,NULL);
567:   BVSetRandomColumn(eps->W,0);
568:   BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL);
569:   BVCopyVec(eps->V,0,v);
570:   BVCopyVec(eps->W,0,w);
571:   STGetShift(eps->st,&sigma);    /* original shift */
572:   rho = sigma;

574:   while (eps->reason == EPS_CONVERGED_ITERATING) {
575:     eps->its++;
576:     k = eps->nconv;

578:     /* y = OP v, z = OP' w */
579:     BVGetColumn(eps->V,k,&y);
580:     STApply(eps->st,v,y);
581:     BVRestoreColumn(eps->V,k,&y);
582:     BVGetColumn(eps->W,k,&z);
583:     VecConjugate(w);
584:     STApplyTranspose(eps->st,w,z);
585:     VecConjugate(w);
586:     VecConjugate(z);
587:     BVRestoreColumn(eps->W,k,&z);

589:     /* purge previously converged eigenvectors */
590:     BVBiorthogonalizeColumn(eps->V,eps->W,k);

592:     /* theta = (w,y)_B */
593:     BVSetActiveColumns(eps->V,k,k+1);
594:     BVDotVec(eps->V,w,&theta);
595:     theta = PetscConj(theta);

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

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

602:       /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
603:       BVCopyVec(eps->V,k,e);
604:       VecAXPY(e,-theta,v);
605:       VecNorm(e,NORM_2,&relerr);
606:       BVCopyVec(eps->W,k,e);
607:       VecAXPY(e,-PetscConj(theta),w);
608:       VecNorm(e,NORM_2,&relerrl);
609:       relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
610:     }

612:     /* normalize */
613:     BVSetActiveColumns(eps->V,k,k+1);
614:     BVGetColumn(eps->W,k,&z);
615:     BVDotVec(eps->V,z,&alpha);
616:     BVRestoreColumn(eps->W,k,&z);
617:     delta = PetscSqrtReal(PetscAbsScalar(alpha));
618:     if (delta==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Breakdown in two-sided Power/RQI");
619:     BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta));
620:     BVScaleColumn(eps->W,k,1.0/delta);
621:     BVCopyVec(eps->V,k,v);
622:     BVCopyVec(eps->W,k,w);

624:     if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */

626:       /* compute relative error */
627:       if (rho == 0.0) relerr = PETSC_MAX_REAL;
628:       else relerr = 1.0 / PetscAbsScalar(delta*rho);

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

633:       /* compute new shift */
634:       if (relerr<eps->tol) {
635:         rho = sigma;  /* if converged, restore original shift */
636:         STSetShift(eps->st,rho);
637:       } else {
638:         rho = rho + PetscConj(theta)/(delta*delta);  /* Rayleigh quotient R(v) */
639:         /* update operator according to new shift */
640:         KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
641:         STSetShift(eps->st,rho);
642:         KSPGetConvergedReason(ksp,&reason);
643:         if (reason) {
644:           PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
645:           rho *= 1+10*PETSC_MACHINE_EPSILON;
646:           STSetShift(eps->st,rho);
647:           KSPGetConvergedReason(ksp,&reason);
648:           if (reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
649:         }
650:         KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
651:       }
652:     }
653:     eps->errest[eps->nconv] = relerr;

655:     /* if relerr<tol, accept eigenpair */
656:     if (relerr<eps->tol) {
657:       eps->nconv = eps->nconv + 1;
658:       if (eps->nconv<eps->nev) {
659:         EPSGetStartVector(eps,eps->nconv,&breakdown);
660:         if (breakdown) {
661:           eps->reason = EPS_DIVERGED_BREAKDOWN;
662:           PetscInfo(eps,"Unable to generate more start vectors\n");
663:           break;
664:         }
665:         BVSetRandomColumn(eps->W,eps->nconv);
666:         BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL);
667:       }
668:     }
669:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
670:     (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
671:   }

673:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
674:   DSSetState(eps->ds,DS_STATE_RAW);
675:   return(0);
676: }

678: PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
679: {
681:   EPS_POWER      *power = (EPS_POWER*)eps->data;
682:   SNESConvergedReason snesreason;

685:   if (power->update) {
686:     SNESGetConvergedReason(power->snes,&snesreason);
687:     if (snesreason < 0) {
688:       *reason = EPS_DIVERGED_BREAKDOWN;
689:       return(0);
690:     }
691:   }
692:   EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx);
693:   return(0);
694: }

696: PetscErrorCode EPSBackTransform_Power(EPS eps)
697: {
699:   EPS_POWER      *power = (EPS_POWER*)eps->data;

702:   if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
703:   else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
704:     EPSBackTransform_Default(eps);
705:   }
706:   return(0);
707: }

709: PetscErrorCode EPSSetFromOptions_Power(PetscOptionItems *PetscOptionsObject,EPS eps)
710: {
711:   PetscErrorCode    ierr;
712:   EPS_POWER         *power = (EPS_POWER*)eps->data;
713:   PetscBool         flg,val;
714:   EPSPowerShiftType shift;

717:   PetscOptionsHead(PetscOptionsObject,"EPS Power Options");

719:     PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
720:     if (flg) { EPSPowerSetShiftType(eps,shift); }

722:     PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg);
723:     if (flg) { EPSPowerSetNonlinear(eps,val); }

725:     PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg);
726:     if (flg) { EPSPowerSetUpdate(eps,val); }

728:   PetscOptionsTail();
729:   return(0);
730: }

732: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
733: {
734:   EPS_POWER *power = (EPS_POWER*)eps->data;

737:   switch (shift) {
738:     case EPS_POWER_SHIFT_CONSTANT:
739:     case EPS_POWER_SHIFT_RAYLEIGH:
740:     case EPS_POWER_SHIFT_WILKINSON:
741:       if (power->shift_type != shift) {
742:         power->shift_type = shift;
743:         eps->state = EPS_STATE_INITIAL;
744:       }
745:       break;
746:     default:
747:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
748:   }
749:   return(0);
750: }

752: /*@
753:    EPSPowerSetShiftType - Sets the type of shifts used during the power
754:    iteration. This can be used to emulate the Rayleigh Quotient Iteration
755:    (RQI) method.

757:    Logically Collective on eps

759:    Input Parameters:
760: +  eps - the eigenproblem solver context
761: -  shift - the type of shift

763:    Options Database Key:
764: .  -eps_power_shift_type - Sets the shift type (either 'constant' or
765:                            'rayleigh' or 'wilkinson')

767:    Notes:
768:    By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
769:    is the simple power method (or inverse iteration if a shift-and-invert
770:    transformation is being used).

772:    A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
773:    EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
774:    a cubic converging method such as RQI.

776:    Level: advanced

778: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
779: @*/
780: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
781: {

787:   PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
788:   return(0);
789: }

791: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
792: {
793:   EPS_POWER *power = (EPS_POWER*)eps->data;

796:   *shift = power->shift_type;
797:   return(0);
798: }

800: /*@
801:    EPSPowerGetShiftType - Gets the type of shifts used during the power
802:    iteration.

804:    Not Collective

806:    Input Parameter:
807: .  eps - the eigenproblem solver context

809:    Input Parameter:
810: .  shift - the type of shift

812:    Level: advanced

814: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
815: @*/
816: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
817: {

823:   PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
824:   return(0);
825: }

827: static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
828: {
829:   EPS_POWER *power = (EPS_POWER*)eps->data;

832:   if (power->nonlinear != nonlinear) {
833:     power->nonlinear = nonlinear;
834:     eps->useds = PetscNot(nonlinear);
835:     eps->state = EPS_STATE_INITIAL;
836:   }
837:   return(0);
838: }

840: /*@
841:    EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.

843:    Logically Collective on eps

845:    Input Parameters:
846: +  eps - the eigenproblem solver context
847: -  nonlinear - whether the problem is nonlinear or not

849:    Options Database Key:
850: .  -eps_power_nonlinear - Sets the nonlinear flag

852:    Notes:
853:    If this flag is set, the solver assumes that the problem is nonlinear,
854:    that is, the operators that define the eigenproblem are not constant
855:    matrices, but depend on the eigenvector: A(x)*x=lambda*B(x)*x. This is
856:    different from the case of nonlinearity with respect to the eigenvalue
857:    (use the NEP solver class for this kind of problems).

859:    The way in which nonlinear operators are specified is very similar to
860:    the case of PETSc's SNES solver. The difference is that the callback
861:    functions are provided via composed functions "formFunction" and
862:    "formJacobian" in each of the matrix objects passed as arguments of
863:    EPSSetOperators(). The application context required for these functions
864:    can be attached via a composed PetscContainer.

866:    Level: advanced

868: .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
869: @*/
870: PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
871: {

877:   PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
878:   return(0);
879: }

881: static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
882: {
883:   EPS_POWER *power = (EPS_POWER*)eps->data;

886:   *nonlinear = power->nonlinear;
887:   return(0);
888: }

890: /*@
891:    EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.

893:    Not Collective

895:    Input Parameter:
896: .  eps - the eigenproblem solver context

898:    Input Parameter:
899: .  nonlinear - the nonlinear flag

901:    Level: advanced

903: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
904: @*/
905: PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
906: {

912:   PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
913:   return(0);
914: }

916: static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
917: {
918:   EPS_POWER *power = (EPS_POWER*)eps->data;

921:   if (!power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems");
922:   power->update = update;
923:   eps->state = EPS_STATE_INITIAL;
924:   return(0);
925: }

927: /*@
928:    EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
929:    for nonlinear problems. This potentially has a better convergence.

931:    Logically Collective on eps

933:    Input Parameters:
934: +  eps - the eigenproblem solver context
935: -  update - whether the residual is updated monolithically or not

937:    Options Database Key:
938: .  -eps_power_update - Sets the update flag

940:    Level: advanced

942: .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
943: @*/
944: PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
945: {

951:   PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
952:   return(0);
953: }

955: static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
956: {
957:   EPS_POWER *power = (EPS_POWER*)eps->data;

960:   *update = power->update;
961:   return(0);
962: }

964: /*@
965:    EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
966:    for nonlinear problems.

968:    Not Collective

970:    Input Parameter:
971: .  eps - the eigenproblem solver context

973:    Input Parameter:
974: .  update - the update flag

976:    Level: advanced

978: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
979: @*/
980: PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
981: {

987:   PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
988:   return(0);
989: }

991: static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
992: {
994:   EPS_POWER      *power = (EPS_POWER*)eps->data;

997:   PetscObjectReference((PetscObject)snes);
998:   SNESDestroy(&power->snes);
999:   power->snes = snes;
1000:   PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
1001:   eps->state = EPS_STATE_INITIAL;
1002:   return(0);
1003: }

1005: /*@
1006:    EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
1007:    eigenvalue solver (to be used in nonlinear inverse iteration).

1009:    Collective on eps

1011:    Input Parameters:
1012: +  eps  - the eigenvalue solver
1013: -  snes - the nonlinear solver object

1015:    Level: advanced

1017: .seealso: EPSPowerGetSNES()
1018: @*/
1019: PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
1020: {

1027:   PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
1028:   return(0);
1029: }

1031: static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
1032: {
1034:   EPS_POWER      *power = (EPS_POWER*)eps->data;

1037:   if (!power->snes) {
1038:     SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes);
1039:     PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1);
1040:     SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix);
1041:     SNESAppendOptionsPrefix(power->snes,"eps_power_");
1042:     PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
1043:     PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options);
1044:   }
1045:   *snes = power->snes;
1046:   return(0);
1047: }

1049: /*@
1050:    EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
1051:    with the eigenvalue solver.

1053:    Not Collective

1055:    Input Parameter:
1056: .  eps - the eigenvalue solver

1058:    Output Parameter:
1059: .  snes - the nonlinear solver object

1061:    Level: advanced

1063: .seealso: EPSPowerSetSNES()
1064: @*/
1065: PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
1066: {

1072:   PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
1073:   return(0);
1074: }

1076: PetscErrorCode EPSReset_Power(EPS eps)
1077: {
1079:   EPS_POWER      *power = (EPS_POWER*)eps->data;

1082:   if (power->snes) { SNESReset(power->snes); }
1083:   return(0);
1084: }

1086: PetscErrorCode EPSDestroy_Power(EPS eps)
1087: {
1089:   EPS_POWER      *power = (EPS_POWER*)eps->data;

1092:   if (power->nonlinear) {
1093:     SNESDestroy(&power->snes);
1094:   }
1095:   PetscFree(eps->data);
1096:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
1097:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
1098:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL);
1099:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL);
1100:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL);
1101:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL);
1102:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL);
1103:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL);
1104:   return(0);
1105: }

1107: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
1108: {
1110:   EPS_POWER      *power = (EPS_POWER*)eps->data;
1111:   PetscBool      isascii;

1114:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1115:   if (isascii) {
1116:     if (power->nonlinear) {
1117:       PetscViewerASCIIPrintf(viewer,"  using nonlinear inverse iteration\n");
1118:       if (power->update) {
1119:         PetscViewerASCIIPrintf(viewer,"  updating the residual monolithically\n");
1120:       }
1121:       if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
1122:       PetscViewerASCIIPushTab(viewer);
1123:       SNESView(power->snes,viewer);
1124:       PetscViewerASCIIPopTab(viewer);
1125:     } else {
1126:       PetscViewerASCIIPrintf(viewer,"  %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
1127:     }
1128:   }
1129:   return(0);
1130: }

1132: PetscErrorCode EPSComputeVectors_Power(EPS eps)
1133: {
1135:   EPS_POWER      *power = (EPS_POWER*)eps->data;
1136:   PetscReal      norm;
1137:   PetscInt       i;

1140:   if (eps->twosided) {
1141:     EPSComputeVectors_Twosided(eps);
1142:     /* normalize (no need to take care of 2x2 blocks */
1143:     for (i=0;i<eps->nconv;i++) {
1144:       BVNormColumn(eps->V,i,NORM_2,&norm);
1145:       BVScaleColumn(eps->V,i,1.0/norm);
1146:       BVNormColumn(eps->W,i,NORM_2,&norm);
1147:       BVScaleColumn(eps->W,i,1.0/norm);
1148:     }
1149:   } else if (!power->nonlinear) {
1150:     EPSComputeVectors_Schur(eps);
1151:   }
1152:   return(0);
1153: }

1155: PetscErrorCode EPSSetDefaultST_Power(EPS eps)
1156: {
1158:   EPS_POWER      *power = (EPS_POWER*)eps->data;
1159:   KSP            ksp;
1160:   PC             pc;

1163:   if (power->nonlinear) {
1164:     STGetKSP(eps->st,&ksp);
1165:     KSPGetPC(ksp,&pc);
1166:     PCSetType(pc,PCNONE);
1167:   }
1168:   return(0);
1169: }

1171: SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
1172: {
1173:   EPS_POWER      *ctx;

1177:   PetscNewLog(eps,&ctx);
1178:   eps->data = (void*)ctx;

1180:   eps->useds = PETSC_TRUE;
1181:   eps->hasts = PETSC_TRUE;
1182:   eps->categ = EPS_CATEGORY_OTHER;

1184:   eps->ops->setup          = EPSSetUp_Power;
1185:   eps->ops->setfromoptions = EPSSetFromOptions_Power;
1186:   eps->ops->reset          = EPSReset_Power;
1187:   eps->ops->destroy        = EPSDestroy_Power;
1188:   eps->ops->view           = EPSView_Power;
1189:   eps->ops->backtransform  = EPSBackTransform_Power;
1190:   eps->ops->computevectors = EPSComputeVectors_Power;
1191:   eps->ops->setdefaultst   = EPSSetDefaultST_Power;
1192:   eps->stopping            = EPSStopping_Power;

1194:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
1195:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
1196:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power);
1197:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power);
1198:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power);
1199:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power);
1200:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power);
1201:   PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power);
1202:   return(0);
1203: }