Actual source code: slp.c
1: /*
3: SLEPc nonlinear eigensolver: "slp"
5: Method: Succesive linear problems
7: Algorithm:
9: Newton-type iteration based on first order Taylor approximation.
11: References:
13: [1] A. Ruhe, "Algorithms for the nonlinear eigenvalue problem", SIAM J.
14: Numer. Anal. 10(4):674-689, 1973.
16: Last update: Feb 2013
18: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
22: This file is part of SLEPc.
24: SLEPc is free software: you can redistribute it and/or modify it under the
25: terms of version 3 of the GNU Lesser General Public License as published by
26: the Free Software Foundation.
28: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
29: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
30: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
31: more details.
33: You should have received a copy of the GNU Lesser General Public License
34: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: */
38: #include <slepc-private/nepimpl.h> /*I "slepcnep.h" I*/
39: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
41: typedef struct {
42: EPS eps; /* linear eigensolver for T*z = mu*Tp*z */
43: PetscBool setfromoptionscalled;
44: } NEP_SLP;
48: PetscErrorCode NEPSetUp_SLP(NEP nep)
49: {
51: NEP_SLP *ctx = (NEP_SLP*)nep->data;
52: ST st;
55: if (nep->ncv) { /* ncv set */
56: if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
57: } else if (nep->mpd) { /* mpd set */
58: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
59: } else { /* neither set: defaults depend on nev being small or large */
60: if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
61: else {
62: nep->mpd = 500;
63: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
64: }
65: }
66: if (!nep->mpd) nep->mpd = nep->ncv;
67: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
68: if (nep->nev>1) { PetscInfo(nep,"Warning: requested more than one eigenpair but SLP can only compute one\n"); }
69: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
70: if (!nep->max_funcs) nep->max_funcs = nep->max_it;
72: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
73: EPSSetWhichEigenpairs(ctx->eps,EPS_TARGET_MAGNITUDE);
74: EPSSetTarget(ctx->eps,0.0);
75: EPSGetST(ctx->eps,&st);
76: STSetType(st,STSINVERT);
77: EPSSetDimensions(ctx->eps,1,nep->ncv,nep->mpd);
78: EPSSetTolerances(ctx->eps,nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:nep->rtol/10.0,nep->max_it);
79: if (ctx->setfromoptionscalled) {
80: EPSSetFromOptions(ctx->eps);
81: ctx->setfromoptionscalled = PETSC_FALSE;
82: }
84: NEPAllocateSolution(nep);
85: NEPSetWorkVecs(nep,1);
86: return(0);
87: }
91: PetscErrorCode NEPSolve_SLP(NEP nep)
92: {
94: NEP_SLP *ctx = (NEP_SLP*)nep->data;
95: Mat T=nep->function,Tp=nep->jacobian;
96: Vec u=nep->V[0],r=nep->work[0];
97: PetscScalar lambda,mu,im;
98: PetscReal relerr;
99: PetscInt nconv;
100: MatStructure mats;
103: /* get initial approximation of eigenvalue and eigenvector */
104: NEPGetDefaultShift(nep,&lambda);
105: if (!nep->nini) {
106: SlepcVecSetRandom(u,nep->rand);
107: }
109: /* Restart loop */
110: while (nep->reason == NEP_CONVERGED_ITERATING) {
111: nep->its++;
113: /* evaluate T(lambda) and T'(lambda) */
114: NEPComputeFunction(nep,lambda,&T,&T,&mats);
115: NEPComputeJacobian(nep,lambda,&Tp,&mats);
117: /* form residual, r = T(lambda)*u (used in convergence test only) */
118: MatMult(T,u,r);
120: /* convergence test */
121: VecNorm(r,NORM_2,&relerr);
122: nep->errest[nep->nconv] = relerr;
123: nep->eig[nep->nconv] = lambda;
124: if (relerr<=nep->rtol) {
125: nep->nconv = nep->nconv + 1;
126: nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
127: }
128: NEPMonitor(nep,nep->its,nep->nconv,nep->eig,nep->errest,1);
130: if (!nep->nconv) {
131: /* compute eigenvalue correction mu and eigenvector approximation u */
132: EPSSetOperators(ctx->eps,T,Tp);
133: EPSSetInitialSpace(ctx->eps,1,&u);
134: EPSSolve(ctx->eps);
135: EPSGetConverged(ctx->eps,&nconv);
136: if (!nconv) {
137: PetscInfo1(nep,"iter=%D, inner iteration failed, stopping solve\n",nep->its);
138: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
139: break;
140: }
141: EPSGetEigenpair(ctx->eps,0,&mu,&im,u,NULL);
142: if (PetscAbsScalar(im)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Complex eigenvalue approximation - not implemented in real scalars");
144: /* correct eigenvalue */
145: lambda = lambda - mu;
146: }
147: if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
148: }
149: return(0);
150: }
154: PetscErrorCode NEPSetFromOptions_SLP(NEP nep)
155: {
156: NEP_SLP *ctx = (NEP_SLP*)nep->data;
159: ctx->setfromoptionscalled = PETSC_TRUE;
160: return(0);
161: }
165: static PetscErrorCode NEPSLPSetEPS_SLP(NEP nep,EPS eps)
166: {
168: NEP_SLP *ctx = (NEP_SLP*)nep->data;
171: PetscObjectReference((PetscObject)eps);
172: EPSDestroy(&ctx->eps);
173: ctx->eps = eps;
174: PetscLogObjectParent(nep,ctx->eps);
175: nep->setupcalled = 0;
176: return(0);
177: }
181: /*@
182: NEPSLPSetEPS - Associate a linear eigensolver object (EPS) to the
183: nonlinear eigenvalue solver.
185: Collective on NEP
187: Input Parameters:
188: + nep - nonlinear eigenvalue solver
189: - eps - the eigensolver object
191: Level: advanced
193: .seealso: NEPSLPGetEPS()
194: @*/
195: PetscErrorCode NEPSLPSetEPS(NEP nep,EPS eps)
196: {
203: PetscTryMethod(nep,"NEPSLPSetEPS_C",(NEP,EPS),(nep,eps));
204: return(0);
205: }
209: static PetscErrorCode NEPSLPGetEPS_SLP(NEP nep,EPS *eps)
210: {
212: NEP_SLP *ctx = (NEP_SLP*)nep->data;
215: if (!ctx->eps) {
216: EPSCreate(PetscObjectComm((PetscObject)nep),&ctx->eps);
217: EPSSetOptionsPrefix(ctx->eps,((PetscObject)nep)->prefix);
218: EPSAppendOptionsPrefix(ctx->eps,"nep_");
219: STSetOptionsPrefix(ctx->eps->st,((PetscObject)ctx->eps)->prefix);
220: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)nep,1);
221: PetscLogObjectParent(nep,ctx->eps);
222: if (!nep->ip) { NEPGetIP(nep,&nep->ip); }
223: EPSSetIP(ctx->eps,nep->ip);
224: }
225: *eps = ctx->eps;
226: return(0);
227: }
231: /*@
232: NEPSLPGetEPS - Retrieve the linear eigensolver object (EPS) associated
233: to the nonlinear eigenvalue solver.
235: Not Collective
237: Input Parameter:
238: . nep - nonlinear eigenvalue solver
240: Output Parameter:
241: . eps - the eigensolver object
243: Level: advanced
245: .seealso: NEPSLPSetEPS()
246: @*/
247: PetscErrorCode NEPSLPGetEPS(NEP nep,EPS *eps)
248: {
254: PetscTryMethod(nep,"NEPSLPGetEPS_C",(NEP,EPS*),(nep,eps));
255: return(0);
256: }
260: PetscErrorCode NEPView_SLP(NEP nep,PetscViewer viewer)
261: {
263: NEP_SLP *ctx = (NEP_SLP*)nep->data;
266: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
267: PetscViewerASCIIPushTab(viewer);
268: EPSView(ctx->eps,viewer);
269: PetscViewerASCIIPopTab(viewer);
270: return(0);
271: }
275: PetscErrorCode NEPReset_SLP(NEP nep)
276: {
278: NEP_SLP *ctx = (NEP_SLP*)nep->data;
281: if (!ctx->eps) { EPSReset(ctx->eps); }
282: NEPReset_Default(nep);
283: return(0);
284: }
288: PetscErrorCode NEPDestroy_SLP(NEP nep)
289: {
291: NEP_SLP *ctx = (NEP_SLP*)nep->data;
294: EPSDestroy(&ctx->eps);
295: PetscFree(nep->data);
296: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NULL);
297: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NULL);
298: return(0);
299: }
303: PETSC_EXTERN PetscErrorCode NEPCreate_SLP(NEP nep)
304: {
306: NEP_SLP *ctx;
309: PetscNewLog(nep,NEP_SLP,&ctx);
310: nep->data = (void*)ctx;
311: nep->ops->solve = NEPSolve_SLP;
312: nep->ops->setup = NEPSetUp_SLP;
313: nep->ops->setfromoptions = NEPSetFromOptions_SLP;
314: nep->ops->reset = NEPReset_SLP;
315: nep->ops->destroy = NEPDestroy_SLP;
316: nep->ops->view = NEPView_SLP;
317: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NEPSLPSetEPS_SLP);
318: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NEPSLPGetEPS_SLP);
319: return(0);
320: }