Actual source code: slp.c
slepc-3.8.0 2017-10-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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 nonlinear eigensolver: "slp"
13: Method: Succesive linear problems
15: Algorithm:
17: Newton-type iteration based on first order Taylor approximation.
19: References:
21: [1] A. Ruhe, "Algorithms for the nonlinear eigenvalue problem", SIAM J.
22: Numer. Anal. 10(4):674-689, 1973.
23: */
25: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
27: typedef struct {
28: EPS eps; /* linear eigensolver for T*z = mu*Tp*z */
29: } NEP_SLP;
31: PetscErrorCode NEPSetUp_SLP(NEP nep)
32: {
34: NEP_SLP *ctx = (NEP_SLP*)nep->data;
35: ST st;
36: PetscBool istrivial;
39: if (nep->nev>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Requested several eigenpairs but this solver can compute only one");
40: if (nep->ncv) { PetscInfo(nep,"Setting ncv = 1, ignoring user-provided value\n"); }
41: nep->ncv = 1;
42: if (nep->mpd) { PetscInfo(nep,"Setting mpd = 1, ignoring user-provided value\n"); }
43: nep->mpd = 1;
44: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
45: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
46: if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");
48: RGIsTrivial(nep->rg,&istrivial);
49: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");
51: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
52: EPSSetWhichEigenpairs(ctx->eps,EPS_TARGET_MAGNITUDE);
53: EPSSetTarget(ctx->eps,0.0);
54: EPSGetST(ctx->eps,&st);
55: STSetType(st,STSINVERT);
56: EPSSetTolerances(ctx->eps,nep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:nep->tol/10.0,nep->max_it?nep->max_it:PETSC_DEFAULT);
58: NEPAllocateSolution(nep,0);
59: NEPSetWorkVecs(nep,1);
60: return(0);
61: }
63: PetscErrorCode NEPSolve_SLP(NEP nep)
64: {
66: NEP_SLP *ctx = (NEP_SLP*)nep->data;
67: Mat T=nep->function,Tp=nep->jacobian;
68: Vec u,r=nep->work[0];
69: PetscScalar lambda,mu,im;
70: PetscReal resnorm;
71: PetscInt nconv;
74: /* get initial approximation of eigenvalue and eigenvector */
75: NEPGetDefaultShift(nep,&lambda);
76: if (!nep->nini) {
77: BVSetRandomColumn(nep->V,0);
78: }
79: BVGetColumn(nep->V,0,&u);
81: /* Restart loop */
82: while (nep->reason == NEP_CONVERGED_ITERATING) {
83: nep->its++;
85: /* evaluate T(lambda) and T'(lambda) */
86: NEPComputeFunction(nep,lambda,T,T);
87: NEPComputeJacobian(nep,lambda,Tp);
89: /* form residual, r = T(lambda)*u (used in convergence test only) */
90: MatMult(T,u,r);
92: /* convergence test */
93: VecNorm(r,NORM_2,&resnorm);
94: (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
95: nep->eigr[nep->nconv] = lambda;
96: if (nep->errest[nep->nconv]<=nep->tol) {
97: nep->nconv = nep->nconv + 1;
98: }
99: (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
100: NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,1);
102: if (nep->reason == NEP_CONVERGED_ITERATING) {
103: /* compute eigenvalue correction mu and eigenvector approximation u */
104: EPSSetOperators(ctx->eps,T,Tp);
105: EPSSetInitialSpace(ctx->eps,1,&u);
106: EPSSolve(ctx->eps);
107: EPSGetConverged(ctx->eps,&nconv);
108: if (!nconv) {
109: PetscInfo1(nep,"iter=%D, inner iteration failed, stopping solve\n",nep->its);
110: nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
111: break;
112: }
113: EPSGetEigenpair(ctx->eps,0,&mu,&im,u,NULL);
114: if (PetscAbsScalar(im)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Complex eigenvalue approximation - not implemented in real scalars");
116: /* correct eigenvalue */
117: lambda = lambda - mu;
118: }
119: }
120: BVRestoreColumn(nep->V,0,&u);
121: return(0);
122: }
124: PetscErrorCode NEPSetFromOptions_SLP(PetscOptionItems *PetscOptionsObject,NEP nep)
125: {
127: NEP_SLP *ctx = (NEP_SLP*)nep->data;
130: PetscOptionsHead(PetscOptionsObject,"NEP SLP Options");
131: PetscOptionsTail();
133: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
134: EPSSetFromOptions(ctx->eps);
135: return(0);
136: }
138: static PetscErrorCode NEPSLPSetEPS_SLP(NEP nep,EPS eps)
139: {
141: NEP_SLP *ctx = (NEP_SLP*)nep->data;
144: PetscObjectReference((PetscObject)eps);
145: EPSDestroy(&ctx->eps);
146: ctx->eps = eps;
147: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
148: nep->state = NEP_STATE_INITIAL;
149: return(0);
150: }
152: /*@
153: NEPSLPSetEPS - Associate a linear eigensolver object (EPS) to the
154: nonlinear eigenvalue solver.
156: Collective on NEP
158: Input Parameters:
159: + nep - nonlinear eigenvalue solver
160: - eps - the eigensolver object
162: Level: advanced
164: .seealso: NEPSLPGetEPS()
165: @*/
166: PetscErrorCode NEPSLPSetEPS(NEP nep,EPS eps)
167: {
174: PetscTryMethod(nep,"NEPSLPSetEPS_C",(NEP,EPS),(nep,eps));
175: return(0);
176: }
178: static PetscErrorCode NEPSLPGetEPS_SLP(NEP nep,EPS *eps)
179: {
181: NEP_SLP *ctx = (NEP_SLP*)nep->data;
184: if (!ctx->eps) {
185: EPSCreate(PetscObjectComm((PetscObject)nep),&ctx->eps);
186: EPSSetOptionsPrefix(ctx->eps,((PetscObject)nep)->prefix);
187: EPSAppendOptionsPrefix(ctx->eps,"nep_slp_");
188: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)nep,1);
189: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
190: }
191: *eps = ctx->eps;
192: return(0);
193: }
195: /*@
196: NEPSLPGetEPS - Retrieve the linear eigensolver object (EPS) associated
197: to the nonlinear eigenvalue solver.
199: Not Collective
201: Input Parameter:
202: . nep - nonlinear eigenvalue solver
204: Output Parameter:
205: . eps - the eigensolver object
207: Level: advanced
209: .seealso: NEPSLPSetEPS()
210: @*/
211: PetscErrorCode NEPSLPGetEPS(NEP nep,EPS *eps)
212: {
218: PetscUseMethod(nep,"NEPSLPGetEPS_C",(NEP,EPS*),(nep,eps));
219: return(0);
220: }
222: PetscErrorCode NEPView_SLP(NEP nep,PetscViewer viewer)
223: {
225: NEP_SLP *ctx = (NEP_SLP*)nep->data;
226: PetscBool isascii;
229: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
230: if (isascii) {
231: if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
232: PetscViewerASCIIPushTab(viewer);
233: EPSView(ctx->eps,viewer);
234: PetscViewerASCIIPopTab(viewer);
235: }
236: return(0);
237: }
239: PetscErrorCode NEPReset_SLP(NEP nep)
240: {
242: NEP_SLP *ctx = (NEP_SLP*)nep->data;
245: EPSReset(ctx->eps);
246: return(0);
247: }
249: PetscErrorCode NEPDestroy_SLP(NEP nep)
250: {
252: NEP_SLP *ctx = (NEP_SLP*)nep->data;
255: EPSDestroy(&ctx->eps);
256: PetscFree(nep->data);
257: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NULL);
258: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NULL);
259: return(0);
260: }
262: PETSC_EXTERN PetscErrorCode NEPCreate_SLP(NEP nep)
263: {
265: NEP_SLP *ctx;
268: PetscNewLog(nep,&ctx);
269: nep->data = (void*)ctx;
271: nep->ops->solve = NEPSolve_SLP;
272: nep->ops->setup = NEPSetUp_SLP;
273: nep->ops->setfromoptions = NEPSetFromOptions_SLP;
274: nep->ops->reset = NEPReset_SLP;
275: nep->ops->destroy = NEPDestroy_SLP;
276: nep->ops->view = NEPView_SLP;
278: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NEPSLPSetEPS_SLP);
279: PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NEPSLPGetEPS_SLP);
280: return(0);
281: }