Actual source code: slp.c
slepc-3.6.1 2015-09-03
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: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17: SLEPc - Scalable Library for Eigenvalue Problem Computations
18: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
20: This file is part of SLEPc.
22: SLEPc is free software: you can redistribute it and/or modify it under the
23: terms of version 3 of the GNU Lesser General Public License as published by
24: the Free Software Foundation.
26: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
27: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
28: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
29: more details.
31: You should have received a copy of the GNU Lesser General Public License
32: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: */
36: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
38: typedef struct {
39: EPS eps; /* linear eigensolver for T*z = mu*Tp*z */
40: } NEP_SLP;
44: PetscErrorCode NEPSetUp_SLP(NEP nep)
45: {
47: NEP_SLP *ctx = (NEP_SLP*)nep->data;
48: ST st;
49: PetscBool istrivial;
52: if (nep->ncv) { /* ncv set */
53: if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
54: } else if (nep->mpd) { /* mpd set */
55: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
56: } else { /* neither set: defaults depend on nev being small or large */
57: if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
58: else {
59: nep->mpd = 500;
60: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
61: }
62: }
63: if (!nep->mpd) nep->mpd = nep->ncv;
64: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
65: if (nep->nev>1) { PetscInfo(nep,"Warning: requested more than one eigenpair but SLP can only compute one\n"); }
66: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
67: if (!nep->max_funcs) nep->max_funcs = nep->max_it;
69: RGIsTrivial(nep->rg,&istrivial);
70: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");
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->ncv:PETSC_DEFAULT,nep->mpd?nep->mpd:PETSC_DEFAULT);
78: EPSSetTolerances(ctx->eps,nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:nep->rtol/10.0,nep->max_it?nep->max_it:PETSC_DEFAULT);
80: NEPAllocateSolution(nep,0);
81: NEPSetWorkVecs(nep,1);
82: return(0);
83: }
87: PetscErrorCode NEPSolve_SLP(NEP nep)
88: {
90: NEP_SLP *ctx = (NEP_SLP*)nep->data;
91: Mat T=nep->function,Tp=nep->jacobian;
92: Vec u,r=nep->work[0];
93: PetscScalar lambda,mu,im;
94: PetscReal relerr;
95: PetscInt nconv;
98: /* get initial approximation of eigenvalue and eigenvector */
99: NEPGetDefaultShift(nep,&lambda);
100: if (!nep->nini) {
101: BVSetRandomColumn(nep->V,0,nep->rand);
102: }
103: BVGetColumn(nep->V,0,&u);
105: /* Restart loop */
106: while (nep->reason == NEP_CONVERGED_ITERATING) {
107: nep->its++;
109: /* evaluate T(lambda) and T'(lambda) */
110: NEPComputeFunction(nep,lambda,T,T);
111: NEPComputeJacobian(nep,lambda,Tp);
113: /* form residual, r = T(lambda)*u (used in convergence test only) */
114: MatMult(T,u,r);
116: /* convergence test */
117: VecNorm(r,NORM_2,&relerr);
118: nep->errest[nep->nconv] = relerr;
119: nep->eigr[nep->nconv] = lambda;
120: if (relerr<=nep->rtol) {
121: nep->nconv = nep->nconv + 1;
122: nep->reason = NEP_CONVERGED_FNORM_RELATIVE;
123: }
124: NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->errest,1);
126: if (!nep->nconv) {
127: /* compute eigenvalue correction mu and eigenvector approximation u */
128: EPSSetOperators(ctx->eps,T,Tp);
129: EPSSetInitialSpace(ctx->eps,1,&u);
130: EPSSolve(ctx->eps);
131: EPSGetConverged(ctx->eps,&nconv);
132: if (!nconv) {
133: PetscInfo1(nep,"iter=%D, inner iteration failed, stopping solve\n",nep->its);
134: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
135: break;
136: }
137: EPSGetEigenpair(ctx->eps,0,&mu,&im,u,NULL);
138: if (PetscAbsScalar(im)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Complex eigenvalue approximation - not implemented in real scalars");
140: /* correct eigenvalue */
141: lambda = lambda - mu;
142: }
143: if (nep->its >= nep->max_it) nep->reason = NEP_DIVERGED_MAX_IT;
144: }
145: BVRestoreColumn(nep->V,0,&u);
146: return(0);
147: }
151: PetscErrorCode NEPSetFromOptions_SLP(PetscOptions *PetscOptionsObject,NEP nep)
152: {
154: NEP_SLP *ctx = (NEP_SLP*)nep->data;
157: PetscOptionsHead(PetscOptionsObject,"NEP SLP Options");
158: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
159: EPSSetFromOptions(ctx->eps);
160: PetscOptionsTail();
161: return(0);
162: }
166: static PetscErrorCode NEPSLPSetEPS_SLP(NEP nep,EPS eps)
167: {
169: NEP_SLP *ctx = (NEP_SLP*)nep->data;
172: PetscObjectReference((PetscObject)eps);
173: EPSDestroy(&ctx->eps);
174: ctx->eps = eps;
175: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
176: nep->state = NEP_STATE_INITIAL;
177: return(0);
178: }
182: /*@
183: NEPSLPSetEPS - Associate a linear eigensolver object (EPS) to the
184: nonlinear eigenvalue solver.
186: Collective on NEP
188: Input Parameters:
189: + nep - nonlinear eigenvalue solver
190: - eps - the eigensolver object
192: Level: advanced
194: .seealso: NEPSLPGetEPS()
195: @*/
196: PetscErrorCode NEPSLPSetEPS(NEP nep,EPS eps)
197: {
204: PetscTryMethod(nep,"NEPSLPSetEPS_C",(NEP,EPS),(nep,eps));
205: return(0);
206: }
210: static PetscErrorCode NEPSLPGetEPS_SLP(NEP nep,EPS *eps)
211: {
213: NEP_SLP *ctx = (NEP_SLP*)nep->data;
214: ST st;
217: if (!ctx->eps) {
218: EPSCreate(PetscObjectComm((PetscObject)nep),&ctx->eps);
219: EPSSetOptionsPrefix(ctx->eps,((PetscObject)nep)->prefix);
220: EPSAppendOptionsPrefix(ctx->eps,"nep_");
221: EPSGetST(ctx->eps,&st);
222: STSetOptionsPrefix(st,((PetscObject)ctx->eps)->prefix);
223: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)nep,1);
224: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
225: }
226: *eps = ctx->eps;
227: return(0);
228: }
232: /*@
233: NEPSLPGetEPS - Retrieve the linear eigensolver object (EPS) associated
234: to the nonlinear eigenvalue solver.
236: Not Collective
238: Input Parameter:
239: . nep - nonlinear eigenvalue solver
241: Output Parameter:
242: . eps - the eigensolver object
244: Level: advanced
246: .seealso: NEPSLPSetEPS()
247: @*/
248: PetscErrorCode NEPSLPGetEPS(NEP nep,EPS *eps)
249: {
255: PetscTryMethod(nep,"NEPSLPGetEPS_C",(NEP,EPS*),(nep,eps));
256: return(0);
257: }
261: PetscErrorCode NEPView_SLP(NEP nep,PetscViewer viewer)
262: {
264: NEP_SLP *ctx = (NEP_SLP*)nep->data;
267: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
268: PetscViewerASCIIPushTab(viewer);
269: EPSView(ctx->eps,viewer);
270: PetscViewerASCIIPopTab(viewer);
271: return(0);
272: }
276: PetscErrorCode NEPReset_SLP(NEP nep)
277: {
279: NEP_SLP *ctx = (NEP_SLP*)nep->data;
282: if (!ctx->eps) { EPSReset(ctx->eps); }
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,&ctx);
310: nep->data = (void*)ctx;
312: nep->ops->solve = NEPSolve_SLP;
313: nep->ops->setup = NEPSetUp_SLP;
314: nep->ops->setfromoptions = NEPSetFromOptions_SLP;
315: nep->ops->reset = NEPReset_SLP;
316: nep->ops->destroy = NEPDestroy_SLP;
317: nep->ops->view = NEPView_SLP;
318: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NEPSLPSetEPS_SLP);
319: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NEPSLPGetEPS_SLP);
320: return(0);
321: }