Actual source code: slp.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  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-2016, 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->nev>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Requested several eigenpairs but this solver can compute only one");
 53:   if (nep->ncv) { PetscInfo(nep,"Setting ncv = 1, ignoring user-provided value\n"); }
 54:   nep->ncv = 1;
 55:   if (nep->mpd) { PetscInfo(nep,"Setting mpd = 1, ignoring user-provided value\n"); }
 56:   nep->mpd = 1;
 57:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 58:   if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 59:   if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");

 61:   RGIsTrivial(nep->rg,&istrivial);
 62:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");

 64:   if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
 65:   EPSSetWhichEigenpairs(ctx->eps,EPS_TARGET_MAGNITUDE);
 66:   EPSSetTarget(ctx->eps,0.0);
 67:   EPSGetST(ctx->eps,&st);
 68:   STSetType(st,STSINVERT);
 69:   EPSSetTolerances(ctx->eps,nep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:nep->tol/10.0,nep->max_it?nep->max_it:PETSC_DEFAULT);

 71:   NEPAllocateSolution(nep,0);
 72:   NEPSetWorkVecs(nep,1);
 73:   return(0);
 74: }

 78: PetscErrorCode NEPSolve_SLP(NEP nep)
 79: {
 81:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;
 82:   Mat            T=nep->function,Tp=nep->jacobian;
 83:   Vec            u,r=nep->work[0];
 84:   PetscScalar    lambda,mu,im;
 85:   PetscReal      resnorm;
 86:   PetscInt       nconv;

 89:   /* get initial approximation of eigenvalue and eigenvector */
 90:   NEPGetDefaultShift(nep,&lambda);
 91:   if (!nep->nini) {
 92:     BVSetRandomColumn(nep->V,0);
 93:   }
 94:   BVGetColumn(nep->V,0,&u);

 96:   /* Restart loop */
 97:   while (nep->reason == NEP_CONVERGED_ITERATING) {
 98:     nep->its++;

100:     /* evaluate T(lambda) and T'(lambda) */
101:     NEPComputeFunction(nep,lambda,T,T);
102:     NEPComputeJacobian(nep,lambda,Tp);

104:     /* form residual,  r = T(lambda)*u (used in convergence test only) */
105:     MatMult(T,u,r);

107:     /* convergence test */
108:     VecNorm(r,NORM_2,&resnorm);
109:     (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
110:     nep->eigr[nep->nconv] = lambda;
111:     if (nep->errest[nep->nconv]<=nep->tol) {
112:       nep->nconv = nep->nconv + 1;
113:     }
114:     (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
115:     NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,1);

117:     if (nep->reason == NEP_CONVERGED_ITERATING) {
118:       /* compute eigenvalue correction mu and eigenvector approximation u */
119:       EPSSetOperators(ctx->eps,T,Tp);
120:       EPSSetInitialSpace(ctx->eps,1,&u);
121:       EPSSolve(ctx->eps);
122:       EPSGetConverged(ctx->eps,&nconv);
123:       if (!nconv) {
124:         PetscInfo1(nep,"iter=%D, inner iteration failed, stopping solve\n",nep->its);
125:         nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
126:         break;
127:       }
128:       EPSGetEigenpair(ctx->eps,0,&mu,&im,u,NULL);
129:       if (PetscAbsScalar(im)>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Complex eigenvalue approximation - not implemented in real scalars");

131:       /* correct eigenvalue */
132:       lambda = lambda - mu;
133:     }
134:   }
135:   BVRestoreColumn(nep->V,0,&u);
136:   return(0);
137: }

141: PetscErrorCode NEPSetFromOptions_SLP(PetscOptionItems *PetscOptionsObject,NEP nep)
142: {
144:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;

147:   PetscOptionsHead(PetscOptionsObject,"NEP SLP Options");
148:   if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
149:   EPSSetFromOptions(ctx->eps);
150:   PetscOptionsTail();
151:   return(0);
152: }

156: static PetscErrorCode NEPSLPSetEPS_SLP(NEP nep,EPS eps)
157: {
159:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;

162:   PetscObjectReference((PetscObject)eps);
163:   EPSDestroy(&ctx->eps);
164:   ctx->eps = eps;
165:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
166:   nep->state = NEP_STATE_INITIAL;
167:   return(0);
168: }

172: /*@
173:    NEPSLPSetEPS - Associate a linear eigensolver object (EPS) to the
174:    nonlinear eigenvalue solver.

176:    Collective on NEP

178:    Input Parameters:
179: +  nep - nonlinear eigenvalue solver
180: -  eps - the eigensolver object

182:    Level: advanced

184: .seealso: NEPSLPGetEPS()
185: @*/
186: PetscErrorCode NEPSLPSetEPS(NEP nep,EPS eps)
187: {

194:   PetscTryMethod(nep,"NEPSLPSetEPS_C",(NEP,EPS),(nep,eps));
195:   return(0);
196: }

200: static PetscErrorCode NEPSLPGetEPS_SLP(NEP nep,EPS *eps)
201: {
203:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;
204:   ST             st;

207:   if (!ctx->eps) {
208:     EPSCreate(PetscObjectComm((PetscObject)nep),&ctx->eps);
209:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)nep)->prefix);
210:     EPSAppendOptionsPrefix(ctx->eps,"nep_slp_");
211:     EPSGetST(ctx->eps,&st);
212:     STSetOptionsPrefix(st,((PetscObject)ctx->eps)->prefix);
213:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)nep,1);
214:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->eps);
215:   }
216:   *eps = ctx->eps;
217:   return(0);
218: }

