Actual source code: power.c
slepc-3.13.3 2020-06-14
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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 https://slepc.upv.es.
31: */
33: #include <slepc/private/epsimpl.h>
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: PetscMPIInt p; /* process id of the owner of idx */
56: PetscReal norm0; /* norm of initial vector */
57: } EPS_POWER;
59: PetscErrorCode EPSSetUp_Power(EPS eps)
60: {
62: EPS_POWER *power = (EPS_POWER*)eps->data;
63: PetscBool flg,istrivial;
64: STMatMode mode;
65: Mat A,B,P;
66: Vec res;
67: PetscContainer container;
68: PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
69: PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
70: void *ctx;
71: SNESLineSearch linesearch;
74: if (eps->ncv!=PETSC_DEFAULT) {
75: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
76: } else eps->ncv = eps->nev;
77: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
78: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
79: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
80: if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
81: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
82: if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
83: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
84: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
85: STGetMatMode(eps->st,&mode);
86: if (mode == ST_MATMODE_INPLACE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
87: }
88: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
89: if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
90: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
91: RGIsTrivial(eps->rg,&istrivial);
92: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
93: EPSAllocateSolution(eps,0);
94: EPS_SetInnerProduct(eps);
96: if (power->nonlinear) {
97: if (eps->nev>1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
98: EPSSetWorkVecs(eps,4);
100: /* Set up SNES solver */
101: /* If SNES was setup, we need to reset it so that it will be consistent with the current EPS */
102: if (power->snes) { SNESReset(power->snes); }
103: else { EPSPowerGetSNES(eps,&power->snes); }
105: /* For nonlinear solver (Newton), we should scale the initial vector back.
106: The initial vector will be scaled in EPSSetUp. */
107: if (eps->IS) {
108: VecNorm((eps->IS)[0],NORM_2,&power->norm0);
109: }
111: EPSGetOperators(eps,&A,&B);
113: PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA);
114: if (!formFunctionA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
115: PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container);
116: if (container) {
117: PetscContainerGetPointer(container,&ctx);
118: } else ctx = NULL;
119: MatCreateVecs(A,&res,NULL);
120: power->formFunctionA = formFunctionA;
121: power->formFunctionActx = ctx;
122: if (power->update) {
123: SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx);
124: PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB);
125: }
126: else { SNESSetFunction(power->snes,res,formFunctionA,ctx); }
127: VecDestroy(&res);
129: PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA);
130: if (!formJacobianA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
131: PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container);
132: if (container) {
133: PetscContainerGetPointer(container,&ctx);
134: } else ctx = NULL;
135: /* This allows users to compute a different preconditioning matrix than A.
136: * It is useful when A and B are shell matrices.
137: */
138: STPrecondGetMatForPC(eps->st,&P);
139: /* If users do not provide a matrix, we simply use A */
140: SNESSetJacobian(power->snes,A,P? P:A,formJacobianA,ctx);
141: SNESSetFromOptions(power->snes);
142: SNESGetLineSearch(power->snes,&linesearch);
143: if (power->update) {
144: SNESLineSearchSetPostCheck(linesearch,SNESLineSearchPostheckFunction,ctx);
145: }
146: SNESSetUp(power->snes);
147: if (B) {
148: PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB);
149: PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container);
150: if (power->formFunctionB && container) {
151: PetscContainerGetPointer(container,&power->formFunctionBctx);
152: } else power->formFunctionBctx = NULL;
153: }
154: } else {
155: if (eps->twosided) {
156: EPSSetWorkVecs(eps,3);
157: } else {
158: EPSSetWorkVecs(eps,2);
159: }
160: DSSetType(eps->ds,DSNHEP);
161: DSAllocate(eps->ds,eps->nev);
162: }
163: /* dispatch solve method */
164: if (eps->twosided) {
165: if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration does not have two-sided variant");
166: if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided variant does not support Wilkinson shifts");
167: eps->ops->solve = EPSSolve_TS_Power;
168: } else eps->ops->solve = EPSSolve_Power;
169: return(0);
170: }
172: /*
173: Find the index of the first nonzero entry of x, and the process that owns it.
174: */
175: static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscMPIInt *p)
176: {
177: PetscErrorCode ierr;
178: PetscInt i,first,last,N;
179: PetscLayout map;
180: const PetscScalar *xx;
183: VecGetSize(x,&N);
184: VecGetOwnershipRange(x,&first,&last);
185: VecGetArrayRead(x,&xx);
186: for (i=first;i<last;i++) {
187: if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
188: }
189: VecRestoreArrayRead(x,&xx);
190: if (i==last) i=N;
191: MPI_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x));
192: if (*idx==N) SETERRQ(PetscObjectComm((PetscObject)x),1,"Zero vector found");
193: VecGetLayout(x,&map);
194: PetscLayoutFindOwner(map,*idx,p);
195: return(0);
196: }
198: /*
199: Normalize a vector x with respect to a given norm as well as the
200: sign of the first nonzero entry (assumed to be idx in process p).
201: */
202: static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscMPIInt p,PetscScalar *sign)
203: {
204: PetscErrorCode ierr;
205: PetscScalar alpha=1.0;
206: PetscInt first,last;
207: PetscReal absv;
208: const PetscScalar *xx;
211: VecGetOwnershipRange(x,&first,&last);
212: if (idx>=first && idx<last) {
213: VecGetArrayRead(x,&xx);
214: absv = PetscAbsScalar(xx[idx-first]);
215: if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
216: VecRestoreArrayRead(x,&xx);
217: }
218: MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x));
219: if (sign) *sign = alpha;
220: alpha *= norm;
221: VecScale(x,1.0/alpha);
222: return(0);
223: }
225: static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
226: {
228: EPS_POWER *power = (EPS_POWER*)eps->data;
229: Mat B;
232: STResetMatrixState(eps->st);
233: EPSGetOperators(eps,NULL,&B);
234: if (B) {
235: if (power->formFunctionB) {
236: (*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx);
237: } else {
238: MatMult(B,x,Bx);
239: }
240: } else {
241: VecCopy(x,Bx);
242: }
243: return(0);
244: }
246: static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
247: {
249: EPS_POWER *power = (EPS_POWER*)eps->data;
250: Mat A;
253: STResetMatrixState(eps->st);
254: EPSGetOperators(eps,&A,NULL);
255: if (A) {
256: if (power->formFunctionA) {
257: (*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx);
258: } else {
259: MatMult(A,x,Ax);
260: }
261: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem");
262: return(0);
263: }
265: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
266: {
268: SNES snes;
269: EPS eps;
270: Vec oldx;
273: SNESLineSearchGetSNES(linesearch,&snes);
274: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
275: if (!eps) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
276: oldx = eps->work[3];
277: VecCopy(x,oldx);
278: return(0);
279: }
281: static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
282: {
284: EPS eps;
285: PetscReal bx;
286: Vec Bx;
287: EPS_POWER *power;
290: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
291: if (!eps) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
292: power = (EPS_POWER*)eps->data;
293: Bx = eps->work[2];
294: if (power->formFunctionAB) {
295: (*power->formFunctionAB)(snes,x,y,Bx,ctx);
296: } else {
297: EPSPowerUpdateFunctionA(eps,x,y);
298: EPSPowerUpdateFunctionB(eps,x,Bx);
299: }
300: VecNorm(Bx,NORM_2,&bx);
301: Normalize(Bx,bx,power->idx,power->p,NULL);
302: VecAXPY(y,-1.0,Bx);
303: return(0);
304: }
306: /*
307: Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
308: */
309: static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
310: {
312: EPS_POWER *power = (EPS_POWER*)eps->data;
313: Vec Bx;
316: VecCopy(x,y);
317: if (power->update) {
318: SNESSolve(power->snes,NULL,y);
319: } else {
320: Bx = eps->work[2];
321: SNESSolve(power->snes,Bx,y);
322: }
323: return(0);
324: }
326: /*
327: Use nonlinear inverse power to compute an initial guess
328: */
329: static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
330: {
331: EPS powereps;
332: Mat A,B,P;
333: Vec v1,v2;
334: SNES snes;
335: DM dm,newdm;
339: EPSCreate(PetscObjectComm((PetscObject)eps),&powereps);
340: EPSGetOperators(eps,&A,&B);
341: EPSSetType(powereps,EPSPOWER);
342: EPSSetOperators(powereps,A,B);
343: EPSSetTolerances(powereps,1e-6,4);
344: EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix);
345: EPSAppendOptionsPrefix(powereps,"init_");
346: EPSSetProblemType(powereps,EPS_GNHEP);
347: EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE);
348: EPSPowerSetNonlinear(powereps,PETSC_TRUE);
349: STPrecondGetMatForPC(eps->st,&P);
350: /* attach dm to initial solve */
351: EPSPowerGetSNES(eps,&snes);
352: SNESGetDM(snes,&dm);
353: /* use dmshell to temporarily store snes context */
354: DMCreate(PetscObjectComm((PetscObject)eps),&newdm);
355: DMSetType(newdm,DMSHELL);
356: DMSetUp(newdm);
357: DMCopyDMSNES(dm,newdm);
358: EPSPowerGetSNES(powereps,&snes);
359: SNESSetDM(snes,dm);
360: EPSSetFromOptions(powereps);
361: if (P) {
362: STPrecondSetMatForPC(powereps->st,P);
363: }
364: EPSSolve(powereps);
365: BVGetColumn(eps->V,0,&v2);
366: BVGetColumn(powereps->V,0,&v1);
367: VecCopy(v1,v2);
368: BVRestoreColumn(powereps->V,0,&v1);
369: BVRestoreColumn(eps->V,0,&v2);
370: EPSDestroy(&powereps);
371: /* restore context back to the old nonlinear solver */
372: DMCopyDMSNES(newdm,dm);
373: DMDestroy(&newdm);
374: return(0);
375: }
377: PetscErrorCode EPSSolve_Power(EPS eps)
378: {
379: PetscErrorCode ierr;
380: EPS_POWER *power = (EPS_POWER*)eps->data;
381: PetscInt k,ld;
382: Vec v,y,e,Bx;
383: Mat A;
384: KSP ksp;
385: PetscReal relerr,norm,norm1,rt1,rt2,cs1;
386: PetscScalar theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
387: PetscBool breakdown;
388: KSPConvergedReason reason;
389: SNESConvergedReason snesreason;
392: e = eps->work[0];
393: y = eps->work[1];
394: if (power->nonlinear) Bx = eps->work[2];
395: else Bx = NULL;
397: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
398: if (power->nonlinear) {
399: PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps);
400: /* Compute an intial guess only when users do not provide one */
401: if (power->update && !eps->nini) {
402: EPSPowerComputeInitialGuess_Update(eps);
403: }
404: } else {
405: DSGetLeadingDimension(eps->ds,&ld);
406: }
407: if (!power->update) {
408: EPSGetStartVector(eps,0,NULL);
409: }
410: if (power->nonlinear) {
411: BVGetColumn(eps->V,0,&v);
412: if (eps->nini) {
413: /* We scale the initial vector back if the initial vector was provided by users */
414: VecScale(v,power->norm0);
415: }
416: EPSPowerUpdateFunctionB(eps,v,Bx);
417: VecNorm(Bx,NORM_2,&norm);
418: FirstNonzeroIdx(Bx,&power->idx,&power->p);
419: Normalize(Bx,norm,power->idx,power->p,NULL);
420: BVRestoreColumn(eps->V,0,&v);
421: }
423: STGetShift(eps->st,&sigma); /* original shift */
424: rho = sigma;
426: while (eps->reason == EPS_CONVERGED_ITERATING) {
427: eps->its++;
428: k = eps->nconv;
430: /* y = OP v */
431: BVGetColumn(eps->V,k,&v);
432: if (power->nonlinear) {
433: VecCopy(v,eps->work[3]);
434: EPSPowerApply_SNES(eps,v,y);
435: VecCopy(eps->work[3],v);
436: } else {
437: STApply(eps->st,v,y);
438: }
439: BVRestoreColumn(eps->V,k,&v);
441: /* purge previously converged eigenvectors */
442: if (power->nonlinear) {
443: EPSPowerUpdateFunctionB(eps,y,Bx);
444: VecNorm(Bx,NORM_2,&norm);
445: Normalize(Bx,norm,power->idx,power->p,&sign);
446: } else {
447: DSGetArray(eps->ds,DS_MAT_A,&T);
448: BVSetActiveColumns(eps->V,0,k);
449: BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL);
450: }
452: /* theta = (v,y)_B */
453: BVSetActiveColumns(eps->V,k,k+1);
454: BVDotVec(eps->V,y,&theta);
455: if (!power->nonlinear) {
456: T[k+k*ld] = theta;
457: DSRestoreArray(eps->ds,DS_MAT_A,&T);
458: }
460: /* Eigenvalue: 1/|Bx| */
461: if (power->nonlinear) theta = 1.0/norm*sign;
463: if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
465: /* approximate eigenvalue is the Rayleigh quotient */
466: eps->eigr[eps->nconv] = theta;
468: /* compute relative error as ||y-theta v||_2/|theta| */
469: VecCopy(y,e);
470: BVGetColumn(eps->V,k,&v);
471: VecAXPY(e,power->nonlinear?-1.0:-theta,v);
472: BVRestoreColumn(eps->V,k,&v);
473: VecNorm(e,NORM_2,&relerr);
474: if (power->nonlinear)
475: relerr *= PetscAbsScalar(theta);
476: else
477: relerr /= PetscAbsScalar(theta);
479: } else { /* RQI */
481: /* delta = ||y||_B */
482: delta = norm;
484: /* compute relative error */
485: if (rho == 0.0) relerr = PETSC_MAX_REAL;
486: else relerr = 1.0 / (norm*PetscAbsScalar(rho));
488: /* approximate eigenvalue is the shift */
489: eps->eigr[eps->nconv] = rho;
491: /* compute new shift */
492: if (relerr<eps->tol) {
493: rho = sigma; /* if converged, restore original shift */
494: STSetShift(eps->st,rho);
495: } else {
496: rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
497: if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
498: /* beta1 is the norm of the residual associated with R(v) */
499: BVGetColumn(eps->V,k,&v);
500: VecAXPY(v,-PetscConj(theta)/(delta*delta),y);
501: BVRestoreColumn(eps->V,k,&v);
502: BVScaleColumn(eps->V,k,1.0/delta);
503: BVNormColumn(eps->V,k,NORM_2,&norm1);
504: beta1 = norm1;
506: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
507: STGetMatrix(eps->st,0,&A);
508: BVGetColumn(eps->V,k,&v);
509: MatMult(A,v,e);
510: VecDot(v,e,&alpha2);
511: BVRestoreColumn(eps->V,k,&v);
512: alpha2 = alpha2 / (beta1 * beta1);
514: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
515: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
516: PetscStackCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
517: PetscFPTrapPop();
518: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
519: else rho = rt2;
520: }
521: /* update operator according to new shift */
522: KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
523: STSetShift(eps->st,rho);
524: KSPGetConvergedReason(ksp,&reason);
525: if (reason) {
526: PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
527: rho *= 1+10*PETSC_MACHINE_EPSILON;
528: STSetShift(eps->st,rho);
529: KSPGetConvergedReason(ksp,&reason);
530: if (reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
531: }
532: KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
533: }
534: }
535: eps->errest[eps->nconv] = relerr;
537: /* normalize */
538: if (!power->nonlinear) { Normalize(y,norm,power->idx,power->p,NULL); }
539: BVInsertVec(eps->V,k,y);
541: if (power->update) {
542: SNESGetConvergedReason(power->snes,&snesreason);
543: }
544: /* if relerr<tol, accept eigenpair */
545: if (relerr<eps->tol || (power->update && snesreason > 0)) {
546: eps->nconv = eps->nconv + 1;
547: if (eps->nconv<eps->nev) {
548: EPSGetStartVector(eps,eps->nconv,&breakdown);
549: if (breakdown) {
550: eps->reason = EPS_DIVERGED_BREAKDOWN;
551: PetscInfo(eps,"Unable to generate more start vectors\n");
552: break;
553: }
554: }
555: }
556: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
557: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
558: }
560: if (power->nonlinear) {
561: PetscObjectCompose((PetscObject)power->snes,"eps",NULL);
562: } else {
563: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
564: DSSetState(eps->ds,DS_STATE_RAW);
565: }
566: return(0);
567: }
569: PetscErrorCode EPSSolve_TS_Power(EPS eps)
570: {
571: PetscErrorCode ierr;
572: EPS_POWER *power = (EPS_POWER*)eps->data;
573: PetscInt k,ld;
574: Vec v,w,y,e,z;
575: KSP ksp;
576: PetscReal relerr=1.0,relerrl,delta;
577: PetscScalar theta,rho,alpha,sigma;
578: PetscBool breakdown,breakdownl;
579: KSPConvergedReason reason;
582: e = eps->work[0];
583: v = eps->work[1];
584: w = eps->work[2];
586: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
587: DSGetLeadingDimension(eps->ds,&ld);
588: EPSGetStartVector(eps,0,NULL);
589: EPSGetLeftStartVector(eps,0,NULL);
590: BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL);
591: BVCopyVec(eps->V,0,v);
592: BVCopyVec(eps->W,0,w);
593: STGetShift(eps->st,&sigma); /* original shift */
594: rho = sigma;
596: while (eps->reason == EPS_CONVERGED_ITERATING) {
597: eps->its++;
598: k = eps->nconv;
600: /* y = OP v, z = OP' w */
601: BVGetColumn(eps->V,k,&y);
602: STApply(eps->st,v,y);
603: BVRestoreColumn(eps->V,k,&y);
604: BVGetColumn(eps->W,k,&z);
605: STApplyHermitianTranspose(eps->st,w,z);
606: BVRestoreColumn(eps->W,k,&z);
608: /* purge previously converged eigenvectors */
609: BVBiorthogonalizeColumn(eps->V,eps->W,k);
611: /* theta = (w,y)_B */
612: BVSetActiveColumns(eps->V,k,k+1);
613: BVDotVec(eps->V,w,&theta);
614: theta = PetscConj(theta);
616: if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
618: /* approximate eigenvalue is the Rayleigh quotient */
619: eps->eigr[eps->nconv] = theta;
621: /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
622: BVCopyVec(eps->V,k,e);
623: VecAXPY(e,-theta,v);
624: VecNorm(e,NORM_2,&relerr);
625: BVCopyVec(eps->W,k,e);
626: VecAXPY(e,-PetscConj(theta),w);
627: VecNorm(e,NORM_2,&relerrl);
628: relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
629: }
631: /* normalize */
632: BVSetActiveColumns(eps->V,k,k+1);
633: BVGetColumn(eps->W,k,&z);
634: BVDotVec(eps->V,z,&alpha);
635: BVRestoreColumn(eps->W,k,&z);
636: delta = PetscSqrtReal(PetscAbsScalar(alpha));
637: if (delta==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Breakdown in two-sided Power/RQI");
638: BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta));
639: BVScaleColumn(eps->W,k,1.0/delta);
640: BVCopyVec(eps->V,k,v);
641: BVCopyVec(eps->W,k,w);
643: if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */
645: /* compute relative error */
646: if (rho == 0.0) relerr = PETSC_MAX_REAL;
647: else relerr = 1.0 / PetscAbsScalar(delta*rho);
649: /* approximate eigenvalue is the shift */
650: eps->eigr[eps->nconv] = rho;
652: /* compute new shift */
653: if (relerr<eps->tol) {
654: rho = sigma; /* if converged, restore original shift */
655: STSetShift(eps->st,rho);
656: } else {
657: rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
658: /* update operator according to new shift */
659: KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
660: STSetShift(eps->st,rho);
661: KSPGetConvergedReason(ksp,&reason);
662: if (reason) {
663: PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
664: rho *= 1+10*PETSC_MACHINE_EPSILON;
665: STSetShift(eps->st,rho);
666: KSPGetConvergedReason(ksp,&reason);
667: if (reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
668: }
669: KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
670: }
671: }
672: eps->errest[eps->nconv] = relerr;
674: /* if relerr<tol, accept eigenpair */
675: if (relerr<eps->tol) {
676: eps->nconv = eps->nconv + 1;
677: if (eps->nconv<eps->nev) {
678: EPSGetStartVector(eps,eps->nconv,&breakdown);
679: EPSGetLeftStartVector(eps,eps->nconv,&breakdownl);
680: if (breakdown || breakdownl) {
681: eps->reason = EPS_DIVERGED_BREAKDOWN;
682: PetscInfo(eps,"Unable to generate more start vectors\n");
683: break;
684: }
685: BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL);
686: }
687: }
688: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
689: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
690: }
692: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
693: DSSetState(eps->ds,DS_STATE_RAW);
694: return(0);
695: }
697: PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
698: {
700: EPS_POWER *power = (EPS_POWER*)eps->data;
701: SNESConvergedReason snesreason;
704: if (power->update) {
705: SNESGetConvergedReason(power->snes,&snesreason);
706: if (snesreason < 0) {
707: *reason = EPS_DIVERGED_BREAKDOWN;
708: return(0);
709: }
710: }
711: EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx);
712: return(0);
713: }
715: PetscErrorCode EPSBackTransform_Power(EPS eps)
716: {
718: EPS_POWER *power = (EPS_POWER*)eps->data;
721: if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
722: else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
723: EPSBackTransform_Default(eps);
724: }
725: return(0);
726: }
728: PetscErrorCode EPSSetFromOptions_Power(PetscOptionItems *PetscOptionsObject,EPS eps)
729: {
730: PetscErrorCode ierr;
731: EPS_POWER *power = (EPS_POWER*)eps->data;
732: PetscBool flg,val;
733: EPSPowerShiftType shift;
736: PetscOptionsHead(PetscOptionsObject,"EPS Power Options");
738: PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
739: if (flg) { EPSPowerSetShiftType(eps,shift); }
741: PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg);
742: if (flg) { EPSPowerSetNonlinear(eps,val); }
744: PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg);
745: if (flg) { EPSPowerSetUpdate(eps,val); }
747: PetscOptionsTail();
748: return(0);
749: }
751: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
752: {
753: EPS_POWER *power = (EPS_POWER*)eps->data;
756: switch (shift) {
757: case EPS_POWER_SHIFT_CONSTANT:
758: case EPS_POWER_SHIFT_RAYLEIGH:
759: case EPS_POWER_SHIFT_WILKINSON:
760: if (power->shift_type != shift) {
761: power->shift_type = shift;
762: eps->state = EPS_STATE_INITIAL;
763: }
764: break;
765: default:
766: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
767: }
768: return(0);
769: }
771: /*@
772: EPSPowerSetShiftType - Sets the type of shifts used during the power
773: iteration. This can be used to emulate the Rayleigh Quotient Iteration
774: (RQI) method.
776: Logically Collective on eps
778: Input Parameters:
779: + eps - the eigenproblem solver context
780: - shift - the type of shift
782: Options Database Key:
783: . -eps_power_shift_type - Sets the shift type (either 'constant' or
784: 'rayleigh' or 'wilkinson')
786: Notes:
787: By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
788: is the simple power method (or inverse iteration if a shift-and-invert
789: transformation is being used).
791: A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
792: EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
793: a cubic converging method such as RQI.
795: Level: advanced
797: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
798: @*/
799: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
800: {
806: PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
807: return(0);
808: }
810: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
811: {
812: EPS_POWER *power = (EPS_POWER*)eps->data;
815: *shift = power->shift_type;
816: return(0);
817: }
819: /*@
820: EPSPowerGetShiftType - Gets the type of shifts used during the power
821: iteration.
823: Not Collective
825: Input Parameter:
826: . eps - the eigenproblem solver context
828: Input Parameter:
829: . shift - the type of shift
831: Level: advanced
833: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
834: @*/
835: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
836: {
842: PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
843: return(0);
844: }
846: static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
847: {
848: EPS_POWER *power = (EPS_POWER*)eps->data;
851: if (power->nonlinear != nonlinear) {
852: power->nonlinear = nonlinear;
853: eps->useds = PetscNot(nonlinear);
854: eps->state = EPS_STATE_INITIAL;
855: }
856: return(0);
857: }
859: /*@
860: EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
862: Logically Collective on eps
864: Input Parameters:
865: + eps - the eigenproblem solver context
866: - nonlinear - whether the problem is nonlinear or not
868: Options Database Key:
869: . -eps_power_nonlinear - Sets the nonlinear flag
871: Notes:
872: If this flag is set, the solver assumes that the problem is nonlinear,
873: that is, the operators that define the eigenproblem are not constant
874: matrices, but depend on the eigenvector: A(x)*x=lambda*B(x)*x. This is
875: different from the case of nonlinearity with respect to the eigenvalue
876: (use the NEP solver class for this kind of problems).
878: The way in which nonlinear operators are specified is very similar to
879: the case of PETSc's SNES solver. The difference is that the callback
880: functions are provided via composed functions "formFunction" and
881: "formJacobian" in each of the matrix objects passed as arguments of
882: EPSSetOperators(). The application context required for these functions
883: can be attached via a composed PetscContainer.
885: Level: advanced
887: .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
888: @*/
889: PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
890: {
896: PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
897: return(0);
898: }
900: static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
901: {
902: EPS_POWER *power = (EPS_POWER*)eps->data;
905: *nonlinear = power->nonlinear;
906: return(0);
907: }
909: /*@
910: EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
912: Not Collective
914: Input Parameter:
915: . eps - the eigenproblem solver context
917: Input Parameter:
918: . nonlinear - the nonlinear flag
920: Level: advanced
922: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
923: @*/
924: PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
925: {
931: PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
932: return(0);
933: }
935: static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
936: {
937: EPS_POWER *power = (EPS_POWER*)eps->data;
940: if (!power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems");
941: power->update = update;
942: eps->state = EPS_STATE_INITIAL;
943: return(0);
944: }
946: /*@
947: EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
948: for nonlinear problems. This potentially has a better convergence.
950: Logically Collective on eps
952: Input Parameters:
953: + eps - the eigenproblem solver context
954: - update - whether the residual is updated monolithically or not
956: Options Database Key:
957: . -eps_power_update - Sets the update flag
959: Level: advanced
961: .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
962: @*/
963: PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
964: {
970: PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
971: return(0);
972: }
974: static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
975: {
976: EPS_POWER *power = (EPS_POWER*)eps->data;
979: *update = power->update;
980: return(0);
981: }
983: /*@
984: EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
985: for nonlinear problems.
987: Not Collective
989: Input Parameter:
990: . eps - the eigenproblem solver context
992: Input Parameter:
993: . update - the update flag
995: Level: advanced
997: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
998: @*/
999: PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
1000: {
1006: PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
1007: return(0);
1008: }
1010: static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
1011: {
1013: EPS_POWER *power = (EPS_POWER*)eps->data;
1016: PetscObjectReference((PetscObject)snes);
1017: SNESDestroy(&power->snes);
1018: power->snes = snes;
1019: PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
1020: eps->state = EPS_STATE_INITIAL;
1021: return(0);
1022: }
1024: /*@
1025: EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
1026: eigenvalue solver (to be used in nonlinear inverse iteration).
1028: Collective on eps
1030: Input Parameters:
1031: + eps - the eigenvalue solver
1032: - snes - the nonlinear solver object
1034: Level: advanced
1036: .seealso: EPSPowerGetSNES()
1037: @*/
1038: PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
1039: {
1046: PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
1047: return(0);
1048: }
1050: static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
1051: {
1053: EPS_POWER *power = (EPS_POWER*)eps->data;
1056: if (!power->snes) {
1057: SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes);
1058: PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1);
1059: SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix);
1060: SNESAppendOptionsPrefix(power->snes,"eps_power_");
1061: PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
1062: PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options);
1063: }
1064: *snes = power->snes;
1065: return(0);
1066: }
1068: /*@
1069: EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
1070: with the eigenvalue solver.
1072: Not Collective
1074: Input Parameter:
1075: . eps - the eigenvalue solver
1077: Output Parameter:
1078: . snes - the nonlinear solver object
1080: Level: advanced
1082: .seealso: EPSPowerSetSNES()
1083: @*/
1084: PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
1085: {
1091: PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
1092: return(0);
1093: }
1095: PetscErrorCode EPSReset_Power(EPS eps)
1096: {
1098: EPS_POWER *power = (EPS_POWER*)eps->data;
1101: if (power->snes) { SNESReset(power->snes); }
1102: return(0);
1103: }
1105: PetscErrorCode EPSDestroy_Power(EPS eps)
1106: {
1108: EPS_POWER *power = (EPS_POWER*)eps->data;
1111: if (power->nonlinear) {
1112: SNESDestroy(&power->snes);
1113: }
1114: PetscFree(eps->data);
1115: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
1116: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
1117: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL);
1118: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL);
1119: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL);
1120: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL);
1121: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL);
1122: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL);
1123: return(0);
1124: }
1126: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
1127: {
1129: EPS_POWER *power = (EPS_POWER*)eps->data;
1130: PetscBool isascii;
1133: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1134: if (isascii) {
1135: if (power->nonlinear) {
1136: PetscViewerASCIIPrintf(viewer," using nonlinear inverse iteration\n");
1137: if (power->update) {
1138: PetscViewerASCIIPrintf(viewer," updating the residual monolithically\n");
1139: }
1140: if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
1141: PetscViewerASCIIPushTab(viewer);
1142: SNESView(power->snes,viewer);
1143: PetscViewerASCIIPopTab(viewer);
1144: } else {
1145: PetscViewerASCIIPrintf(viewer," %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
1146: }
1147: }
1148: return(0);
1149: }
1151: PetscErrorCode EPSComputeVectors_Power(EPS eps)
1152: {
1154: EPS_POWER *power = (EPS_POWER*)eps->data;
1155: PetscReal norm;
1156: PetscInt i;
1159: if (eps->twosided) {
1160: EPSComputeVectors_Twosided(eps);
1161: /* normalize (no need to take care of 2x2 blocks */
1162: for (i=0;i<eps->nconv;i++) {
1163: BVNormColumn(eps->V,i,NORM_2,&norm);
1164: BVScaleColumn(eps->V,i,1.0/norm);
1165: BVNormColumn(eps->W,i,NORM_2,&norm);
1166: BVScaleColumn(eps->W,i,1.0/norm);
1167: }
1168: } else if (!power->nonlinear) {
1169: EPSComputeVectors_Schur(eps);
1170: }
1171: return(0);
1172: }
1174: PetscErrorCode EPSSetDefaultST_Power(EPS eps)
1175: {
1177: EPS_POWER *power = (EPS_POWER*)eps->data;
1178: KSP ksp;
1179: PC pc;
1182: if (power->nonlinear) {
1183: eps->categ=EPS_CATEGORY_PRECOND;
1184: STGetKSP(eps->st,&ksp);
1185: /* Set ST as STPRECOND so it can carry one preconditioning matrix
1186: * It is useful when A and B are shell matrices
1187: */
1188: STSetType(eps->st,STPRECOND);
1189: KSPGetPC(ksp,&pc);
1190: PCSetType(pc,PCNONE);
1191: }
1192: return(0);
1193: }
1195: SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
1196: {
1197: EPS_POWER *ctx;
1201: PetscNewLog(eps,&ctx);
1202: eps->data = (void*)ctx;
1204: eps->useds = PETSC_TRUE;
1205: eps->hasts = PETSC_TRUE;
1206: eps->categ = EPS_CATEGORY_OTHER;
1208: eps->ops->setup = EPSSetUp_Power;
1209: eps->ops->setfromoptions = EPSSetFromOptions_Power;
1210: eps->ops->reset = EPSReset_Power;
1211: eps->ops->destroy = EPSDestroy_Power;
1212: eps->ops->view = EPSView_Power;
1213: eps->ops->backtransform = EPSBackTransform_Power;
1214: eps->ops->computevectors = EPSComputeVectors_Power;
1215: eps->ops->setdefaultst = EPSSetDefaultST_Power;
1216: eps->stopping = EPSStopping_Power;
1218: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
1219: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
1220: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power);
1221: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power);
1222: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power);
1223: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power);
1224: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power);
1225: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power);
1226: return(0);
1227: }