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: June 2005
23: */
24: #include src/eps/epsimpl.h
25: #include slepcblaslapack.h
27: typedef struct {
28: EPSPowerShiftType shift_type;
29: } EPS_POWER;
33: PetscErrorCode EPSSetUp_POWER(EPS eps)
34: {
36: EPS_POWER *power = (EPS_POWER *)eps->data;
37: PetscInt N;
38: PetscTruth flg;
39: STMatMode mode;
42: VecGetSize(eps->vec_initial,&N);
43: if (eps->ncv) {
44: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
45: }
46: else eps->ncv = eps->nev;
47: if (!eps->max_it) eps->max_it = PetscMax(2000,100*N);
48: if (!eps->tol) eps->tol = 1.e-7;
49: if (eps->which!=EPS_LARGEST_MAGNITUDE)
50: SETERRQ(1,"Wrong value of eps->which");
51: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
52: PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
53: if (!flg)
54: SETERRQ(PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert ST");
55: STGetMatMode(eps->OP,&mode);
56: if (mode == STMATMODE_INPLACE)
57: SETERRQ(PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
58: }
59: EPSAllocateSolution(eps);
60: EPSDefaultGetWork(eps,1);
61: return(0);
62: }
66: PetscErrorCode EPSSolve_POWER(EPS eps)
67: {
69: EPS_POWER *power = (EPS_POWER *)eps->data;
70: int i, nsv;
71: Vec v, y, e, *SV;
72: Mat A;
73: PetscReal relerr, norm, rt1, rt2, cs1, anorm;
74: PetscScalar theta, rho, delta, sigma, alpha2, beta1, sn1;
75: PetscTruth breakdown;
78: v = eps->V[0];
79: y = eps->AV[0];
80: e = eps->work[0];
82: /* prepare for selective orthogonalization of converged vectors */
83: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
84: PetscMalloc(eps->nev*sizeof(Vec),&SV);
85: for (i=0;i<eps->nds;i++) SV[i]=eps->DS[i];
86: if (eps->nev>1) {
87: STGetOperators(eps->OP,&A,PETSC_NULL);
88: MatNorm(A,NORM_INFINITY,&anorm);
89: }
90: }
92: EPSGetStartVector(eps,0,v,PETSC_NULL);
93: STGetShift(eps->OP,&sigma); /* original shift */
94: rho = sigma;
96: eps->nconv = 0;
97: eps->its = 0;
99: for (i=0;i<eps->ncv;i++) eps->eigi[i]=0.0;
101: while (eps->reason == EPS_CONVERGED_ITERATING) {
103: /* y = OP v */
104: STApply(eps->OP,v,y);
106: /* theta = (v,y)_B */
107: STInnerProduct(eps->OP,v,y,&theta);
109: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
111: /* approximate eigenvalue is the Rayleigh quotient */
112: eps->eigr[eps->nconv] = theta;
114: /* compute relative error as ||y-theta v||_2/|theta| */
115: VecCopy(y,e);
116: VecAXPY(e,-theta,v);
117: VecNorm(e,NORM_2,&norm);
118: relerr = norm / PetscAbsScalar(theta);
120: } else { /* RQI */
122: /* delta = ||y||_B */
123: STNorm(eps->OP,y,&norm);
124: delta = norm;
126: /* compute relative error */
127: if (rho == 0.0) relerr = PETSC_MAX;
128: else relerr = 1.0 / (norm*PetscAbsScalar(rho));
130: /* approximate eigenvalue is the shift */
131: eps->eigr[eps->nconv] = rho;
133: /* compute new shift */
134: if (relerr<eps->tol) {
135: rho = sigma; /* if converged, restore original shift */
136: STSetShift(eps->OP,rho);
137: } else {
138: rho = rho + theta/(delta*delta); /* Rayleigh quotient R(v) */
139: if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
140: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
141: SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
142: #else
143: /* beta1 is the norm of the residual associated to R(v) */
144: VecAXPY(v,-theta/(delta*delta),y);
145: VecScale(v,1.0/delta);
146: STNorm(eps->OP,v,&norm);
147: beta1 = norm;
148:
149: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
150: STGetOperators(eps->OP,&A,PETSC_NULL);
151: MatMult(A,v,e);
152: VecDot(v,e,&alpha2);
153: alpha2 = alpha2 / (beta1 * beta1);
155: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
156: LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
157: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
158: else rho = rt2;
159: #endif
160: }
161: /* update operator according to new shift */
162: PetscPushErrorHandler(SlepcQuietErrorHandler,PETSC_NULL);
163: STSetShift(eps->OP,rho);
164: PetscPopErrorHandler();
165: if (ierr) {
166: eps->eigr[eps->nconv] = rho;
167: relerr = PETSC_MACHINE_EPSILON;
168: rho = sigma;
169: STSetShift(eps->OP,rho);
170: }
171: }
172: }
174: eps->errest[eps->nconv] = relerr;
175: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
177: /* purge previously converged eigenvectors */
178: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
179: nsv = eps->nds;
180: for (i=0;i<eps->nconv;i++) {
181: if(PetscAbsScalar(rho-eps->eigr[i])>(eps->its+1)*anorm/1000) SV[nsv++]=eps->V[i];
182: }
183: EPSOrthogonalize(eps,nsv,SV,y,PETSC_NULL,&norm,PETSC_NULL);
184: } else {
185: EPSOrthogonalize(eps,eps->nds+eps->nconv,eps->DSV,y,PETSC_NULL,&norm,PETSC_NULL);
186: }
188: /* v = y/||y||_B */
189: VecCopy(y,v);
190: VecScale(v,1.0/norm);
192: /* if relerr<tol, accept eigenpair */
193: if (relerr<eps->tol) {
194: eps->nconv = eps->nconv + 1;
195: if (eps->nconv==eps->nev) break;
196: v = eps->V[eps->nconv];
197: EPSGetStartVector(eps,eps->nconv,v,&breakdown);
198: if (breakdown) {
199: eps->reason = EPS_DIVERGED_BREAKDOWN;
200: PetscInfo(eps,"Unable to generate more start vectors\n");
201: }
202: }
204: eps->its = eps->its + 1;
205: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
206: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
207: }
209: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
210: PetscFree(SV);
211: }
213: return(0);
214: }
218: PetscErrorCode EPSSolve_TS_POWER(EPS eps)
219: {
221: EPS_POWER *power = (EPS_POWER *)eps->data;
222: int i;
223: Vec v, w, y, z, e;
224: Mat A;
225: PetscReal relerr, norm, rt1, rt2, cs1;
226: PetscScalar theta, alpha, beta, rho, delta, sigma, alpha2, beta1, sn1;
229: v = eps->V[0];
230: y = eps->AV[0];
231: e = eps->work[0];
232: w = eps->W[0];
233: z = eps->AW[0];
235: EPSGetStartVector(eps,0,v,PETSC_NULL);
236: EPSGetLeftStartVector(eps,0,w);
237: STGetShift(eps->OP,&sigma); /* original shift */
238: rho = sigma;
240: eps->nconv = 0;
241: eps->its = 0;
243: for (i=0;i<eps->ncv;i++) eps->eigi[i]=0.0;
245: while (eps->its<eps->max_it) {
247: /* y = OP v, z = OP' w */
248: STApply(eps->OP,v,y);
249: STApplyTranspose(eps->OP,w,z);
251: /* theta = (v,z)_B */
252: STInnerProduct(eps->OP,v,z,&theta);
254: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
256: /* approximate eigenvalue is the Rayleigh quotient */
257: eps->eigr[eps->nconv] = theta;
259: /* compute relative errors (right and left) */
260: VecCopy(y,e);
261: VecAXPY(e,-theta,v);
262: VecNorm(e,NORM_2,&norm);
263: relerr = norm / PetscAbsScalar(theta);
264: eps->errest[eps->nconv] = relerr;
265: VecCopy(z,e);
266: VecAXPY(e,-theta,w);
267: VecNorm(e,NORM_2,&norm);
268: relerr = norm / PetscAbsScalar(theta);
269: eps->errest_left[eps->nconv] = relerr;
271: } else { /* RQI */
273: /* delta = sqrt(y,z)_B */
274: STInnerProduct(eps->OP,y,z,&alpha);
275: if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
276: delta = PetscSqrtScalar(alpha);
278: /* compute relative error */
279: if (rho == 0.0) relerr = PETSC_MAX;
280: else relerr = 1.0 / (PetscAbsScalar(delta*rho));
281: eps->errest[eps->nconv] = relerr;
282: eps->errest_left[eps->nconv] = relerr;
284: /* approximate eigenvalue is the shift */
285: eps->eigr[eps->nconv] = rho;
287: /* compute new shift */
288: if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
289: rho = sigma; /* if converged, restore original shift */
290: STSetShift(eps->OP,rho);
291: } else {
292: rho = rho + theta/(delta*delta); /* Rayleigh quotient R(v,w) */
293: if (power->shift_type == EPSPOWER_SHIFT_WILKINSON) {
294: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
295: SETERRQ(PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable.");
296: #else
297: /* beta1 is the norm of the residual associated to R(v,w) */
298: VecAXPY(v,-theta/(delta*delta),y);
299: VecScale(v,1.0/delta);
300: STNorm(eps->OP,v,&norm);
301: beta1 = norm;
302:
303: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
304: STGetOperators(eps->OP,&A,PETSC_NULL);
305: MatMult(A,v,e);
306: VecDot(v,e,&alpha2);
307: alpha2 = alpha2 / (beta1 * beta1);
309: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
310: LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1);
311: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
312: else rho = rt2;
313: #endif
314: }
315: /* update operator according to new shift */
316: PetscPushErrorHandler(SlepcQuietErrorHandler,PETSC_NULL);
317: STSetShift(eps->OP,rho);
318: PetscPopErrorHandler();
319: if (ierr) {
320: eps->eigr[eps->nconv] = rho;
321: eps->errest[eps->nconv] = PETSC_MACHINE_EPSILON;
322: eps->errest_left[eps->nconv] = PETSC_MACHINE_EPSILON;
323: rho = sigma;
324: STSetShift(eps->OP,rho);
325: }
326: }
327: }
329: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nconv+1);
330: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest_left,eps->nconv+1);
332: /* purge previously converged eigenvectors */
333: EPSBiOrthogonalize(eps,eps->nconv,eps->V,eps->W,z,PETSC_NULL,PETSC_NULL);
334: EPSBiOrthogonalize(eps,eps->nconv,eps->W,eps->V,y,PETSC_NULL,PETSC_NULL);
336: /* normalize so that (y,z)_B=1 */
337: VecCopy(y,v);
338: VecCopy(z,w);
339: STInnerProduct(eps->OP,y,z,&alpha);
340: if (alpha==0.0) SETERRQ(1,"Breakdown in two-sided Power/RQI");
341: delta = PetscSqrtScalar(PetscAbsScalar(alpha));
342: beta = 1.0/PetscConj(alpha/delta);
343: delta = 1.0/delta;
344: VecScale(w,beta);
345: VecScale(v,delta);
347: /* if relerr<tol (both right and left), accept eigenpair */
348: if (eps->errest[eps->nconv]<eps->tol && eps->errest_left[eps->nconv]<eps->tol) {
349: eps->nconv = eps->nconv + 1;
350: if (eps->nconv==eps->nev) break;
351: v = eps->V[eps->nconv];
352: EPSGetStartVector(eps,eps->nconv,v,PETSC_NULL);
353: w = eps->W[eps->nconv];
354: EPSGetLeftStartVector(eps,eps->nconv,w);
355: }
357: eps->its = eps->its + 1;
358: }
360: if( eps->nconv == eps->nev ) eps->reason = EPS_CONVERGED_TOL;
361: else eps->reason = EPS_DIVERGED_ITS;
363: return(0);
364: }
368: PetscErrorCode EPSBackTransform_POWER(EPS eps)
369: {
371: EPS_POWER *power = (EPS_POWER *)eps->data;
374: if (power->shift_type == EPSPOWER_SHIFT_CONSTANT) {
375: EPSBackTransform_Default(eps);
376: }
377: return(0);
378: }
382: PetscErrorCode EPSSetFromOptions_POWER(EPS eps)
383: {
385: EPS_POWER *power = (EPS_POWER *)eps->data;
386: PetscTruth flg;
387: PetscInt i;
388: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
391: PetscOptionsHead("POWER options");
392: PetscOptionsEList("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",shift_list,3,shift_list[power->shift_type],&i,&flg);
393: if (flg ) power->shift_type = (EPSPowerShiftType)i;
394: if (power->shift_type != EPSPOWER_SHIFT_CONSTANT) {
395: STSetType(eps->OP,STSINV);
396: }
397: PetscOptionsTail();
398: return(0);
399: }
404: PetscErrorCode EPSPowerSetShiftType_POWER(EPS eps,EPSPowerShiftType shift)
405: {
406: EPS_POWER *power = (EPS_POWER *)eps->data;
409: switch (shift) {
410: case EPSPOWER_SHIFT_CONSTANT:
411: case EPSPOWER_SHIFT_RAYLEIGH:
412: case EPSPOWER_SHIFT_WILKINSON:
413: power->shift_type = shift;
414: break;
415: default:
416: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
417: }
418: return(0);
419: }
424: /*@
425: EPSPowerSetShiftType - Sets the type of shifts used during the power
426: iteration. This can be used to emulate the Rayleigh Quotient Iteration
427: (RQI) method.
429: Collective on EPS
431: Input Parameters:
432: + eps - the eigenproblem solver context
433: - shift - the type of shift
435: Options Database Key:
436: . -eps_power_shift_type - Sets the shift type (either 'constant' or
437: 'rayleigh' or 'wilkinson')
439: Notes:
440: By default, shifts are constant (EPSPOWER_SHIFT_CONSTANT) and the iteration
441: is the simple power method (or inverse iteration if a shift-and-invert
442: transformation is being used).
444: A variable shift can be specified (EPSPOWER_SHIFT_RAYLEIGH or
445: EPSPOWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
446: a cubic converging method as RQI. See the users manual for details.
447:
448: Level: advanced
450: .seealso: EPSGetShiftType(), STSetShift()
451: @*/
452: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
453: {
454: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType);
458: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerSetShiftType_C",(void (**)())&f);
459: if (f) {
460: (*f)(eps,shift);
461: }
462: return(0);
463: }
468: PetscErrorCode EPSPowerGetShiftType_POWER(EPS eps,EPSPowerShiftType *shift)
469: {
470: EPS_POWER *power = (EPS_POWER *)eps->data;
472: *shift = power->shift_type;
473: return(0);
474: }
479: /*@C
480: EPSPowerGetShiftType - Gets the type of shifts used during the power
481: iteration.
483: Collective on EPS
485: Input Parameter:
486: . eps - the eigenproblem solver context
488: Input Parameter:
489: . shift - the type of shift
491: Level: advanced
493: .seealso: EPSSetShiftType()
494: @*/
495: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
496: {
497: PetscErrorCode ierr, (*f)(EPS,EPSPowerShiftType*);
501: PetscObjectQueryFunction((PetscObject)eps,"EPSPowerGetShiftType_C",(void (**)())&f);
502: if (f) {
503: (*f)(eps,shift);
504: }
505: return(0);
506: }
510: PetscErrorCode EPSView_POWER(EPS eps,PetscViewer viewer)
511: {
513: EPS_POWER *power = (EPS_POWER *)eps->data;
514: PetscTruth isascii;
515: const char *shift_list[3] = { "constant", "rayleigh", "wilkinson" };
518: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
519: if (!isascii) {
520: SETERRQ1(1,"Viewer type %s not supported for EPSPOWER",((PetscObject)viewer)->type_name);
521: }
522: PetscViewerASCIIPrintf(viewer,"shift type: %s\n",shift_list[power->shift_type]);
523: return(0);
524: }
529: PetscErrorCode EPSCreate_POWER(EPS eps)
530: {
532: EPS_POWER *power;
535: PetscNew(EPS_POWER,&power);
536: PetscLogObjectMemory(eps,sizeof(EPS_POWER));
537: eps->data = (void *) power;
538: eps->ops->solve = EPSSolve_POWER;
539: eps->ops->solvets = EPSSolve_TS_POWER;
540: eps->ops->setup = EPSSetUp_POWER;
541: eps->ops->setfromoptions = EPSSetFromOptions_POWER;
542: eps->ops->destroy = EPSDestroy_Default;
543: eps->ops->view = EPSView_POWER;
544: eps->ops->backtransform = EPSBackTransform_POWER;
545: eps->ops->computevectors = EPSComputeVectors_Default;
546: power->shift_type = EPSPOWER_SHIFT_CONSTANT;
547: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerSetShiftType_C","EPSPowerSetShiftType_POWER",EPSPowerSetShiftType_POWER);
548: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPowerGetShiftType_C","EPSPowerGetShiftType_POWER",EPSPowerGetShiftType_POWER);
549: return(0);
550: }