222: /*@
223:    NEPSLPGetEPS - Retrieve the linear eigensolver object (EPS) associated
224:    to the nonlinear eigenvalue solver.

226:    Not Collective

228:    Input Parameter:
229: .  nep - nonlinear eigenvalue solver

231:    Output Parameter:
232: .  eps - the eigensolver object

234:    Level: advanced

236: .seealso: NEPSLPSetEPS()
237: @*/
238: PetscErrorCode NEPSLPGetEPS(NEP nep,EPS *eps)
239: {

245:   PetscUseMethod(nep,"NEPSLPGetEPS_C",(NEP,EPS*),(nep,eps));
246:   return(0);
247: }

251: PetscErrorCode NEPView_SLP(NEP nep,PetscViewer viewer)
252: {
254:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;
255:   PetscBool      isascii;

258:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
259:   if (isascii) {
260:     if (!ctx->eps) { NEPSLPGetEPS(nep,&ctx->eps); }
261:     PetscViewerASCIIPushTab(viewer);
262:     EPSView(ctx->eps,viewer);
263:     PetscViewerASCIIPopTab(viewer);
264:   }
265:   return(0);
266: }

270: PetscErrorCode NEPReset_SLP(NEP nep)
271: {
273:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;

276:   if (!ctx->eps) { EPSReset(ctx->eps); }
277:   return(0);
278: }

282: PetscErrorCode NEPDestroy_SLP(NEP nep)
283: {
285:   NEP_SLP        *ctx = (NEP_SLP*)nep->data;

288:   EPSDestroy(&ctx->eps);
289:   PetscFree(nep->data);
290:   PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NULL);
291:   PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NULL);
292:   return(0);
293: }

297: PETSC_EXTERN PetscErrorCode NEPCreate_SLP(NEP nep)
298: {
300:   NEP_SLP        *ctx;

303:   PetscNewLog(nep,&ctx);
304:   nep->data = (void*)ctx;

306:   nep->ops->solve          = NEPSolve_SLP;
307:   nep->ops->setup          = NEPSetUp_SLP;
308:   nep->ops->setfromoptions = NEPSetFromOptions_SLP;
309:   nep->ops->reset          = NEPReset_SLP;
310:   nep->ops->destroy        = NEPDestroy_SLP;
311:   nep->ops->view           = NEPView_SLP;
312:   PetscObjectComposeFunction((PetscObject)nep,"NEPSLPSetEPS_C",NEPSLPSetEPS_SLP);
313:   PetscObjectComposeFunction((PetscObject)nep,"NEPSLPGetEPS_C",NEPSLPGetEPS_SLP);
314:   return(0);
315: }