Actual source code: power.c
slepc-3.18.1 2022-11-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: PetscErrorCode EPSSolve_Power(EPS);
41: PetscErrorCode EPSSolve_TS_Power(EPS);
43: typedef struct {
44: EPSPowerShiftType shift_type;
45: PetscBool nonlinear;
46: PetscBool update;
47: SNES snes;
48: PetscErrorCode (*formFunctionB)(SNES,Vec,Vec,void*);
49: void *formFunctionBctx;
50: PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
51: void *formFunctionActx;
52: PetscErrorCode (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
53: PetscInt idx; /* index of the first nonzero entry in the iteration vector */
54: PetscMPIInt p; /* process id of the owner of idx */
55: PetscReal norm0; /* norm of initial vector */
56: } EPS_POWER;
58: static PetscErrorCode SNESMonitor_PowerUpdate(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
59: {
60: EPS eps;
62: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
64: /* Call EPS monitor on each SNES iteration */
65: EPSMonitor(eps,its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
66: return 0;
67: }
69: PetscErrorCode EPSSetUp_Power(EPS eps)
70: {
71: EPS_POWER *power = (EPS_POWER*)eps->data;
72: STMatMode mode;
73: Mat A,B,P;
74: Vec res;
75: PetscContainer container;
76: PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
77: PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
78: void *ctx;
80: if (eps->ncv!=PETSC_DEFAULT) {
82: } else eps->ncv = eps->nev;
83: if (eps->mpd!=PETSC_DEFAULT) PetscInfo(eps,"Warning: parameter mpd ignored\n");
84: if (eps->max_it==PETSC_DEFAULT) {
85: /* SNES will directly return the solution for us, and we need to do only one iteration */
86: if (power->nonlinear && power->update) eps->max_it = 1;
87: else eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
88: }
89: if (!eps->which) EPSSetWhichEigenpairs_Default(eps);
91: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
93: EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," (with variable shifts)");
94: STGetMatMode(eps->st,&mode);
96: }
97: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE);
98: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
99: EPSAllocateSolution(eps,0);
100: EPS_SetInnerProduct(eps);
102: if (power->nonlinear) {
104: EPSSetWorkVecs(eps,3);
107: /* Set up SNES solver */
108: /* If SNES was setup, we need to reset it so that it will be consistent with the current EPS */
109: if (power->snes) SNESReset(power->snes);
110: else EPSPowerGetSNES(eps,&power->snes);
112: /* For nonlinear solver (Newton), we should scale the initial vector back.
113: The initial vector will be scaled in EPSSetUp. */
114: if (eps->IS) VecNorm((eps->IS)[0],NORM_2,&power->norm0);
116: EPSGetOperators(eps,&A,&B);
118: PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA);
120: PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container);
121: if (container) PetscContainerGetPointer(container,&ctx);
122: else ctx = NULL;
123: MatCreateVecs(A,&res,NULL);
124: power->formFunctionA = formFunctionA;
125: power->formFunctionActx = ctx;
126: if (power->update) {
127: SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx);
128: PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB);
129: SNESMonitorSet(power->snes,SNESMonitor_PowerUpdate,NULL,NULL);
130: }
131: else SNESSetFunction(power->snes,res,formFunctionA,ctx);
132: VecDestroy(&res);
134: PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA);
136: PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container);
137: if (container) PetscContainerGetPointer(container,&ctx);
138: else ctx = NULL;
139: /* This allows users to compute a different preconditioning matrix than A.
140: * It is useful when A and B are shell matrices.
141: */
142: STGetPreconditionerMat(eps->st,&P);
143: /* If users do not provide a matrix, we simply use A */
144: SNESSetJacobian(power->snes,A,P? P:A,formJacobianA,ctx);
145: SNESSetFromOptions(power->snes);
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) PetscContainerGetPointer(container,&power->formFunctionBctx);
151: else power->formFunctionBctx = NULL;
152: }
153: } else {
154: if (eps->twosided) EPSSetWorkVecs(eps,3);
155: else EPSSetWorkVecs(eps,2);
156: DSSetType(eps->ds,DSNHEP);
157: DSAllocate(eps->ds,eps->nev);
158: }
159: /* dispatch solve method */
160: if (eps->twosided) {
163: eps->ops->solve = EPSSolve_TS_Power;
164: } else eps->ops->solve = EPSSolve_Power;
165: return 0;
166: }
168: /*
169: Find the index of the first nonzero entry of x, and the process that owns it.
170: */
171: static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscMPIInt *p)
172: {
173: PetscInt i,first,last,N;
174: PetscLayout map;
175: const PetscScalar *xx;
177: VecGetSize(x,&N);
178: VecGetOwnershipRange(x,&first,&last);
179: VecGetArrayRead(x,&xx);
180: for (i=first;i<last;i++) {
181: if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
182: }
183: VecRestoreArrayRead(x,&xx);
184: if (i==last) i=N;
185: MPIU_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x));
187: VecGetLayout(x,&map);
188: PetscLayoutFindOwner(map,*idx,p);
189: return 0;
190: }
192: /*
193: Normalize a vector x with respect to a given norm as well as the
194: sign of the first nonzero entry (assumed to be idx in process p).
195: */
196: static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscMPIInt p,PetscScalar *sign)
197: {
198: PetscScalar alpha=1.0;
199: PetscInt first,last;
200: PetscReal absv;
201: const PetscScalar *xx;
203: VecGetOwnershipRange(x,&first,&last);
204: if (idx>=first && idx<last) {
205: VecGetArrayRead(x,&xx);
206: absv = PetscAbsScalar(xx[idx-first]);
207: if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
208: VecRestoreArrayRead(x,&xx);
209: }
210: MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x));
211: if (sign) *sign = alpha;
212: alpha *= norm;
213: VecScale(x,1.0/alpha);
214: return 0;
215: }
217: static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
218: {
219: EPS_POWER *power = (EPS_POWER*)eps->data;
220: Mat B;
222: STResetMatrixState(eps->st);
223: EPSGetOperators(eps,NULL,&B);
224: if (B) {
225: if (power->formFunctionB) (*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx);
226: else MatMult(B,x,Bx);
227: } else VecCopy(x,Bx);
228: return 0;
229: }
231: static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
232: {
233: EPS_POWER *power = (EPS_POWER*)eps->data;
234: Mat A;
236: STResetMatrixState(eps->st);
237: EPSGetOperators(eps,&A,NULL);
239: if (power->formFunctionA) (*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx);
240: else MatMult(A,x,Ax);
241: return 0;
242: }
244: static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
245: {
246: EPS eps;
247: PetscReal bx;
248: Vec Bx;
249: PetscScalar sign;
250: EPS_POWER *power;
252: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
254: STResetMatrixState(eps->st);
255: power = (EPS_POWER*)eps->data;
256: Bx = eps->work[2];
257: if (power->formFunctionAB) (*power->formFunctionAB)(snes,x,y,Bx,ctx);
258: else {
259: EPSPowerUpdateFunctionA(eps,x,y);
260: EPSPowerUpdateFunctionB(eps,x,Bx);
261: }
262: VecNorm(Bx,NORM_2,&bx);
263: Normalize(Bx,bx,power->idx,power->p,&sign);
264: VecAXPY(y,-1.0,Bx);
265: /* Keep tracking eigenvalue update. It would be useful when we want to monitor solver progress via snes monitor. */
266: eps->eigr[(eps->nconv < eps->nev)? eps->nconv:(eps->nconv-1)] = 1.0/(bx*sign);
267: return 0;
268: }
270: /*
271: Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
272: */
273: static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
274: {
275: EPS_POWER *power = (EPS_POWER*)eps->data;
276: Vec Bx;
278: VecCopy(x,y);
279: if (power->update) SNESSolve(power->snes,NULL,y);
280: else {
281: Bx = eps->work[2];
282: SNESSolve(power->snes,Bx,y);
283: }
284: return 0;
285: }
287: /*
288: Use nonlinear inverse power to compute an initial guess
289: */
290: static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
291: {
292: EPS powereps;
293: Mat A,B,P;
294: Vec v1,v2;
295: SNES snes;
296: DM dm,newdm;
298: EPSCreate(PetscObjectComm((PetscObject)eps),&powereps);
299: EPSGetOperators(eps,&A,&B);
300: EPSSetType(powereps,EPSPOWER);
301: EPSSetOperators(powereps,A,B);
302: EPSSetTolerances(powereps,1e-6,4);
303: EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix);
304: EPSAppendOptionsPrefix(powereps,"init_");
305: EPSSetProblemType(powereps,EPS_GNHEP);
306: EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE);
307: EPSPowerSetNonlinear(powereps,PETSC_TRUE);
308: STGetPreconditionerMat(eps->st,&P);
309: /* attach dm to initial solve */
310: EPSPowerGetSNES(eps,&snes);
311: SNESGetDM(snes,&dm);
312: /* use dmshell to temporarily store snes context */
313: DMCreate(PetscObjectComm((PetscObject)eps),&newdm);
314: DMSetType(newdm,DMSHELL);
315: DMSetUp(newdm);
316: DMCopyDMSNES(dm,newdm);
317: EPSPowerGetSNES(powereps,&snes);
318: SNESSetDM(snes,dm);
319: EPSSetFromOptions(powereps);
320: if (P) STSetPreconditionerMat(powereps->st,P);
321: EPSSolve(powereps);
322: BVGetColumn(eps->V,0,&v2);
323: BVGetColumn(powereps->V,0,&v1);
324: VecCopy(v1,v2);
325: BVRestoreColumn(powereps->V,0,&v1);
326: BVRestoreColumn(eps->V,0,&v2);
327: EPSDestroy(&powereps);
328: /* restore context back to the old nonlinear solver */
329: DMCopyDMSNES(newdm,dm);
330: DMDestroy(&newdm);
331: return 0;
332: }
334: PetscErrorCode EPSSolve_Power(EPS eps)
335: {
336: EPS_POWER *power = (EPS_POWER*)eps->data;
337: PetscInt k,ld;
338: Vec v,y,e,Bx;
339: Mat A;
340: KSP ksp;
341: PetscReal relerr,norm,norm1,rt1,rt2,cs1;
342: PetscScalar theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
343: PetscBool breakdown;
344: KSPConvergedReason reason;
345: SNESConvergedReason snesreason;
347: e = eps->work[0];
348: y = eps->work[1];
349: if (power->nonlinear) Bx = eps->work[2];
350: else Bx = NULL;
352: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) STGetKSP(eps->st,&ksp);
353: if (power->nonlinear) {
354: PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps);
355: /* Compute an initial guess only when users do not provide one */
356: if (power->update && !eps->nini) EPSPowerComputeInitialGuess_Update(eps);
357: } else DSGetLeadingDimension(eps->ds,&ld);
358: if (!power->update) EPSGetStartVector(eps,0,NULL);
359: if (power->nonlinear) {
360: BVGetColumn(eps->V,0,&v);
361: if (eps->nini) {
362: /* We scale the initial vector back if the initial vector was provided by users */
363: VecScale(v,power->norm0);
364: }
365: EPSPowerUpdateFunctionB(eps,v,Bx);
366: VecNorm(Bx,NORM_2,&norm);
367: FirstNonzeroIdx(Bx,&power->idx,&power->p);
368: Normalize(Bx,norm,power->idx,power->p,NULL);
369: BVRestoreColumn(eps->V,0,&v);
370: }
372: STGetShift(eps->st,&sigma); /* original shift */
373: rho = sigma;
375: while (eps->reason == EPS_CONVERGED_ITERATING) {
376: eps->its++;
377: k = eps->nconv;
379: /* y = OP v */
380: BVGetColumn(eps->V,k,&v);
381: if (power->nonlinear) EPSPowerApply_SNES(eps,v,y);
382: else STApply(eps->st,v,y);
383: BVRestoreColumn(eps->V,k,&v);
385: /* purge previously converged eigenvectors */
386: if (PetscUnlikely(power->nonlinear)) {
387: /* We do not need to call this for Newton eigensolver since eigenvalue is
388: * updated in function evaluations.
389: */
390: if (!power->update) {
391: EPSPowerUpdateFunctionB(eps,y,Bx);
392: VecNorm(Bx,NORM_2,&norm);
393: Normalize(Bx,norm,power->idx,power->p,&sign);
394: }
395: } else {
396: DSGetArray(eps->ds,DS_MAT_A,&T);
397: BVSetActiveColumns(eps->V,0,k);
398: BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL);
399: }
401: /* theta = (v,y)_B */
402: BVSetActiveColumns(eps->V,k,k+1);
403: BVDotVec(eps->V,y,&theta);
404: if (!power->nonlinear) {
405: T[k+k*ld] = theta;
406: DSRestoreArray(eps->ds,DS_MAT_A,&T);
407: }
409: /* Eigenvalue is already stored in function evaluations.
410: * Assign eigenvalue to theta to make the rest of the code consistent
411: */
412: if (power->update) theta = eps->eigr[eps->nconv];
413: else if (power->nonlinear) theta = 1.0/norm*sign; /* Eigenvalue: 1/|Bx| */
415: if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
417: /* approximate eigenvalue is the Rayleigh quotient */
418: eps->eigr[eps->nconv] = theta;
420: /**
421: * If the Newton method (update, SNES) is used, we do not compute "relerr"
422: * since SNES determines the convergence.
423: */
424: if (PetscUnlikely(power->update)) relerr = 0.;
425: else {
426: /* compute relative error as ||y-theta v||_2/|theta| */
427: VecCopy(y,e);
428: BVGetColumn(eps->V,k,&v);
429: VecAXPY(e,power->nonlinear?-1.0:-theta,v);
430: BVRestoreColumn(eps->V,k,&v);
431: VecNorm(e,NORM_2,&relerr);
432: if (PetscUnlikely(power->nonlinear)) relerr *= PetscAbsScalar(theta);
433: else relerr /= PetscAbsScalar(theta);
434: }
436: } else { /* RQI */
438: /* delta = ||y||_B */
439: delta = norm;
441: /* compute relative error */
442: if (rho == 0.0) relerr = PETSC_MAX_REAL;
443: else relerr = 1.0 / (norm*PetscAbsScalar(rho));
445: /* approximate eigenvalue is the shift */
446: eps->eigr[eps->nconv] = rho;
448: /* compute new shift */
449: if (relerr<eps->tol) {
450: rho = sigma; /* if converged, restore original shift */
451: STSetShift(eps->st,rho);
452: } else {
453: rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
454: if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
455: /* beta1 is the norm of the residual associated with R(v) */
456: BVGetColumn(eps->V,k,&v);
457: VecAXPY(v,-PetscConj(theta)/(delta*delta),y);
458: BVRestoreColumn(eps->V,k,&v);
459: BVScaleColumn(eps->V,k,1.0/delta);
460: BVNormColumn(eps->V,k,NORM_2,&norm1);
461: beta1 = norm1;
463: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
464: STGetMatrix(eps->st,0,&A);
465: BVGetColumn(eps->V,k,&v);
466: MatMult(A,v,e);
467: VecDot(v,e,&alpha2);
468: BVRestoreColumn(eps->V,k,&v);
469: alpha2 = alpha2 / (beta1 * beta1);
471: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
472: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
473: PetscCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
474: PetscFPTrapPop();
475: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
476: else rho = rt2;
477: }
478: /* update operator according to new shift */
479: KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
480: STSetShift(eps->st,rho);
481: KSPGetConvergedReason(ksp,&reason);
482: if (reason) {
483: PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
484: rho *= 1+10*PETSC_MACHINE_EPSILON;
485: STSetShift(eps->st,rho);
486: KSPGetConvergedReason(ksp,&reason);
488: }
489: KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
490: }
491: }
492: eps->errest[eps->nconv] = relerr;
494: /* normalize */
495: if (!power->nonlinear) Normalize(y,norm,power->idx,power->p,NULL);
496: BVInsertVec(eps->V,k,y);
498: if (PetscUnlikely(power->update)) {
499: SNESGetConvergedReason(power->snes,&snesreason);
500: /* For Newton eigensolver, we are ready to return once SNES converged. */
501: if (snesreason>0) eps->nconv = 1;
502: } else if (PetscUnlikely(relerr<eps->tol)) { /* accept eigenpair */
503: eps->nconv = eps->nconv + 1;
504: if (eps->nconv<eps->nev) {
505: EPSGetStartVector(eps,eps->nconv,&breakdown);
506: if (breakdown) {
507: eps->reason = EPS_DIVERGED_BREAKDOWN;
508: PetscInfo(eps,"Unable to generate more start vectors\n");
509: break;
510: }
511: }
512: }
513: /* For Newton eigensolver, monitor will be called from SNES monitor */
514: if (!power->update) EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
515: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
517: /**
518: * When a customized stopping test is used, and reason can be set to be converged (EPS_CONVERGED_USER).
519: * In this case, we need to increase eps->nconv to "1" so users can retrieve the solution.
520: */
521: if (PetscUnlikely(power->nonlinear && eps->reason>0)) eps->nconv = 1;
522: }
524: if (power->nonlinear) PetscObjectCompose((PetscObject)power->snes,"eps",NULL);
525: else {
526: DSSetDimensions(eps->ds,eps->nconv,0,0);
527: DSSetState(eps->ds,DS_STATE_RAW);
528: }
529: return 0;
530: }
532: PetscErrorCode EPSSolve_TS_Power(EPS eps)
533: {
534: EPS_POWER *power = (EPS_POWER*)eps->data;
535: PetscInt k,ld;
536: Vec v,w,y,e,z;
537: KSP ksp;
538: PetscReal relerr=1.0,relerrl,delta;
539: PetscScalar theta,rho,alpha,sigma;
540: PetscBool breakdown,breakdownl;
541: KSPConvergedReason reason;
543: e = eps->work[0];
544: v = eps->work[1];
545: w = eps->work[2];
547: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) STGetKSP(eps->st,&ksp);
548: DSGetLeadingDimension(eps->ds,&ld);
549: EPSGetStartVector(eps,0,NULL);
550: EPSGetLeftStartVector(eps,0,NULL);
551: BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL);
552: BVCopyVec(eps->V,0,v);
553: BVCopyVec(eps->W,0,w);
554: STGetShift(eps->st,&sigma); /* original shift */
555: rho = sigma;
557: while (eps->reason == EPS_CONVERGED_ITERATING) {
558: eps->its++;
559: k = eps->nconv;
561: /* y = OP v, z = OP' w */
562: BVGetColumn(eps->V,k,&y);
563: STApply(eps->st,v,y);
564: BVRestoreColumn(eps->V,k,&y);
565: BVGetColumn(eps->W,k,&z);
566: STApplyHermitianTranspose(eps->st,w,z);
567: BVRestoreColumn(eps->W,k,&z);
569: /* purge previously converged eigenvectors */
570: BVBiorthogonalizeColumn(eps->V,eps->W,k);
572: /* theta = (w,y)_B */
573: BVSetActiveColumns(eps->V,k,k+1);
574: BVDotVec(eps->V,w,&theta);
575: theta = PetscConj(theta);
577: if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
579: /* approximate eigenvalue is the Rayleigh quotient */
580: eps->eigr[eps->nconv] = theta;
582: /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
583: BVCopyVec(eps->V,k,e);
584: VecAXPY(e,-theta,v);
585: VecNorm(e,NORM_2,&relerr);
586: BVCopyVec(eps->W,k,e);
587: VecAXPY(e,-PetscConj(theta),w);
588: VecNorm(e,NORM_2,&relerrl);
589: relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
590: }
592: /* normalize */
593: BVSetActiveColumns(eps->V,k,k+1);
594: BVGetColumn(eps->W,k,&z);
595: BVDotVec(eps->V,z,&alpha);
596: BVRestoreColumn(eps->W,k,&z);
597: delta = PetscSqrtReal(PetscAbsScalar(alpha));
599: BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta));
600: BVScaleColumn(eps->W,k,1.0/delta);
601: BVCopyVec(eps->V,k,v);
602: BVCopyVec(eps->W,k,w);
604: if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */
606: /* compute relative error */
607: if (rho == 0.0) relerr = PETSC_MAX_REAL;
608: else relerr = 1.0 / PetscAbsScalar(delta*rho);
610: /* approximate eigenvalue is the shift */
611: eps->eigr[eps->nconv] = rho;
613: /* compute new shift */
614: if (relerr<eps->tol) {
615: rho = sigma; /* if converged, restore original shift */
616: STSetShift(eps->st,rho);
617: } else {
618: rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
619: /* update operator according to new shift */
620: KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
621: STSetShift(eps->st,rho);
622: KSPGetConvergedReason(ksp,&reason);
623: if (reason) {
624: PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
625: rho *= 1+10*PETSC_MACHINE_EPSILON;
626: STSetShift(eps->st,rho);
627: KSPGetConvergedReason(ksp,&reason);
629: }
630: KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
631: }
632: }
633: eps->errest[eps->nconv] = relerr;
635: /* if relerr<tol, accept eigenpair */
636: if (relerr<eps->tol) {
637: eps->nconv = eps->nconv + 1;
638: if (eps->nconv<eps->nev) {
639: EPSGetStartVector(eps,eps->nconv,&breakdown);
640: EPSGetLeftStartVector(eps,eps->nconv,&breakdownl);
641: if (breakdown || breakdownl) {
642: eps->reason = EPS_DIVERGED_BREAKDOWN;
643: PetscInfo(eps,"Unable to generate more start vectors\n");
644: break;
645: }
646: BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL);
647: }
648: }
649: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
650: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
651: }
653: DSSetDimensions(eps->ds,eps->nconv,0,0);
654: DSSetState(eps->ds,DS_STATE_RAW);
655: return 0;
656: }
658: PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
659: {
660: EPS_POWER *power = (EPS_POWER*)eps->data;
661: SNESConvergedReason snesreason;
663: if (PetscUnlikely(power->update)) {
664: SNESGetConvergedReason(power->snes,&snesreason);
665: if (snesreason < 0) {
666: *reason = EPS_DIVERGED_BREAKDOWN;
667: return 0;
668: }
669: }
670: EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx);
671: return 0;
672: }
674: PetscErrorCode EPSBackTransform_Power(EPS eps)
675: {
676: EPS_POWER *power = (EPS_POWER*)eps->data;
678: if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
679: else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) EPSBackTransform_Default(eps);
680: return 0;
681: }
683: PetscErrorCode EPSSetFromOptions_Power(EPS eps,PetscOptionItems *PetscOptionsObject)
684: {
685: EPS_POWER *power = (EPS_POWER*)eps->data;
686: PetscBool flg,val;
687: EPSPowerShiftType shift;
689: PetscOptionsHeadBegin(PetscOptionsObject,"EPS Power Options");
691: PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
692: if (flg) EPSPowerSetShiftType(eps,shift);
694: PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg);
695: if (flg) EPSPowerSetNonlinear(eps,val);
697: PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg);
698: if (flg) EPSPowerSetUpdate(eps,val);
700: PetscOptionsHeadEnd();
701: return 0;
702: }
704: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
705: {
706: EPS_POWER *power = (EPS_POWER*)eps->data;
708: switch (shift) {
709: case EPS_POWER_SHIFT_CONSTANT:
710: case EPS_POWER_SHIFT_RAYLEIGH:
711: case EPS_POWER_SHIFT_WILKINSON:
712: if (power->shift_type != shift) {
713: power->shift_type = shift;
714: eps->state = EPS_STATE_INITIAL;
715: }
716: break;
717: default:
718: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
719: }
720: return 0;
721: }
723: /*@
724: EPSPowerSetShiftType - Sets the type of shifts used during the power
725: iteration. This can be used to emulate the Rayleigh Quotient Iteration
726: (RQI) method.
728: Logically Collective on eps
730: Input Parameters:
731: + eps - the eigenproblem solver context
732: - shift - the type of shift
734: Options Database Key:
735: . -eps_power_shift_type - Sets the shift type (either 'constant' or
736: 'rayleigh' or 'wilkinson')
738: Notes:
739: By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
740: is the simple power method (or inverse iteration if a shift-and-invert
741: transformation is being used).
743: A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
744: EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
745: a cubic converging method such as RQI.
747: Level: advanced
749: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
750: @*/
751: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
752: {
755: PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
756: return 0;
757: }
759: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
760: {
761: EPS_POWER *power = (EPS_POWER*)eps->data;
763: *shift = power->shift_type;
764: return 0;
765: }
767: /*@
768: EPSPowerGetShiftType - Gets the type of shifts used during the power
769: iteration.
771: Not Collective
773: Input Parameter:
774: . eps - the eigenproblem solver context
776: Output Parameter:
777: . shift - the type of shift
779: Level: advanced
781: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
782: @*/
783: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
784: {
787: PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
788: return 0;
789: }
791: static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
792: {
793: EPS_POWER *power = (EPS_POWER*)eps->data;
795: if (power->nonlinear != nonlinear) {
796: power->nonlinear = nonlinear;
797: eps->useds = PetscNot(nonlinear);
798: eps->ops->setupsort = nonlinear? NULL: EPSSetUpSort_Default;
799: eps->state = EPS_STATE_INITIAL;
800: }
801: return 0;
802: }
804: /*@
805: EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
807: Logically Collective on eps
809: Input Parameters:
810: + eps - the eigenproblem solver context
811: - nonlinear - whether the problem is nonlinear or not
813: Options Database Key:
814: . -eps_power_nonlinear - Sets the nonlinear flag
816: Notes:
817: If this flag is set, the solver assumes that the problem is nonlinear,
818: that is, the operators that define the eigenproblem are not constant
819: matrices, but depend on the eigenvector, A(x)*x=lambda*B(x)*x. This is
820: different from the case of nonlinearity with respect to the eigenvalue
821: (use the NEP solver class for this kind of problems).
823: The way in which nonlinear operators are specified is very similar to
824: the case of PETSc's SNES solver. The difference is that the callback
825: functions are provided via composed functions "formFunction" and
826: "formJacobian" in each of the matrix objects passed as arguments of
827: EPSSetOperators(). The application context required for these functions
828: can be attached via a composed PetscContainer.
830: Level: advanced
832: .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
833: @*/
834: PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
835: {
838: PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
839: return 0;
840: }
842: static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
843: {
844: EPS_POWER *power = (EPS_POWER*)eps->data;
846: *nonlinear = power->nonlinear;
847: return 0;
848: }
850: /*@
851: EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
853: Not Collective
855: Input Parameter:
856: . eps - the eigenproblem solver context
858: Output Parameter:
859: . nonlinear - the nonlinear flag
861: Level: advanced
863: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
864: @*/
865: PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
866: {
869: PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
870: return 0;
871: }
873: static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
874: {
875: EPS_POWER *power = (EPS_POWER*)eps->data;
878: power->update = update;
879: eps->state = EPS_STATE_INITIAL;
880: return 0;
881: }
883: /*@
884: EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
885: for nonlinear problems. This potentially has a better convergence.
887: Logically Collective on eps
889: Input Parameters:
890: + eps - the eigenproblem solver context
891: - update - whether the residual is updated monolithically or not
893: Options Database Key:
894: . -eps_power_update - Sets the update flag
896: Level: advanced
898: .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
899: @*/
900: PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
901: {
904: PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
905: return 0;
906: }
908: static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
909: {
910: EPS_POWER *power = (EPS_POWER*)eps->data;
912: *update = power->update;
913: return 0;
914: }
916: /*@
917: EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
918: for nonlinear problems.
920: Not Collective
922: Input Parameter:
923: . eps - the eigenproblem solver context
925: Output Parameter:
926: . update - the update flag
928: Level: advanced
930: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
931: @*/
932: PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
933: {
936: PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
937: return 0;
938: }
940: static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
941: {
942: EPS_POWER *power = (EPS_POWER*)eps->data;
944: PetscObjectReference((PetscObject)snes);
945: SNESDestroy(&power->snes);
946: power->snes = snes;
947: eps->state = EPS_STATE_INITIAL;
948: return 0;
949: }
951: /*@
952: EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
953: eigenvalue solver (to be used in nonlinear inverse iteration).
955: Collective on eps
957: Input Parameters:
958: + eps - the eigenvalue solver
959: - snes - the nonlinear solver object
961: Level: advanced
963: .seealso: EPSPowerGetSNES()
964: @*/
965: PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
966: {
970: PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
971: return 0;
972: }
974: static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
975: {
976: EPS_POWER *power = (EPS_POWER*)eps->data;
978: if (!power->snes) {
979: SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes);
980: PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1);
981: SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix);
982: SNESAppendOptionsPrefix(power->snes,"eps_power_");
983: PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options);
984: }
985: *snes = power->snes;
986: return 0;
987: }
989: /*@
990: EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
991: with the eigenvalue solver.
993: Not Collective
995: Input Parameter:
996: . eps - the eigenvalue solver
998: Output Parameter:
999: . snes - the nonlinear solver object
1001: Level: advanced
1003: .seealso: EPSPowerSetSNES()
1004: @*/
1005: PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
1006: {
1009: PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
1010: return 0;
1011: }
1013: PetscErrorCode EPSReset_Power(EPS eps)
1014: {
1015: EPS_POWER *power = (EPS_POWER*)eps->data;
1017: if (power->snes) SNESReset(power->snes);
1018: return 0;
1019: }
1021: PetscErrorCode EPSDestroy_Power(EPS eps)
1022: {
1023: EPS_POWER *power = (EPS_POWER*)eps->data;
1025: if (power->nonlinear) SNESDestroy(&power->snes);
1026: PetscFree(eps->data);
1027: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
1028: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
1029: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL);
1030: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL);
1031: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL);
1032: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL);
1033: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL);
1034: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL);
1035: return 0;
1036: }
1038: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
1039: {
1040: EPS_POWER *power = (EPS_POWER*)eps->data;
1041: PetscBool isascii;
1043: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1044: if (isascii) {
1045: if (power->nonlinear) {
1046: PetscViewerASCIIPrintf(viewer," using nonlinear inverse iteration\n");
1047: if (power->update) PetscViewerASCIIPrintf(viewer," updating the residual monolithically\n");
1048: if (!power->snes) EPSPowerGetSNES(eps,&power->snes);
1049: PetscViewerASCIIPushTab(viewer);
1050: SNESView(power->snes,viewer);
1051: PetscViewerASCIIPopTab(viewer);
1052: } else PetscViewerASCIIPrintf(viewer," %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
1053: }
1054: return 0;
1055: }
1057: PetscErrorCode EPSComputeVectors_Power(EPS eps)
1058: {
1059: EPS_POWER *power = (EPS_POWER*)eps->data;
1061: if (eps->twosided) {
1062: EPSComputeVectors_Twosided(eps);
1063: BVNormalize(eps->V,NULL);
1064: BVNormalize(eps->W,NULL);
1065: } else if (!power->nonlinear) EPSComputeVectors_Schur(eps);
1066: return 0;
1067: }
1069: PetscErrorCode EPSSetDefaultST_Power(EPS eps)
1070: {
1071: EPS_POWER *power = (EPS_POWER*)eps->data;
1072: KSP ksp;
1073: PC pc;
1075: if (power->nonlinear) {
1076: eps->categ=EPS_CATEGORY_PRECOND;
1077: STGetKSP(eps->st,&ksp);
1078: /* Set ST as STPRECOND so it can carry one preconditioning matrix
1079: * It is useful when A and B are shell matrices
1080: */
1081: STSetType(eps->st,STPRECOND);
1082: KSPGetPC(ksp,&pc);
1083: PCSetType(pc,PCNONE);
1084: }
1085: return 0;
1086: }
1088: SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
1089: {
1090: EPS_POWER *ctx;
1092: PetscNew(&ctx);
1093: eps->data = (void*)ctx;
1095: eps->useds = PETSC_TRUE;
1096: eps->categ = EPS_CATEGORY_OTHER;
1098: eps->ops->setup = EPSSetUp_Power;
1099: eps->ops->setupsort = EPSSetUpSort_Default;
1100: eps->ops->setfromoptions = EPSSetFromOptions_Power;
1101: eps->ops->reset = EPSReset_Power;
1102: eps->ops->destroy = EPSDestroy_Power;
1103: eps->ops->view = EPSView_Power;
1104: eps->ops->backtransform = EPSBackTransform_Power;
1105: eps->ops->computevectors = EPSComputeVectors_Power;
1106: eps->ops->setdefaultst = EPSSetDefaultST_Power;
1107: eps->stopping = EPSStopping_Power;
1109: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
1110: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
1111: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power);
1112: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power);
1113: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power);
1114: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power);
1115: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power);
1116: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power);
1117: return 0;
1118: }