Actual source code: interpol.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*

  3:    SLEPc nonlinear eigensolver: "interpol"

  5:    Method: Polynomial interpolation

  7:    Algorithm:

  9:        Uses a PEP object to solve the interpolated NEP. Currently supports
 10:        only Chebyshev interpolation on an interval.

 12:    References:

 14:        [1] C. Effenberger and D. Kresser, "Chebyshev interpolation for
 15:            nonlinear eigenvalue problems", BIT 52:933-951, 2012.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 21:    This file is part of SLEPc.

 23:    SLEPc is free software: you can redistribute it and/or modify it under  the
 24:    terms of version 3 of the GNU Lesser General Public License as published by
 25:    the Free Software Foundation.

 27:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 28:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 29:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 30:    more details.

 32:    You  should have received a copy of the GNU Lesser General  Public  License
 33:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 34:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35: */

 37: #include <slepc/private/nepimpl.h>         /*I "slepcnep.h" I*/
 38: #include <slepc/private/pepimpl.h>         /*I "slepcpep.h" I*/

 40: typedef struct {
 41:   PEP       pep;
 42:   PetscInt  deg;
 43: } NEP_INTERPOL;

 47: PetscErrorCode NEPSetUp_Interpol(NEP nep)
 48: {
 50:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
 51:   ST             st;
 52:   RG             rg;
 53:   PetscReal      a,b,c,d,s,tol;
 54:   PetscScalar    zero=0.0;
 55:   PetscBool      flg,istrivial,trackall;
 56:   PetscInt       its,in;

 59:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
 60:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 61:   if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 62:   if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL only available for split operator");
 63:   if (nep->stopping!=NEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support user-defined stopping test");

 65:   /* transfer PEP options */
 66:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
 67:   PEPSetBV(ctx->pep,nep->V);
 68:   PEPSetBasis(ctx->pep,PEP_BASIS_CHEBYSHEV1);
 69:   PEPSetWhichEigenpairs(ctx->pep,PEP_TARGET_MAGNITUDE);
 70:   PEPGetST(ctx->pep,&st);
 71:   STSetType(st,STSINVERT);
 72:   PEPSetDimensions(ctx->pep,nep->nev,nep->ncv?nep->ncv:PETSC_DEFAULT,nep->mpd?nep->mpd:PETSC_DEFAULT);
 73:   tol=ctx->pep->tol;
 74:   if (tol==PETSC_DEFAULT) tol = (nep->tol==PETSC_DEFAULT)?SLEPC_DEFAULT_TOL/10.0:nep->tol/10.0;
 75:   its=ctx->pep->max_it;
 76:   if (!its) its = nep->max_it?nep->max_it:PETSC_DEFAULT;
 77:   PEPSetTolerances(ctx->pep,tol,its);
 78:   NEPGetTrackAll(nep,&trackall);
 79:   PEPSetTrackAll(ctx->pep,trackall);

 81:   /* transfer region options */
 82:   RGIsTrivial(nep->rg,&istrivial);
 83:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL requires a nontrivial region");
 84:   PetscObjectTypeCompare((PetscObject)nep->rg,RGINTERVAL,&flg);
 85:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for interval regions");
 86:   RGIntervalGetEndpoints(nep->rg,&a,&b,&c,&d);
 87:   if (a<=-PETSC_MAX_REAL || b>=PETSC_MAX_REAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for bounded intervals");
 88:   PEPGetRG(ctx->pep,&rg);
 89:   RGSetType(rg,RGINTERVAL);
 90:   if (a==b) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for intervals on the real axis");
 91:   s = 2.0/(b-a);
 92:   c = c*s;
 93:   d = d*s;
 94:   RGIntervalSetEndpoints(rg,-1.0,1.0,c,d);
 95:   RGCheckInside(nep->rg,1,&nep->target,&zero,&in);
 96:   if (in<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The target is not inside the target set");
 97:   PEPSetTarget(ctx->pep,(nep->target-(a+b)/2)*s);

 99:   NEPAllocateSolution(nep,0);
100:   return(0);
101: }

105: /*
106:   Input: 
107:     d, number of nodes to compute
108:     a,b, interval extrems
109:   Output:
110:     *x, array containing the d Chebyshev nodes of the interval [a,b]
111:     *dct2, coefficients to compute a discrete cosine transformation (DCT-II)
112: */
113: static PetscErrorCode ChebyshevNodes(PetscInt d,PetscReal a,PetscReal b,PetscScalar *x,PetscReal *dct2)
114: {
115:   PetscInt  j,i;
116:   PetscReal t;

119:   for (j=0;j<d+1;j++) {
120:     t = ((2*j+1)*PETSC_PI)/(2*(d+1));
121:     x[j] = (a+b)/2.0+((b-a)/2.0)*PetscCosReal(t);
122:     for (i=0;i<d+1;i++) dct2[j*(d+1)+i] = PetscCosReal(i*t);
123:   }
124:   return(0);
125: }

129: PetscErrorCode NEPSolve_Interpol(NEP nep)
130: {
132:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
133:   Mat            *A;   /*T=nep->function,Tp=nep->jacobian;*/
134:   PetscScalar    *x,*fx,t;
135:   PetscReal      *cs,a,b,s;
136:   PetscInt       i,j,k,deg=ctx->deg;

139:   PetscMalloc4(deg+1,&A,(deg+1)*(deg+1),&cs,deg+1,&x,(deg+1)*nep->nt,&fx);
140:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
141:   ChebyshevNodes(deg,a,b,x,cs);
142:   for (j=0;j<nep->nt;j++) {
143:     for (i=0;i<=deg;i++) {
144:       FNEvaluateFunction(nep->f[j],x[i],&fx[i+j*(deg+1)]);
145:     }
146:   }

148:   /* Polynomial coefficients */
149:   for (k=0;k<=deg;k++) {
150:     MatDuplicate(nep->A[0],MAT_COPY_VALUES,&A[k]);
151:     t = 0.0;
152:     for (i=0;i<deg+1;i++) t += fx[i]*cs[i*(deg+1)+k];
153:     t *= 2.0/(deg+1); 
154:     if (k==0) t /= 2.0;
155:     MatScale(A[k],t);
156:     for (j=1;j<nep->nt;j++) {
157:       t = 0.0;
158:       for (i=0;i<deg+1;i++) t += fx[i+j*(deg+1)]*cs[i*(deg+1)+k];
159:       t *= 2.0/(deg+1); 
160:       if (k==0) t /= 2.0;
161:       MatAXPY(A[k],t,nep->A[j],nep->mstr);
162:     }
163:   }

165:   PEPSetOperators(ctx->pep,deg+1,A);
166:   for (k=0;k<=deg;k++) {
167:     MatDestroy(&A[k]);
168:   }
169:   PetscFree4(A,cs,x,fx);

171:   /* Solve polynomial eigenproblem */
172:   PEPSolve(ctx->pep);
173:   PEPGetConverged(ctx->pep,&nep->nconv);
174:   PEPGetIterationNumber(ctx->pep,&nep->its);
175:   PEPGetConvergedReason(ctx->pep,(PEPConvergedReason*)&nep->reason);
176:   s = 2.0/(b-a);
177:   for (i=0;i<nep->nconv;i++) {
178:     PEPGetEigenpair(ctx->pep,i,&nep->eigr[i],&nep->eigi[i],NULL,NULL);
179:     nep->eigr[i] /= s;
180:     nep->eigr[i] += (a+b)/2.0;
181:     nep->eigi[i] /= s;
182:   }
183:   nep->state = NEP_STATE_EIGENVECTORS;
184:   return(0);
185: }

189: static PetscErrorCode PEPMonitor_Interpol(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
190: {
191:   PetscInt       i,n;
192:   NEP            nep = (NEP)ctx;
193:   PetscReal      a,b,s;

197:   n = PetscMin(nest,nep->ncv);
198:   for (i=0;i<n;i++) {
199:     nep->eigr[i]   = eigr[i];
200:     nep->eigi[i]   = eigi[i];
201:     nep->errest[i] = errest[i];
202:   }
203:   STBackTransform(pep->st,n,nep->eigr,nep->eigi);
204:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
205:   s = 2.0/(b-a);
206:   for (i=0;i<n;i++) {
207:     nep->eigr[i] /= s;
208:     nep->eigr[i] += (a+b)/2.0;
209:     nep->eigi[i] /= s;
210:   }  
211:   NEPMonitor(nep,its,nconv,nep->eigr,nep->eigi,nep->errest,nest);
212:   return(0);
213: }

217: PetscErrorCode NEPSetFromOptions_Interpol(PetscOptionItems *PetscOptionsObject,NEP nep)
218: {
220:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

223:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
224:   PEPSetFromOptions(ctx->pep);
225:   PetscOptionsHead(PetscOptionsObject,"NEP Interpol Options");
226:   PetscOptionsInt("-nep_interpol_degree","Degree of interpolation polynomial","NEPInterpolSetDegree",ctx->deg,&ctx->deg,NULL);
227:   PetscOptionsTail();
228:   return(0);
229: }

233: static PetscErrorCode NEPInterpolSetDegree_Interpol(NEP nep,PetscInt deg)
234: {
235:   NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;

238:   ctx->deg = deg;
239:   return(0);
240: }

244: /*@
245:    NEPInterpolSetDegree - Sets the degree of the interpolation polynomial.

247:    Collective on NEP

249:    Input Parameters:
250: +  nep - nonlinear eigenvalue solver
251: -  deg - polynomial degree

253:    Level: advanced

255: .seealso: NEPInterpolGetDegree()
256: @*/
257: PetscErrorCode NEPInterpolSetDegree(NEP nep,PetscInt deg)
258: {

264:   PetscTryMethod(nep,"NEPInterpolSetDegree_C",(NEP,PetscInt),(nep,deg));
265:   return(0);
266: }

270: static PetscErrorCode NEPInterpolGetDegree_Interpol(NEP nep,PetscInt *deg)
271: {
272:   NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;

275:   *deg = ctx->deg;
276:   return(0);
277: }

281: /*@
282:    NEPInterpolGetDegree - Gets the degree of the interpolation polynomial.

284:    Not Collective

286:    Input Parameter:
287: .  nep - nonlinear eigenvalue solver

289:    Output Parameter:
290: .  deg - the polynomial degree

292:    Level: advanced

294: .seealso: NEPInterpolSetDegree()
295: @*/
296: PetscErrorCode NEPInterpolGetDegree(NEP nep,PetscInt *deg)
297: {

303:   PetscUseMethod(nep,"NEPInterpolGetDegree_C",(NEP,PetscInt*),(nep,deg));
304:   return(0);
305: }

309: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
310: {
312:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

315:   PetscObjectReference((PetscObject)pep);
316:   PEPDestroy(&ctx->pep);
317:   ctx->pep = pep;
318:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
319:   nep->state = NEP_STATE_INITIAL;
320:   return(0);
321: }

325: /*@
326:    NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
327:    nonlinear eigenvalue solver.

329:    Collective on NEP

331:    Input Parameters:
332: +  nep - nonlinear eigenvalue solver
333: -  pep - the polynomial eigensolver object

335:    Level: advanced

337: .seealso: NEPInterpolGetPEP()
338: @*/
339: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
340: {

347:   PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
348:   return(0);
349: }

353: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
354: {
356:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
357:   ST             st;

360:   if (!ctx->pep) {
361:     PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
362:     PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
363:     PEPAppendOptionsPrefix(ctx->pep,"nep_interpol_");
364:     PEPGetST(ctx->pep,&st);
365:     STSetOptionsPrefix(st,((PetscObject)ctx->pep)->prefix);
366:     PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
367:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
368:     PEPMonitorSet(ctx->pep,PEPMonitor_Interpol,nep,NULL);
369:   }
370:   *pep = ctx->pep;
371:   return(0);
372: }

376: /*@
377:    NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
378:    associated with the nonlinear eigenvalue solver.

380:    Not Collective

382:    Input Parameter:
383: .  nep - nonlinear eigenvalue solver

385:    Output Parameter:
386: .  pep - the polynomial eigensolver object

388:    Level: advanced

390: .seealso: NEPInterpolSetPEP()
391: @*/
392: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
393: {

399:   PetscUseMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
400:   return(0);
401: }

405: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
406: {
408:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
409:   PetscBool      isascii;

412:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
413:   if (isascii) {
414:     if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
415:     PetscViewerASCIIPrintf(viewer,"  Interpol: polynomial degree %D\n",ctx->deg);
416:     PetscViewerASCIIPushTab(viewer);
417:     PEPView(ctx->pep,viewer);
418:     PetscViewerASCIIPopTab(viewer);
419:   }
420:   return(0);
421: }

425: PetscErrorCode NEPReset_Interpol(NEP nep)
426: {
428:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

431:   if (!ctx->pep) { PEPReset(ctx->pep); }
432:   return(0);
433: }

437: PetscErrorCode NEPDestroy_Interpol(NEP nep)
438: {
440:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

443:   PEPDestroy(&ctx->pep);
444:   PetscFree(nep->data);
445:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NULL);
446:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NULL);
447:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
448:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
449:   return(0);
450: }

454: PETSC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
455: {
457:   NEP_INTERPOL   *ctx;

460:   PetscNewLog(nep,&ctx);
461:   ctx->deg  = 5;
462:   nep->data = (void*)ctx;

464:   nep->ops->solve          = NEPSolve_Interpol;
465:   nep->ops->setup          = NEPSetUp_Interpol;
466:   nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
467:   nep->ops->reset          = NEPReset_Interpol;
468:   nep->ops->destroy        = NEPDestroy_Interpol;
469:   nep->ops->view           = NEPView_Interpol;
470:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NEPInterpolSetDegree_Interpol);
471:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NEPInterpolGetDegree_Interpol);
472:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
473:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
474:   return(0);
475: }