Actual source code: power.c
1: /*
3: SLEPc eigensolver: "power"
5: Method: Power Iteration
7: Algorithm:
9: This solver implements the power iteration for finding dominant
10: eigenpairs. It also includes the following well-known methods:
11: - Inverse Iteration: when used in combination with shift-and-invert
12: spectral transformation.
13: - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
14: a variable shift.
16: References:
18: [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report STR-2,
19: available at http://www.grycap.upv.es/slepc.
21: Last update: Feb 2009
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: SLEPc - Scalable Library for Eigenvalue Problem Computations
25: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
27: This file is part of SLEPc.
28:
29: SLEPc is free software: you can redistribute it and/or modify it under the
30: terms of version 3 of the GNU Lesser General Public License as published by
31: the Free Software Foundation.
33: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
34: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
35: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
36: more details.
38: You should have received a copy of the GNU Lesser General Public License
39: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: */
43: #include private/epsimpl.h
44: #include slepcblaslapack.h
46: typedef struct {
47: EPSPowerShiftType shift_type;
48: } EPS_POWER;
52: PetscErrorCode EPSSetUp_POWER(EPS eps)
53: {
55: EPS_POWER *power = (EPS_POWER *)eps->data;
56: PetscInt N;
57: PetscTruth flg;
58: STMatMode mode;
61: VecGetSize(eps->vec_initial,&N);
62: if (eps->ncv) {
63: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
64: }
65: else eps->ncv = eps->nev;
66: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
67: if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
68: if (eps->which!=EPS_LARGEST_MAGNITUDE)
69: SETERRQ(1,"Wrong value of eps->which");
70: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
71: PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
72: if (!flg)
73: SETERRQ(PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert ST");
74: STGetMatMode(eps->OP,&mode);
75: if (mode == STMATMODE_INPLACE)
76: SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
77: }
78: if (eps->extraction) {
79: PetscInfo(eps,"Warning: extraction type ignored\n");
80: }
81: EPSAllocateSolution(eps);
82: EPSDefaultGetWork(eps,3);
83: return(0);
84: }
88: PetscErrorCode EPSSolve_POWER(EPS eps)
89: {
91: EPS_POWER *power = (EPS_POWER *)eps->data;
92: PetscInt i, nsv;
93: Vec v, y, e, *SV;
94: Mat A;
95: PetscReal relerr, norm, rt1, rt2, cs1, anorm;
96: PetscScalar theta, rho, delta, sigma, alpha2, beta1, sn1;
97: PetscTruth breakdown;
100: v = eps->V[0];
101: y = eps->work[2];
102: e = eps->work[0];
104: /* prepare for selective orthogonalization of converged vectors */
105: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
106: PetscMalloc(eps->nev*sizeof(Vec),&SV);
107: for (i=0;i<eps->nds;i++) SV[i]=eps->DS[i];
108: if (eps->nev>1) {
109: STGetOperators(eps->OP,&A,PETSC_NULL);
110: MatNorm(A,NORM_INFINITY,&anorm);
111: }
112: }
114: EPSGetStartVector(eps,0,v,PETSC_NULL);
115: STGetShift(eps->OP,&sigma); /* original shift */
116: rho = sigma;
118: while (eps->reason == EPS_CONVERGED_ITERATING) {
119: eps->its = eps->its + 1;
121: /* y = OP v */
122: STApply(eps->OP,v,y);
124: /* theta = (v,y)_B */
125: IPInnerProduct(eps->ip,v,y,&theta);
127: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
129: /* approximate eigenvalue is the Rayleigh quotient */
130: eps->eigr[eps->nconv] = theta;
132: /* compute relative error as ||y-theta v||_2/|theta| */
133: VecCopy(y,e);
134: VecAXPY(e,-theta,v);
135: VecNorm(e,NORM_2,&norm);
136: relerr = norm / PetscAbsScalar(theta);
138: } else { /* RQI */
140: /* delta = ||y||_B */
141: IPNorm(eps->ip,y,&norm);
142: delta = norm;
144: /* compute relative error */
145: if (rho == 0.0) relerr = PETSC_MAX;
146: else relerr = 1.0 / (norm*PetscAbsScalar(rho));
148: /* approximate eigenvalue is the shift */
149: eps->eigr[eps->nconv] = rho;
151: /* compute new shift */
152: if (relerr<eps->tol) {
153: rho = sigma; /* if converged, restore original shift */
154: STSetShift(eps->OP,rho);
155: } else {
156: rho = rho + theta/(delta*delta); /* Rayleigh quotient R(v) */
157: if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
158: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
159: SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
160: #else
161: /* beta1 is the norm of the residual associated to R(v) */
162: VecAXPY(v,-theta/(delta*delta),y);
163: VecScale(v,1.0/delta);
164: IPNorm(eps->ip,v,&norm);
165: beta1 = norm;
166:
167: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
168: STGetOperators(eps->OP,&A,PETSC_NULL);
169: MatMult(A,v,e);
170: VecDot(v,e,&alpha2);
171: alpha2 = alpha2 / (beta1 * beta1);
173: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
174: LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
175: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
176: else rho = rt2;
177: #endif
178: }
179: /* update operator according to new shift */
180: PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
181: STSetShift(eps->OP,rho);
182: PetscPopErrorHandler();
183: if (ierr) {
184: eps->eigr[eps->nconv] = rho;
185: relerr = PETSC_MACHINE_EPSILON;
186: rho = sigma;
187: STSetShift(eps->OP,rho);
188: }
189: }
190: }
192: eps->errest[eps->nconv] = relerr;
193: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
195: /* purge previously converged eigenvectors */
196: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
197: nsv = eps->nds;
198: for (i=0;i<eps->nconv;i++) {
199: if(PetscAbsScalar(rho-eps->eigr[i])>eps->its*anorm/1000) SV[nsv++]=eps->V[i];
200: }
201: IPOrthogonalize(eps->ip,nsv,PETSC_NULL,SV,y,PETSC_NULL,&norm,PETSC_NULL,eps->work[1],PETSC_NULL);
202: } else {
203: IPOrthogonalize(eps->ip,eps->nds+eps->nconv,PETSC_NULL,eps->DSV,y,PETSC_NULL,&norm,PETSC_NULL,eps->work[1],PETSC_NULL);
204: }
206: /* v = y/||y||_B */
207: VecCopy(y,v);
208: VecScale(v,1.0/norm);
210: /* if relerr<tol, accept eigenpair */
211: if (relerr<eps->tol) {
212: eps->nconv = eps->nconv + 1;
213: if (eps->nconv==eps->nev) eps->reason = EPS_CONVERGED_TOL;
214: else {
215: v = eps->V[eps->nconv];
216: EPSGetStartVector(eps,eps->nconv,v,&breakdown);
217: if (breakdown) {
218: eps->reason = EPS_DIVERGED_BREAKDOWN;
219: PetscInfo(eps,"Unable to generate more start vectors\n");
220: }
221: }
222: }
224: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
225: }
227: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
228: PetscFree(SV);
229: }
231: return(0);
232: }
236: PetscErrorCode EPSSolve_TS_POWER(EPS eps)
237: {
239: EPS_POWER *power = (EPS_POWER *)eps->data;
240: Vec v, w, y, z, e;
241: Mat A;
242: PetscReal relerr, norm, rt1, rt2, cs1;
243: PetscScalar theta, alpha, beta, rho, delta, sigma, alpha2, beta1, sn1;
246: v = eps->V[0];
247: y = eps->work[1];
248: e = eps->work[0];
249: w = eps->W[0];
250: z = eps->work[2];
252: EPSGetStartVector(eps,0,v,PETSC_NULL);
253: EPSGetLeftStartVector(eps,0,w);
254: STGetShift(eps->OP,&sigma); /* original shift */
255: rho = sigma;
257: while (eps->its<eps->max_it) {
258: eps->its++;
259:
260: /* y = OP v, z = OP' w */
261: STApply(eps->OP,v,y);
262: STApplyTranspose(eps->OP,w,z);
264: /* theta = (v,z)_B */
265: IPInnerProduct(eps->ip,v,z,&theta);
267: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
269: /* approximate eigenvalue is the Rayleigh quotient */
270: eps->eigr[eps->nconv] = theta;
272: /* compute relative errors (right and left) */
273: VecCopy(y,e);
274: VecAXPY(e,-theta,v);
275: VecNorm(e,NORM_2,&norm);
276: relerr = norm / PetscAbsScalar(theta);
277: eps->errest[eps->nconv] = relerr;
278: VecCopy(z,e);
279: VecAXPY(e,-theta,w);
280: VecNorm(e,NORM_2,&norm);
281: relerr = norm / PetscAbsScalar(theta);
282: eps->errest_left[eps->nconv] = relerr;
284: } else { /* RQI */
286: /* delta = sqrt(y,z)_B */
287: IPInnerProduct(eps->ip,y,z,&alpha);
288: if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
289: delta = PetscSqrtScalar(alpha);
291: /* compute relative error */
292: if (rho == 0.0) relerr = PETSC_MAX;
293: else relerr = 1.0 / (PetscAbsScalar(delta*rho));
294: eps->errest[eps->nconv] = relerr;
295: eps->errest_left[eps->nconv] = relerr;
297: /* approximate eigenvalue is the shift */
298: eps->eigr[eps->nconv] = rho;
300: /* compute new shift */
301: if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
302: rho = sigma; /* if converged, restore original shift */
303: STSetShift(eps->OP,rho);
304: } else {
305: rho = rho + theta/(delta*delta); /* Rayleigh quotient R(v,w) */
306: if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
307: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
308: SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
309: #else
310: /* beta1 is the norm of the residual associated to R(v,w) */
311: VecAXPY(v,-theta/(delta*delta),y);
312: VecScale(v,1.0/delta);
313: IPNorm(eps->ip,v,&norm);
314: beta1 = norm;
315:
316: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
317: STGetOperators(eps->OP,&A,PETSC_NULL);
318: MatMult(A,v,e);
319: VecDot(v,e,&alpha2);
320: alpha2 = alpha2 / (beta1 * beta1);
322: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
323: LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
324: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
325: else rho = rt2;
326: #endif
327: }
328: /* update operator according to new shift */
329: PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
330: STSetShift(eps->OP,rho);
331: PetscPopErrorHandler();
332: if (ierr) {
333: eps->eigr[eps->nconv] = rho;
334: eps->errest[eps->nconv] = PETSC_MACHINE_EPSILON;
335: eps->errest_left[eps->nconv] = PETSC_MACHINE_EPSILON;
336: rho = sigma;
337: STSetShift(eps->OP,rho);
338: }
339: }
340: }
342: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
343: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,eps->nconv+1);
345: /* purge previously converged eigenvectors */
346: IPBiOrthogonalize(eps->ip,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
347: IPBiOrthogonalize(eps->ip,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);
349: /* normalize so that (y,z)_B=1 */
350: VecCopy(y,v);
351: VecCopy(z,w);
352: IPInnerProduct(eps->ip,y,z,&alpha);
353: if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
354: delta = PetscSqrtScalar(PetscAbsScalar(alpha));
355: beta = 1.0/PetscConj(alpha/delta);
356: delta = 1.0/delta;
357: VecScale(w,beta);
358: VecScale(v,delta);
360: /* if relerr<tol (both right and left), accept eigenpair */
361: if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
362: eps->nconv = eps->nconv + 1;
363: if (eps->nconv==eps->nev) break;
364: v = eps->V[eps->nconv];
365: EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
366: w = eps->W[eps->nconv];
367: EPSGetLeftStartVector(eps,eps->nconv,w);
368: }
369: }
371: if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
372: else eps->reason = EPS_DIVERGED_ITS;
374: return(0);
375: }
379: PetscErrorCode EPSBackTransform_POWER(EPS eps)
380: {
382: EPS_POWER *power = (EPS_POWER *)eps->data;
385: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) {
386: EPSBackTransform_Default(eps);
387: }
388: return(0);
389: }
393: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
394: {
396: EPS_POWER *power = (EPS_POWER *)eps->data;
397: PetscTruth flg;
398: PetscInt i;
399: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
402: PetscOptionsHead("POWER options");
403: PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],&i,&flg);
404: if (flg ) power->shift_type = (EPSPowerShiftType)i;
405: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
406: STSetType(eps->OP,STSINV);
407: }
408: PetscOptionsTail();
409: return(0);
410: }
415: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
416: {
417: EPS_POWER *power = (EPS_POWER *)eps->data;
420: switch (shift) {
421: case EPSPOWER_SHIFT_CONSTANT:
422: case EPSPOWER_SHIFT_RAYLEIGH:
423: case EPSPOWER_SHIFT_WILKINSON:
424: power->shift_type = shift;
425: break;
426: default:
427: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
428: }
429: return(0);
430: }
435: /*@
436: EPSPowerSetShiftType - Sets the type of shifts used during the power
437: iteration. This can be used to emulate the Rayleigh Quotient Iteration
438: (RQI) method.
440: Collective on EPS
442: Input Parameters:
443: + eps - the eigenproblem solver context
444: - shift - the type of shift
446: Options Database Key:
447: . -eps_power_shift_type - Sets the shift type (either 'constant' or
448: 'rayleigh' or 'wilkinson')
450: Notes:
451: By default, shifts are constant (EPSPOWER_SHIFT_CONSTANT) and the iteration
452: is the simple power method (or inverse iteration if a shift-and-invert
453: transformation is being used).
455: A variable shift can be specified (EPSPOWER_SHIFT_RAYLEIGH or
456: EPSPOWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
457: a cubic converging method as RQI. See the users manual for details.
458:
459: Level: advanced
461: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
462: @*/
463: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
464: {
465: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);
469: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
470: if (f) {
471: (*f)(eps,shift);
472: }
473: return(0);
474: }
479: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
480: {
481: EPS_POWER *power = (EPS_POWER *)eps->data;
483: *shift = power->shift_type;
484: return(0);
485: }
490: /*@C
491: EPSPowerGetShiftType - Gets the type of shifts used during the power
492: iteration.
494: Collective on EPS
496: Input Parameter:
497: . eps - the eigenproblem solver context
499: Input Parameter:
500: . shift - the type of shift
502: Level: advanced
504: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
505: @*/
506: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
507: {
508: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);
512: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
513: if (f) {
514: (*f)(eps,shift);
515: }
516: return(0);
517: }
521: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
522: {
524: EPS_POWER *power = (EPS_POWER *)eps->data;
525: PetscTruth isascii;
526: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
529: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
530: if (!isascii) {
531: SETERRQ1(1,"Viewer type %s not supported for EPSPOWER",((PetscObject)viewer)->type_name);
532: }
533: PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
534: return(0);
535: }
540: PetscErrorCode EPSCreate_POWER(EPS eps)
541: {
543: EPS_POWER *power;
546: PetscNew(EPS_POWER,&power);
547: PetscLogObjectMemory(eps,sizeof(EPS_POWER));
548: eps->data = (void *) power;
549: eps->ops->solve = EPSSolve_POWER;
550: eps->ops->solvets = EPSSolve_TS_POWER;
551: eps->ops->setup = EPSSetUp_POWER;
552: eps->ops->setfromoptions = EPSSetFromOptions_POWER;
553: eps->ops->destroy = EPSDestroy_Default;
554: eps->ops->view = EPSView_POWER;
555: eps->ops->backtransform = EPSBackTransform_POWER;
556: eps->ops->computevectors = EPSComputeVectors_Default;
557: power->shift_type = EPSPOWER_SHIFT_CONSTANT;
558: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
559: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
560: return(0);
561: }