Actual source code: power.c
slepc-3.9.1 2018-05-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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 http://slepc.upv.es.
31: */
33: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
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);
42: typedef struct {
43: EPSPowerShiftType shift_type;
44: PetscBool nonlinear;
45: PetscBool update;
46: SNES snes;
47: PetscErrorCode (*formFunctionB)(SNES,Vec,Vec,void*);
48: void *formFunctionBctx;
49: PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
50: void *formFunctionActx;
51: PetscErrorCode (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
52: } EPS_POWER;
54: PetscErrorCode EPSSetUp_Power(EPS eps)
55: {
57: EPS_POWER *power = (EPS_POWER*)eps->data;
58: PetscBool flg,istrivial;
59: STMatMode mode;
60: Mat A,B;
61: Vec res;
62: PetscContainer container;
63: PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
64: PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
65: void *ctx;
66: SNESLineSearch linesearch;
69: if (eps->ncv) {
70: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
71: } else eps->ncv = eps->nev;
72: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
73: if (!eps->max_it) eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
74: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
75: if (eps->which!=EPS_LARGEST_MAGNITUDE && eps->which !=EPS_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
76: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
77: if (power->nonlinear) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
78: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
79: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts only allowed in shift-and-invert or Cayley ST");
80: STGetMatMode(eps->st,&mode);
81: if (mode == ST_MATMODE_INPLACE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
82: }
83: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
84: if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in this solver");
85: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
86: RGIsTrivial(eps->rg,&istrivial);
87: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
88: EPSAllocateSolution(eps,0);
89: EPS_SetInnerProduct(eps);
91: if (power->nonlinear) {
92: if (eps->nev>1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
93: EPSSetWorkVecs(eps,4);
95: /* set up SNES solver */
96: if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
97: EPSGetOperators(eps,&A,&B);
99: PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA);
100: if (!formFunctionA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
101: PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container);
102: if (container) {
103: PetscContainerGetPointer(container,&ctx);
104: } else ctx = NULL;
105: MatCreateVecs(A,&res,NULL);
106: power->formFunctionA = formFunctionA;
107: power->formFunctionActx = ctx;
108: if (power->update) {
109: SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx);
110: PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB);
111: }
112: else { SNESSetFunction(power->snes,res,formFunctionA,ctx); }
113: VecDestroy(&res);
115: PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA);
116: if (!formJacobianA) SETERRQ(PetscObjectComm((PetscObject)eps),1,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
117: PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container);
118: if (container) {
119: PetscContainerGetPointer(container,&ctx);
120: } else ctx = NULL;
121: SNESSetJacobian(power->snes,A,A,formJacobianA,ctx);
122: SNESSetFromOptions(power->snes);
123: SNESGetLineSearch(power->snes,&linesearch);
124: if (power->update) {
125: SNESLineSearchSetPostCheck(linesearch,SNESLineSearchPostheckFunction,ctx);
126: }
127: SNESSetUp(power->snes);
128: if (B) {
129: PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB);
130: PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container);
131: if (power->formFunctionB && container) {
132: PetscContainerGetPointer(container,&power->formFunctionBctx);
133: } else power->formFunctionBctx = NULL;
134: }
135: } else {
136: EPSSetWorkVecs(eps,2);
137: DSSetType(eps->ds,DSNHEP);
138: DSAllocate(eps->ds,eps->nev);
139: }
140: return(0);
141: }
144: /*
145: Normalize a vector x with respect to a given norm as well as the
146: sign of the first entry.
147: */
148: static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscScalar *sign)
149: {
150: PetscErrorCode ierr;
151: PetscScalar alpha = 1.0;
152: PetscMPIInt rank;
153: PetscReal absv;
154: const PetscScalar *xx;
157: MPI_Comm_rank(PetscObjectComm((PetscObject)x),&rank);
158: if (!rank) {
159: VecGetArrayRead(x,&xx);
160: absv = PetscAbsScalar(*xx);
161: if (absv>10*PETSC_MACHINE_EPSILON) alpha = *xx/absv;
162: VecRestoreArrayRead(x,&xx);
163: }
164: MPI_Bcast(&alpha,1,MPIU_SCALAR,0,PetscObjectComm((PetscObject)x));
165: if (sign) *sign = alpha;
166: alpha *= norm;
167: VecScale(x,1.0/alpha);
168: return(0);
169: }
171: static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
172: {
174: EPS_POWER *power = (EPS_POWER*)eps->data;
175: Mat B;
178: STResetMatrixState(eps->st);
179: EPSGetOperators(eps,NULL,&B);
180: if (B) {
181: if (power->formFunctionB) {
182: (*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx);
183: } else {
184: MatMult(B,x,Bx);
185: }
186: } else {
187: VecCopy(x,Bx);
188: }
189: return(0);
190: }
192: static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
193: {
195: EPS_POWER *power = (EPS_POWER*)eps->data;
196: Mat A;
199: STResetMatrixState(eps->st);
200: EPSGetOperators(eps,&A,NULL);
201: if (A) {
202: if (power->formFunctionA) {
203: (*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx);
204: } else {
205: MatMult(A,x,Ax);
206: }
207: } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem\n");
208: return(0);
209: }
211: static PetscErrorCode SNESLineSearchPostheckFunction(SNESLineSearch linesearch,Vec x,Vec y,Vec w,PetscBool *changed_y,PetscBool *changed_w,void *ctx)
212: {
214: SNES snes;
215: EPS eps;
216: Vec oldx;
219: SNESLineSearchGetSNES(linesearch,&snes);
220: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
221: if (!eps) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"No composed EPS");
222: oldx = eps->work[3];
223: VecCopy(x,oldx);
224: return(0);
225: }
227: static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
228: {
230: EPS eps;
231: PetscReal bx;
232: Vec Bx;
233: EPS_POWER *power;
236: PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps);
237: if (!eps) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"No composed EPS");
238: power = (EPS_POWER*)eps->data;
239: Bx = eps->work[2];
240: if (power->formFunctionAB) {
241: (*power->formFunctionAB)(snes,x,y,Bx,ctx);
242: } else {
243: EPSPowerUpdateFunctionA(eps,x,y);
244: EPSPowerUpdateFunctionB(eps,x,Bx);
245: }
246: VecNorm(Bx,NORM_2,&bx);
247: Normalize(Bx,bx,NULL);
248: VecAXPY(y,-1.0,Bx);
249: return(0);
250: }
252: /*
253: Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
254: */
255: static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
256: {
258: EPS_POWER *power = (EPS_POWER*)eps->data;
259: Vec Bx;
262: VecCopy(x,y);
263: if (power->update) {
264: SNESSolve(power->snes,NULL,y);
265: } else {
266: Bx = eps->work[2];
267: SNESSolve(power->snes,Bx,y);
268: }
269: return(0);
270: }
272: /*
273: Use nonlinear inverse power to compute an initial guess
274: */
275: static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
276: {
277: EPS powereps;
278: Mat A,B;
279: Vec v1,v2;
280: SNES snes;
281: DM dm,newdm;
285: EPSCreate(PetscObjectComm((PetscObject)eps),&powereps);
286: EPSGetOperators(eps,&A,&B);
287: EPSSetType(powereps,EPSPOWER);
288: EPSSetOperators(powereps,A,B);
289: EPSSetTolerances(powereps,1e-6,4);
290: EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix);
291: EPSAppendOptionsPrefix(powereps,"init_");
292: EPSSetProblemType(powereps,EPS_GNHEP);
293: EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE);
294: EPSPowerSetNonlinear(powereps,PETSC_TRUE);
295: STSetType(powereps->st,STSINVERT);
296: /* attach dm to initial solve */
297: EPSPowerGetSNES(eps,&snes);
298: SNESGetDM(snes,&dm);
299: /* use dmshell to temporarily store snes context */
300: DMCreate(PetscObjectComm((PetscObject)eps),&newdm);
301: DMSetType(newdm,DMSHELL);
302: DMSetUp(newdm);
303: DMCopyDMSNES(dm,newdm);
304: EPSPowerGetSNES(powereps,&snes);
305: SNESSetDM(snes,dm);
306: EPSSetFromOptions(powereps);
307: EPSSolve(powereps);
308: BVGetColumn(eps->V,0,&v2);
309: BVGetColumn(powereps->V,0,&v1);
310: VecCopy(v1,v2);
311: BVRestoreColumn(powereps->V,0,&v1);
312: BVRestoreColumn(eps->V,0,&v2);
313: EPSDestroy(&powereps);
314: /* restore context back to the old nonlinear solver */
315: DMCopyDMSNES(newdm,dm);
316: DMDestroy(&newdm);
317: return(0);
318: }
320: PetscErrorCode EPSSolve_Power(EPS eps)
321: {
322: #if defined(SLEPC_MISSING_LAPACK_LAEV2)
324: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"LAEV2 - Lapack routine is unavailable");
325: #else
326: PetscErrorCode ierr;
327: EPS_POWER *power = (EPS_POWER*)eps->data;
328: PetscInt k,ld;
329: Vec v,y,e,Bx;
330: Mat A;
331: KSP ksp;
332: PetscReal relerr,norm,norm1,rt1,rt2,cs1;
333: PetscScalar theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
334: PetscBool breakdown;
335: KSPConvergedReason reason;
336: SNESConvergedReason snesreason;
339: e = eps->work[0];
340: y = eps->work[1];
341: if (power->nonlinear) Bx = eps->work[2];
342: else Bx = NULL;
344: if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) { STGetKSP(eps->st,&ksp); }
345: if (power->nonlinear) {
346: PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps);
347: if (power->update) {
348: EPSPowerComputeInitialGuess_Update(eps);
349: }
350: } else {
351: DSGetLeadingDimension(eps->ds,&ld);
352: }
353: if (!power->update) {
354: EPSGetStartVector(eps,0,NULL);
355: }
356: if (power->nonlinear) {
357: BVGetColumn(eps->V,0,&v);
358: EPSPowerUpdateFunctionB(eps,v,Bx);
359: VecNorm(Bx,NORM_2,&norm);
360: Normalize(Bx,norm,NULL);
361: BVRestoreColumn(eps->V,0,&v);
362: }
364: STGetShift(eps->st,&sigma); /* original shift */
365: rho = sigma;
367: while (eps->reason == EPS_CONVERGED_ITERATING) {
368: eps->its++;
369: k = eps->nconv;
371: /* y = OP v */
372: BVGetColumn(eps->V,k,&v);
373: if (power->nonlinear) {
374: VecCopy(v,eps->work[3]);
375: EPSPowerApply_SNES(eps,v,y);
376: VecCopy(eps->work[3],v);
377: } else {
378: STApply(eps->st,v,y);
379: }
380: BVRestoreColumn(eps->V,k,&v);
382: /* purge previously converged eigenvectors */
383: if (power->nonlinear) {
384: EPSPowerUpdateFunctionB(eps,y,Bx);
385: VecNorm(Bx,NORM_2,&norm);
386: Normalize(Bx,norm,&sign);
387: } else {
388: DSGetArray(eps->ds,DS_MAT_A,&T);
389: BVSetActiveColumns(eps->V,0,k);
390: BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL);
391: }
393: /* theta = (v,y)_B */
394: BVSetActiveColumns(eps->V,k,k+1);
395: BVDotVec(eps->V,y,&theta);
396: if (!power->nonlinear) {
397: T[k+k*ld] = theta;
398: DSRestoreArray(eps->ds,DS_MAT_A,&T);
399: }
401: if (power->nonlinear) theta = norm*sign;
403: if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
405: /* approximate eigenvalue is the Rayleigh quotient */
406: eps->eigr[eps->nconv] = theta;
408: /* compute relative error as ||y-theta v||_2/|theta| */
409: VecCopy(y,e);
410: BVGetColumn(eps->V,k,&v);
411: VecAXPY(e,power->nonlinear?-1.0:-theta,v);
412: BVRestoreColumn(eps->V,k,&v);
413: VecNorm(e,NORM_2,&relerr);
414: relerr /= PetscAbsScalar(theta);
416: } else { /* RQI */
418: /* delta = ||y||_B */
419: delta = norm;
421: /* compute relative error */
422: if (rho == 0.0) relerr = PETSC_MAX_REAL;
423: else relerr = 1.0 / (norm*PetscAbsScalar(rho));
425: /* approximate eigenvalue is the shift */
426: eps->eigr[eps->nconv] = rho;
428: /* compute new shift */
429: if (relerr<eps->tol) {
430: rho = sigma; /* if converged, restore original shift */
431: STSetShift(eps->st,rho);
432: } else {
433: rho = rho + theta/(delta*delta); /* Rayleigh quotient R(v) */
434: if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
435: /* beta1 is the norm of the residual associated with R(v) */
436: BVGetColumn(eps->V,k,&v);
437: VecAXPY(v,-theta/(delta*delta),y);
438: BVRestoreColumn(eps->V,k,&v);
439: BVScaleColumn(eps->V,k,1.0/delta);
440: BVNormColumn(eps->V,k,NORM_2,&norm1);
441: beta1 = norm1;
443: /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
444: STGetMatrix(eps->st,0,&A);
445: BVGetColumn(eps->V,k,&v);
446: MatMult(A,v,e);
447: VecDot(v,e,&alpha2);
448: BVRestoreColumn(eps->V,k,&v);
449: alpha2 = alpha2 / (beta1 * beta1);
451: /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
452: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
453: PetscStackCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
454: PetscFPTrapPop();
455: if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
456: else rho = rt2;
457: }
458: /* update operator according to new shift */
459: KSPSetErrorIfNotConverged(ksp,PETSC_FALSE);
460: STSetShift(eps->st,rho);
461: KSPGetConvergedReason(ksp,&reason);
462: if (reason) {
463: PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n");
464: rho *= 1+PETSC_MACHINE_EPSILON;
465: STSetShift(eps->st,rho);
466: KSPGetConvergedReason(ksp,&reason);
467: if (reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Second factorization failed");
468: }
469: KSPSetErrorIfNotConverged(ksp,PETSC_TRUE);
470: }
471: }
472: eps->errest[eps->nconv] = relerr;
474: /* normalize */
475: if (!power->nonlinear) { Normalize(y,norm,NULL); }
476: BVInsertVec(eps->V,k,y);
478: if (power->update) {
479: SNESGetConvergedReason(power->snes,&snesreason);
480: }
481: /* if relerr<tol, accept eigenpair */
482: if (relerr<eps->tol || (power->update && snesreason > 0)) {
483: eps->nconv = eps->nconv + 1;
484: if (eps->nconv<eps->nev) {
485: EPSGetStartVector(eps,eps->nconv,&breakdown);
486: if (breakdown) {
487: eps->reason = EPS_DIVERGED_BREAKDOWN;
488: PetscInfo(eps,"Unable to generate more start vectors\n");
489: break;
490: }
491: }
492: }
493: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev));
494: (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
495: }
497: if (power->nonlinear) {
498: PetscObjectCompose((PetscObject)power->snes,"eps",NULL);
499: } else {
500: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
501: DSSetState(eps->ds,DS_STATE_RAW);
502: }
503: return(0);
504: #endif
505: }
508: PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
509: {
511: EPS_POWER *power = (EPS_POWER*)eps->data;
512: SNESConvergedReason snesreason;
515: if (power->update) {
516: SNESGetConvergedReason(power->snes,&snesreason);
517: if (snesreason < 0) {
518: *reason = EPS_DIVERGED_BREAKDOWN;
519: return(0);
520: }
521: }
522: EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx);
523: return(0);
524: }
526: PetscErrorCode EPSBackTransform_Power(EPS eps)
527: {
529: EPS_POWER *power = (EPS_POWER*)eps->data;
532: if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
533: else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) {
534: EPSBackTransform_Default(eps);
535: }
536: return(0);
537: }
539: PetscErrorCode EPSSetFromOptions_Power(PetscOptionItems *PetscOptionsObject,EPS eps)
540: {
541: PetscErrorCode ierr;
542: EPS_POWER *power = (EPS_POWER*)eps->data;
543: PetscBool flg,val;
544: EPSPowerShiftType shift;
547: PetscOptionsHead(PetscOptionsObject,"EPS Power Options");
549: PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg);
550: if (flg) { EPSPowerSetShiftType(eps,shift); }
552: PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg);
553: if (flg) { EPSPowerSetNonlinear(eps,val); }
555: PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg);
556: if (flg) { EPSPowerSetUpdate(eps,val); }
558: PetscOptionsTail();
559: return(0);
560: }
562: static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
563: {
564: EPS_POWER *power = (EPS_POWER*)eps->data;
567: switch (shift) {
568: case EPS_POWER_SHIFT_CONSTANT:
569: case EPS_POWER_SHIFT_RAYLEIGH:
570: case EPS_POWER_SHIFT_WILKINSON:
571: if (power->shift_type != shift) {
572: power->shift_type = shift;
573: eps->state = EPS_STATE_INITIAL;
574: }
575: break;
576: default:
577: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
578: }
579: return(0);
580: }
582: /*@
583: EPSPowerSetShiftType - Sets the type of shifts used during the power
584: iteration. This can be used to emulate the Rayleigh Quotient Iteration
585: (RQI) method.
587: Logically Collective on EPS
589: Input Parameters:
590: + eps - the eigenproblem solver context
591: - shift - the type of shift
593: Options Database Key:
594: . -eps_power_shift_type - Sets the shift type (either 'constant' or
595: 'rayleigh' or 'wilkinson')
597: Notes:
598: By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
599: is the simple power method (or inverse iteration if a shift-and-invert
600: transformation is being used).
602: A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
603: EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
604: a cubic converging method such as RQI.
606: Level: advanced
608: .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
609: @*/
610: PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
611: {
617: PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
618: return(0);
619: }
621: static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
622: {
623: EPS_POWER *power = (EPS_POWER*)eps->data;
626: *shift = power->shift_type;
627: return(0);
628: }
630: /*@
631: EPSPowerGetShiftType - Gets the type of shifts used during the power
632: iteration.
634: Not Collective
636: Input Parameter:
637: . eps - the eigenproblem solver context
639: Input Parameter:
640: . shift - the type of shift
642: Level: advanced
644: .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
645: @*/
646: PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
647: {
653: PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
654: return(0);
655: }
657: static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
658: {
659: EPS_POWER *power = (EPS_POWER*)eps->data;
662: if (power->nonlinear != nonlinear) {
663: power->nonlinear = nonlinear;
664: eps->useds = PetscNot(nonlinear);
665: eps->state = EPS_STATE_INITIAL;
666: }
667: return(0);
668: }
670: /*@
671: EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
673: Logically Collective on EPS
675: Input Parameters:
676: + eps - the eigenproblem solver context
677: - nonlinear - whether the problem is nonlinear or not
679: Options Database Key:
680: . -eps_power_nonlinear - Sets the nonlinear flag
682: Notes:
683: If this flag is set, the solver assumes that the problem is nonlinear,
684: that is, the operators that define the eigenproblem are not constant
685: matrices, but depend on the eigenvector: A(x)*x=lambda*B(x)*x. This is
686: different from the case of nonlinearity with respect to the eigenvalue
687: (use the NEP solver class for this kind of problems).
689: The way in which nonlinear operators are specified is very similar to
690: the case of PETSc's SNES solver. The difference is that the callback
691: functions are provided via composed functions "formFunction" and
692: "formJacobian" in each of the matrix objects passed as arguments of
693: EPSSetOperators(). The application context required for these functions
694: can be attached via a composed PetscContainer.
696: Level: advanced
698: .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
699: @*/
700: PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
701: {
707: PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
708: return(0);
709: }
711: static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
712: {
713: EPS_POWER *power = (EPS_POWER*)eps->data;
716: *nonlinear = power->nonlinear;
717: return(0);
718: }
720: /*@
721: EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
723: Not Collective
725: Input Parameter:
726: . eps - the eigenproblem solver context
728: Input Parameter:
729: . nonlinear - the nonlinear flag
731: Level: advanced
733: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
734: @*/
735: PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
736: {
742: PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
743: return(0);
744: }
746: static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
747: {
748: EPS_POWER *power = (EPS_POWER*)eps->data;
751: if (!power->nonlinear) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems \n");
752: power->update = update;
753: eps->state = EPS_STATE_INITIAL;
754: return(0);
755: }
757: /*@
758: EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
759: for nonlinear problems. This potentially has a better convergence.
761: Logically Collective on EPS
763: Input Parameters:
764: + eps - the eigenproblem solver context
765: - update - whether the residual is updated monolithically or not
767: Options Database Key:
768: . -eps_power_update - Sets the update flag
770: Level: advanced
772: .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
773: @*/
774: PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
775: {
781: PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
782: return(0);
783: }
785: static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
786: {
787: EPS_POWER *power = (EPS_POWER*)eps->data;
790: *update = power->update;
791: return(0);
792: }
794: /*@
795: EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
796: for nonlinear problems.
798: Not Collective
800: Input Parameter:
801: . eps - the eigenproblem solver context
803: Input Parameter:
804: . update - the update flag
806: Level: advanced
808: .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
809: @*/
810: PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
811: {
817: PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
818: return(0);
819: }
821: static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
822: {
824: EPS_POWER *power = (EPS_POWER*)eps->data;
827: PetscObjectReference((PetscObject)snes);
828: SNESDestroy(&power->snes);
829: power->snes = snes;
830: PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
831: eps->state = EPS_STATE_INITIAL;
832: return(0);
833: }
835: /*@
836: EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
837: eigenvalue solver (to be used in nonlinear inverse iteration).
839: Collective on EPS
841: Input Parameters:
842: + eps - the eigenvalue solver
843: - snes - the nonlinear solver object
845: Level: advanced
847: .seealso: EPSPowerGetSNES()
848: @*/
849: PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
850: {
857: PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
858: return(0);
859: }
861: static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
862: {
864: EPS_POWER *power = (EPS_POWER*)eps->data;
867: if (!power->snes) {
868: SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes);
869: PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1);
870: SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix);
871: SNESAppendOptionsPrefix(power->snes,"eps_power_");
872: PetscLogObjectParent((PetscObject)eps,(PetscObject)power->snes);
873: }
874: *snes = power->snes;
875: return(0);
876: }
878: /*@
879: EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
880: with the eigenvalue solver.
882: Not Collective
884: Input Parameter:
885: . eps - the eigenvalue solver
887: Output Parameter:
888: . snes - the nonlinear solver object
890: Level: advanced
892: .seealso: EPSPowerSetSNES()
893: @*/
894: PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
895: {
901: PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
902: return(0);
903: }
905: PetscErrorCode EPSReset_Power(EPS eps)
906: {
908: EPS_POWER *power = (EPS_POWER*)eps->data;
911: if (power->snes) { SNESReset(power->snes); }
912: return(0);
913: }
915: PetscErrorCode EPSDestroy_Power(EPS eps)
916: {
918: EPS_POWER *power = (EPS_POWER*)eps->data;
921: if (power->nonlinear) {
922: SNESDestroy(&power->snes);
923: }
924: PetscFree(eps->data);
925: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL);
926: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL);
927: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL);
928: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL);
929: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL);
930: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL);
931: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL);
932: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL);
933: return(0);
934: }
936: PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
937: {
939: EPS_POWER *power = (EPS_POWER*)eps->data;
940: PetscBool isascii;
943: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
944: if (isascii) {
945: if (power->nonlinear) {
946: PetscViewerASCIIPrintf(viewer," using nonlinear inverse iteration\n");
947: if (power->update) {
948: PetscViewerASCIIPrintf(viewer," updating the residual monolithically\n");
949: }
950: if (!power->snes) { EPSPowerGetSNES(eps,&power->snes); }
951: SNESView(power->snes,viewer);
952: } else {
953: PetscViewerASCIIPrintf(viewer," %s shifts\n",EPSPowerShiftTypes[power->shift_type]);
954: }
955: }
956: return(0);
957: }
959: PetscErrorCode EPSComputeVectors_Power(EPS eps)
960: {
962: EPS_POWER *power = (EPS_POWER*)eps->data;
965: if (!power->nonlinear) {
966: EPSComputeVectors_Schur(eps);
967: }
968: return(0);
969: }
971: PetscErrorCode EPSSetDefaultST_Power(EPS eps)
972: {
974: EPS_POWER *power = (EPS_POWER*)eps->data;
975: KSP ksp;
976: PC pc;
979: if (power->nonlinear) {
980: STGetKSP(eps->st,&ksp);
981: KSPGetPC(ksp,&pc);
982: PCSetType(pc,PCNONE);
983: }
984: return(0);
985: }
987: PETSC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
988: {
989: EPS_POWER *ctx;
993: PetscNewLog(eps,&ctx);
994: eps->data = (void*)ctx;
996: eps->useds = PETSC_TRUE;
997: eps->categ = EPS_CATEGORY_OTHER;
999: eps->ops->solve = EPSSolve_Power;
1000: eps->ops->setup = EPSSetUp_Power;
1001: eps->ops->setfromoptions = EPSSetFromOptions_Power;
1002: eps->ops->reset = EPSReset_Power;
1003: eps->ops->destroy = EPSDestroy_Power;
1004: eps->ops->view = EPSView_Power;
1005: eps->ops->backtransform = EPSBackTransform_Power;
1006: eps->ops->computevectors = EPSComputeVectors_Power;
1007: eps->ops->setdefaultst = EPSSetDefaultST_Power;
1008: eps->stopping = EPSStopping_Power;
1010: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power);
1011: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power);
1012: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power);
1013: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power);
1014: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power);
1015: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power);
1016: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power);
1017: PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power);
1018: return(0);
1019: }