Actual source code: power.c
1: /*
3: SLEPc eigensolver: "power"
5: Method: Power Iteration
7: Description:
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: Algorithm:
18: The implemented algorithm is the simple power iteration working with
19: OP, the operator provided by the ST object. Converged eigenpairs are
20: deflated by restriction, so that several eigenpairs can be sought.
21: Symmetry is preserved in symmetric definite pencils. See the SLEPc
22: users guide for details.
24: Variable shifts can be used. There are two possible strategies for
25: updating shift: Rayleigh quotients and Wilkinson shifts.
27: References:
29: [1] B.N. Parlett, "The Symmetric Eigenvalue Problem", SIAM Classics in
30: Applied Mathematics (1998), pp 61-80 and 159-165.
32: Last update: June 2004
34: */
35: #include src/eps/epsimpl.h
36: #include slepcblaslapack.h
38: typedef struct {
39: EPSPowerShiftType shift_type;
40: } EPS_POWER;
44: PetscErrorCode EPSSetUp_POWER(EPS eps)
45: {
46: PetscErrorCode ierr;
47: EPS_POWER *power = (EPS_POWER *)eps->data;
48: int N;
49: PetscTruth flg;
50: STMatMode mode;
53: VecGetSize(eps->vec_initial,&N);
54: if (eps->ncv) {
55: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
56: }
57: else eps->ncv = eps->nev;
58: if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
59: if (!eps->tol) eps->tol = 1.e-7;
60: if (eps->which!=EPS_LARGEST_MAGNITUDE)
61: SETERRQ(1,"Wrong value of eps->which");
62: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
63: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&flg);
64: if (flg)
65: SETERRQ(PETSC_ERR_SUP,"Shift spectral transformation does not work with variable shifts");
66: STGetMatMode(eps->OP,&mode);
67: if (mode == STMATMODE_INPLACE)
68: SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
69: }
70: EPSAllocateSolution(eps);
71: EPSDefaultGetWork(eps,2);
72: return(0);
73: }
77: /*
78: EPSPowerUpdateShift - Computes the new shift to be used in the next
79: iteration of the power method. This function is invoked only when using
80: the option of variable shifts (see EPSPowerSetShiftType).
81: */
82: static PetscErrorCode EPSPowerUpdateShift(EPS eps,Vec v,PetscScalar* shift)
83: {
84: PetscErrorCode ierr;
85: EPS_POWER *power = (EPS_POWER *)eps->data;
86: Vec e, w;
87: Mat A;
88: PetscReal norm, rt1, rt2, cs1;
89: PetscScalar alpha, alpha1, alpha2, beta1, sn1;
92: e = eps->work[0];
93: w = eps->work[1];
94: STGetOperators(eps->OP,&A,PETSC_NULL);
96: /* compute the Rayleigh quotient R(v) assuming v is B-normalized */
97: MatMult(A,v,e);
98: VecDot(v,e,&alpha1);
100: /* in the case of Wilkinson the shift is improved */
101: if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
102: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
103: SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
104: #endif
105: /* beta1 is the norm of the residual associated to R(v) */
106: alpha = -alpha1;
107: VecAXPY(&alpha,v,e);
108: STNorm(eps->OP,e,&norm);
109: beta1 = norm;
110:
111: /* alfa2 = (e'*A*e)/(beta1*beta1), where e is the residual */
112: MatMult(A,e,w);
113: VecDot(e,w,&alpha2);
114: alpha2 = alpha2 / (beta1 * beta1);
116: /* choose the eigenvalue of [alfa1 beta1; beta1 alfa2] closest to alpha1 */
117: LAlaev2_(&alpha1,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
118: if (PetscAbsScalar(rt1-alpha1) < PetscAbsScalar(rt2-alpha1)) {
119: *shift = rt1;
120: } else {
121: *shift = rt2;
122: }
123: }
124: else *shift = alpha1;
126: return(0);
127: }
131: PetscErrorCode EPSSolve_POWER(EPS eps)
132: {
133: PetscErrorCode ierr;
134: EPS_POWER *power = (EPS_POWER *)eps->data;
135: int i;
136: Vec v, y, e;
137: PetscReal relerr, norm;
138: PetscScalar theta, alpha, rho;
141: v = eps->V[0];
142: y = eps->AV[0];
143: e = eps->work[0];
145: VecCopy(eps->vec_initial,y);
147: eps->nconv = 0;
148: eps->its = 0;
150: for (i=0;i<eps->ncv;i++) eps->eigi[i]=0.0;
152: while (eps->its<eps->max_it) {
154: eps->its = eps->its + 1;
156: /* deflation of converged eigenvectors */
157: EPSPurge(eps,y);
159: /* v = y/||y||_B */
160: VecCopy(y,v);
161: STNorm(eps->OP,y,&norm);
162: alpha = 1.0/norm;
163: VecScale(&alpha,v);
165: /* y = OP v */
166: STApply(eps->OP,v,y);
168: /* theta = (y,v)_B */
169: STInnerProduct(eps->OP,y,v,&theta);
171: /* compute residual norm */
172: VecCopy(y,e);
173: alpha = -theta;
174: VecAXPY(&alpha,v,e);
175: VecNorm(e,NORM_2,&relerr);
176: relerr = relerr / PetscAbsScalar(theta);
177: eps->errest[eps->nconv] = relerr;
179: /* update eigenvalue and shift */
180: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
181: EPSPowerUpdateShift(eps,v,&rho);
182: /* change the shift only if rho is not too close to an eigenvalue */
183: if (relerr > 1000*eps->tol) {
184: STSetShift(eps->OP,rho);
185: }
186: eps->eigr[eps->nconv] = rho;
187: } else {
188: eps->eigr[eps->nconv] = theta;
189: }
191: /* if ||y-theta v||_2 / |theta| < tol, accept eigenpair */
192: if (relerr<eps->tol) {
193: eps->nconv = eps->nconv + 1;
194: if (eps->nconv==eps->nev) break;
195: v = eps->V[eps->nconv];
196: }
198: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
199: }
201: if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
202: else eps->reason = EPS_DIVERGED_ITS;
204: return(0);
205: }
209: PetscErrorCode EPSBackTransform_POWER(EPS eps)
210: {
211: PetscErrorCode ierr;
212: EPS_POWER *power = (EPS_POWER *)eps->data;
215: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) {
216: EPSBackTransform_Default(eps);
217: }
218: return(0);
219: }
223: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
224: {
225: PetscErrorCode ierr;
226: EPS_POWER *power = (EPS_POWER *)eps->data;
227: PetscTruth flg;
228: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
231: PetscOptionsHead("POWER options");
232: PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],(int*)&power->shift_type,&flg);
233: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
234: STSetType(eps->OP,STSINV);
235: }
236: PetscOptionsTail();
237: return(0);
238: }
240: EXTERN_C_BEGIN
243: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
244: {
245: EPS_POWER *power = (EPS_POWER *)eps->data;
248: switch (shift) {
249: case EPSPOWER_SHIFT_CONSTANT:
250: case EPSPOWER_SHIFT_RAYLEIGH:
251: case EPSPOWER_SHIFT_WILKINSON:
252: power->shift_type = shift;
253: break;
254: default:
255: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
256: }
257: return(0);
258: }
259: EXTERN_C_END
263: /*@
264: EPSPowerSetShiftType - Sets the type of shifts used during the power
265: iteration. This can be used to emulate the Rayleigh Quotient Iteration
266: (RQI) method.
268: Collective on EPS
270: Input Parameters:
271: + eps - the eigenproblem solver context
272: - shift - the type of shift
274: Options Database Key:
275: . -eps_power_shift_type - Sets the shift type (either 'constant' or
276: 'rayleigh' or 'wilkinson')
278: Notes:
279: By default, shifts are constant (EPSPOWER_SHIFT_CONSTANT) and the iteration
280: is the simple power method (or inverse iteration if a shift-and-invert
281: transformation is being used).
283: A variable shift can be specified (EPSPOWER_SHIFT_RAYLEIGH or
284: EPSPOWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
285: a cubic converging method as RQI. See the users manual for details.
286:
287: Level: advanced
289: .seealso: EPSGetShiftType(), STSetShift()
290: @*/
291: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
292: {
293: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);
297: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
298: if (f) {
299: (*f)(eps,shift);
300: }
301: return(0);
302: }
304: EXTERN_C_BEGIN
307: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
308: {
309: EPS_POWER *power = (EPS_POWER *)eps->data;
311: *shift = power->shift_type;
312: return(0);
313: }
314: EXTERN_C_END
318: /*@
319: EPSPowerGetShiftType - Gets the type of shifts used during the power
320: iteration.
322: Collective on EPS
324: Input Parameter:
325: . eps - the eigenproblem solver context
327: Input Parameter:
328: . shift - the type of shift
330: Level: advanced
332: .seealso: EPSSetShiftType()
333: @*/
334: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
335: {
336: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);
340: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
341: if (f) {
342: (*f)(eps,shift);
343: }
344: return(0);
345: }
349: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
350: {
351: PetscErrorCode ierr;
352: EPS_POWER *power = (EPS_POWER *)eps->data;
353: PetscTruth isascii;
354: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
357: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
358: if (!isascii) {
359: SETERRQ1(1,"Viewer type %s not supported for EPSPOWER",((PetscObject)viewer)->type_name);
360: }
361: PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
362: return(0);
363: }
365: EXTERN_C_BEGIN
368: PetscErrorCode EPSCreate_POWER(EPS eps)
369: {
370: PetscErrorCode ierr;
371: EPS_POWER *power;
374: PetscNew(EPS_POWER,&power);
375: PetscMemzero(power,sizeof(EPS_POWER));
376: PetscLogObjectMemory(eps,sizeof(EPS_POWER));
377: eps->data = (void *) power;
378: eps->ops->solve = EPSSolve_POWER;
379: eps->ops->setup = EPSSetUp_POWER;
380: eps->ops->setfromoptions = EPSSetFromOptions_POWER;
381: eps->ops->destroy = EPSDestroy_Default;
382: eps->ops->view = EPSView_POWER;
383: eps->ops->backtransform = EPSBackTransform_POWER;
384: eps->ops->computevectors = EPSComputeVectors_Default;
385: power->shift_type = EPSPOWER_SHIFT_CONSTANT;
386: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
387: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
388: return(0);
389: }
390: EXTERN_C_END