Actual source code: interpol.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  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: "interpol"

 13:    Method: Polynomial interpolation

 15:    Algorithm:

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

 20:    References:

 22:        [1] C. Effenberger and D. Kresser, "Chebyshev interpolation for
 23:            nonlinear eigenvalue problems", BIT 52:933-951, 2012.
 24: */

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

 29: typedef struct {
 30:   PEP       pep;
 31:   PetscReal tol;       /* tolerance for norm of polynomial coefficients */
 32:   PetscInt  maxdeg;    /* maximum degree of interpolation polynomial */
 33:   PetscInt  deg;       /* actual degree of interpolation polynomial */
 34: } NEP_INTERPOL;

 36: PetscErrorCode NEPSetUp_Interpol(NEP nep)
 37: {
 39:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
 40:   ST             st;
 41:   RG             rg;
 42:   PetscReal      a,b,c,d,s,tol;
 43:   PetscScalar    zero=0.0;
 44:   PetscBool      flg,istrivial,trackall;
 45:   PetscInt       its,in;

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

 54:   /* transfer PEP options */
 55:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
 56:   PEPSetBasis(ctx->pep,PEP_BASIS_CHEBYSHEV1);
 57:   PEPSetWhichEigenpairs(ctx->pep,PEP_TARGET_MAGNITUDE);
 58:   PEPGetST(ctx->pep,&st);
 59:   STSetType(st,STSINVERT);
 60:   PEPSetDimensions(ctx->pep,nep->nev,nep->ncv?nep->ncv:PETSC_DEFAULT,nep->mpd?nep->mpd:PETSC_DEFAULT);
 61:   PEPGetTolerances(ctx->pep,&tol,&its);
 62:   if (tol==PETSC_DEFAULT) tol = (nep->tol==PETSC_DEFAULT)?SLEPC_DEFAULT_TOL:nep->tol;
 63:   if (ctx->tol==PETSC_DEFAULT) ctx->tol = tol;
 64:   if (!its) its = nep->max_it?nep->max_it:PETSC_DEFAULT;
 65:   PEPSetTolerances(ctx->pep,tol,its);
 66:   NEPGetTrackAll(nep,&trackall);
 67:   PEPSetTrackAll(ctx->pep,trackall);

 69:   /* transfer region options */
 70:   RGIsTrivial(nep->rg,&istrivial);
 71:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL requires a nontrivial region");
 72:   PetscObjectTypeCompare((PetscObject)nep->rg,RGINTERVAL,&flg);
 73:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for interval regions");
 74:   RGIntervalGetEndpoints(nep->rg,&a,&b,&c,&d);
 75:   if (a<=-PETSC_MAX_REAL || b>=PETSC_MAX_REAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for bounded intervals");
 76:   PEPGetRG(ctx->pep,&rg);
 77:   RGSetType(rg,RGINTERVAL);
 78:   if (a==b) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for intervals on the real axis");
 79:   s = 2.0/(b-a);
 80:   c = c*s;
 81:   d = d*s;
 82:   RGIntervalSetEndpoints(rg,-1.0,1.0,c,d);
 83:   RGCheckInside(nep->rg,1,&nep->target,&zero,&in);
 84:   if (in<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The target is not inside the target set");
 85:   PEPSetTarget(ctx->pep,(nep->target-(a+b)/2)*s);

 87:   NEPAllocateSolution(nep,0);
 88:   return(0);
 89: }

 91: /*
 92:   Input:
 93:     d, number of nodes to compute
 94:     a,b, interval extrems
 95:   Output:
 96:     *x, array containing the d Chebyshev nodes of the interval [a,b]
 97:     *dct2, coefficients to compute a discrete cosine transformation (DCT-II)
 98: */
 99: static PetscErrorCode ChebyshevNodes(PetscInt d,PetscReal a,PetscReal b,PetscScalar *x,PetscReal *dct2)
100: {
101:   PetscInt  j,i;
102:   PetscReal t;

105:   for (j=0;j<d+1;j++) {
106:     t = ((2*j+1)*PETSC_PI)/(2*(d+1));
107:     x[j] = (a+b)/2.0+((b-a)/2.0)*PetscCosReal(t);
108:     for (i=0;i<d+1;i++) dct2[j*(d+1)+i] = PetscCosReal(i*t);
109:   }
110:   return(0);
111: }

113: PetscErrorCode NEPSolve_Interpol(NEP nep)
114: {
116:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
117:   Mat            *A;
118:   PetscScalar    *x,*fx,t;
119:   PetscReal      *cs,a,b,s,aprox,aprox0=1.0,*matnorm;
120:   PetscInt       i,j,k,deg=ctx->maxdeg;
121:   PetscBool      hasmnorm;
122:   Vec            vr,vi=NULL;

125:   PetscMalloc5(deg+1,&A,(deg+1)*(deg+1),&cs,deg+1,&x,(deg+1)*nep->nt,&fx,nep->nt,&matnorm);
126:   for  (j=0;j<nep->nt;j++) {
127:     MatHasOperation(nep->A[j],MATOP_NORM,&hasmnorm);
128:     if (!hasmnorm) break;
129:     MatNorm(nep->A[j],NORM_INFINITY,matnorm+j);
130:   }
131:   if (!hasmnorm) for (j=0;j<nep->nt;j++) matnorm[j] = 1.0;
132:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
133:   ChebyshevNodes(deg,a,b,x,cs);
134:   for (j=0;j<nep->nt;j++) {
135:     for (i=0;i<=deg;i++) {
136:       FNEvaluateFunction(nep->f[j],x[i],&fx[i+j*(deg+1)]);
137:     }
138:   }
139:   /* Polynomial coefficients */
140:   ctx->deg = deg;
141:   for (k=0;k<=deg;k++) {
142:     MatDuplicate(nep->A[0],MAT_COPY_VALUES,&A[k]);
143:     t = 0.0;
144:     for (i=0;i<deg+1;i++) t += fx[i]*cs[i*(deg+1)+k];
145:     t *= 2.0/(deg+1);
146:     if (k==0) t /= 2.0;
147:     aprox = matnorm[0]*PetscAbsScalar(t);
148:     MatScale(A[k],t);
149:     for (j=1;j<nep->nt;j++) {
150:       t = 0.0;
151:       for (i=0;i<deg+1;i++) t += fx[i+j*(deg+1)]*cs[i*(deg+1)+k];
152:       t *= 2.0/(deg+1);
153:       if (k==0) t /= 2.0;
154:       aprox += matnorm[j]*PetscAbsScalar(t);
155:       MatAXPY(A[k],t,nep->A[j],nep->mstr);
156:     }
157:     if (k==0) aprox0 = aprox;
158:     if (aprox/aprox0<ctx->tol) { ctx->deg = k; deg = k; break; }
159:   }

161:   PEPSetOperators(ctx->pep,deg+1,A);
162:   for (k=0;k<=deg;k++) {
163:     MatDestroy(&A[k]);
164:   }
165:   PetscFree5(A,cs,x,fx,matnorm);

167:   /* Solve polynomial eigenproblem */
168:   PEPSolve(ctx->pep);
169:   PEPGetConverged(ctx->pep,&nep->nconv);
170:   PEPGetIterationNumber(ctx->pep,&nep->its);
171:   PEPGetConvergedReason(ctx->pep,(PEPConvergedReason*)&nep->reason);
172:   BVSetActiveColumns(nep->V,0,nep->nconv);
173:   BVCreateVec(nep->V,&vr);
174: #if !defined(PETSC_USE_COMPLEX)   
175:   VecDuplicate(vr,&vi);
176: #endif
177:   s = 2.0/(b-a);
178:   for (i=0;i<nep->nconv;i++) {
179:     PEPGetEigenpair(ctx->pep,i,&nep->eigr[i],&nep->eigi[i],vr,vi);
180:     nep->eigr[i] /= s;
181:     nep->eigr[i] += (a+b)/2.0;
182:     nep->eigi[i] /= s;
183:     BVInsertVec(nep->V,i,vr);
184: #if !defined(PETSC_USE_COMPLEX)   
185:     if (nep->eigi[i]!=0.0) {
186:       BVInsertVec(nep->V,++i,vi);
187:     }
188: #endif
189:   }
190:   VecDestroy(&vr);
191:   VecDestroy(&vi);

193:   nep->state = NEP_STATE_EIGENVECTORS;
194:   return(0);
195: }

197: static PetscErrorCode PEPMonitor_Interpol(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
198: {
199:   PetscInt       i,n;
200:   NEP            nep = (NEP)ctx;
201:   PetscReal      a,b,s;
202:   ST             st;

206:   n = PetscMin(nest,nep->ncv);
207:   for (i=0;i<n;i++) {
208:     nep->eigr[i]   = eigr[i];
209:     nep->eigi[i]   = eigi[i];
210:     nep->errest[i] = errest[i];
211:   }
212:   PEPGetST(pep,&st);
213:   STBackTransform(st,n,nep->eigr,nep->eigi);
214:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
215:   s = 2.0/(b-a);
216:   for (i=0;i<n;i++) {
217:     nep->eigr[i] /= s;
218:     nep->eigr[i] += (a+b)/2.0;
219:     nep->eigi[i] /= s;
220:   }
221:   NEPMonitor(nep,its,nconv,nep->eigr,nep->eigi,nep->errest,nest);
222:   return(0);
223: }

225: PetscErrorCode NEPSetFromOptions_Interpol(PetscOptionItems *PetscOptionsObject,NEP nep)
226: {
228:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
229:   PetscInt       i;
230:   PetscBool      flg1,flg2;
231:   PetscReal      r;

234:   PetscOptionsHead(PetscOptionsObject,"NEP Interpol Options");

236:     NEPInterpolGetInterpolation(nep,&r,&i);
237:     if (!i) i = PETSC_DEFAULT;
238:     PetscOptionsInt("-nep_interpol_interpolation_degree","Maximum degree of polynomial interpolation","NEPInterpolSetInterpolation",i,&i,&flg1);
239:     PetscOptionsReal("-nep_interpol_interpolation_tol","Tolerance for interpolation coefficients","NEPInterpolSetInterpolation",r,&r,&flg2);
240:     if (flg1 || flg2) { NEPInterpolSetInterpolation(nep,r,i); }

242:   PetscOptionsTail();

244:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
245:   PEPSetFromOptions(ctx->pep);
246:   return(0);
247: }

249: static PetscErrorCode NEPInterpolSetInterpolation_Interpol(NEP nep,PetscReal tol,PetscInt degree)
250: {
252:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

255:   if (tol == PETSC_DEFAULT) {
256:     ctx->tol   = PETSC_DEFAULT;
257:     nep->state = NEP_STATE_INITIAL;
258:   } else {
259:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
260:     ctx->tol = tol;
261:   }
262:   if (degree == PETSC_DEFAULT || degree == PETSC_DECIDE) {
263:     ctx->maxdeg = 0;
264:     if (nep->state) { NEPReset(nep); }
265:     nep->state = NEP_STATE_INITIAL;
266:   } else {
267:     if (degree <= 0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of degree. Must be > 0");
268:     if (ctx->maxdeg != degree) {
269:       ctx->maxdeg = degree;
270:       if (nep->state) { NEPReset(nep); }
271:       nep->state = NEP_STATE_INITIAL;
272:     }
273:   }
274:   return(0);
275: }

277: /*@
278:    NEPInterpolSetInterpolation - Sets the tolerance and maximum degree when building
279:    the interpolation polynomial.

281:    Collective on NEP

283:    Input Parameters:
284: +  nep - nonlinear eigenvalue solver
285: .  tol - tolerance to stop computing polynomial coefficients
286: -  deg - maximum degree of interpolation

288:    Options Database Key:
289: +  -nep_interpol_interpolation_tol <tol> - Sets the tolerance to stop computing polynomial coefficients
290: -  -nep_interpol_interpolation_degree <degree> - Sets the maximum degree of interpolation

292:    Notes:
293:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

295:    Level: advanced

297: .seealso: NEPInterpolGetInterpolation()
298: @*/
299: PetscErrorCode NEPInterpolSetInterpolation(NEP nep,PetscReal tol,PetscInt deg)
300: {

307:   PetscTryMethod(nep,"NEPInterpolSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,deg));
308:   return(0);
309: }

311: static PetscErrorCode NEPInterpolGetInterpolation_Interpol(NEP nep,PetscReal *tol,PetscInt *deg)
312: {
313:   NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;

316:   if (tol) *tol = ctx->tol;
317:   if (deg) *deg = ctx->maxdeg;
318:   return(0);
319: }

321: /*@
322:    NEPInterpolGetInterpolation - Gets the tolerance and maximum degree when building
323:    the interpolation polynomial.

325:    Not Collective

327:    Input Parameter:
328: .  nep - nonlinear eigenvalue solver

330:    Output Parameter:
331: +  tol - tolerance to stop computing polynomial coefficients
332: -  deg - maximum degree of interpolation

334:    Level: advanced

336: .seealso: NEPInterpolSetInterpolation()
337: @*/
338: PetscErrorCode NEPInterpolGetInterpolation(NEP nep,PetscReal *tol,PetscInt *deg)
339: {

344:   PetscUseMethod(nep,"NEPInterpolGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,deg));
345:   return(0);
346: }

348: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
349: {
351:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

354:   PetscObjectReference((PetscObject)pep);
355:   PEPDestroy(&ctx->pep);
356:   ctx->pep = pep;
357:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
358:   nep->state = NEP_STATE_INITIAL;
359:   return(0);
360: }

362: /*@
363:    NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
364:    nonlinear eigenvalue solver.

366:    Collective on NEP

368:    Input Parameters:
369: +  nep - nonlinear eigenvalue solver
370: -  pep - the polynomial eigensolver object

372:    Level: advanced

374: .seealso: NEPInterpolGetPEP()
375: @*/
376: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
377: {

384:   PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
385:   return(0);
386: }

388: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
389: {
391:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

394:   if (!ctx->pep) {
395:     PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
396:     PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
397:     PEPAppendOptionsPrefix(ctx->pep,"nep_interpol_");
398:     PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
399:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
400:     PEPMonitorSet(ctx->pep,PEPMonitor_Interpol,nep,NULL);
401:   }
402:   *pep = ctx->pep;
403:   return(0);
404: }

406: /*@
407:    NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
408:    associated with the nonlinear eigenvalue solver.

410:    Not Collective

412:    Input Parameter:
413: .  nep - nonlinear eigenvalue solver

415:    Output Parameter:
416: .  pep - the polynomial eigensolver object

418:    Level: advanced

420: .seealso: NEPInterpolSetPEP()
421: @*/
422: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
423: {

429:   PetscUseMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
430:   return(0);
431: }

433: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
434: {
436:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
437:   PetscBool      isascii;

440:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
441:   if (isascii) {
442:     if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
443:     PetscViewerASCIIPrintf(viewer,"  polynomial degree %D, max=%D\n",ctx->deg,ctx->maxdeg);
444:     PetscViewerASCIIPrintf(viewer,"  tolerance for norm of polynomial coefficients %g\n",ctx->tol);
445:     PetscViewerASCIIPushTab(viewer);
446:     PEPView(ctx->pep,viewer);
447:     PetscViewerASCIIPopTab(viewer);
448:   }
449:   return(0);
450: }

452: PetscErrorCode NEPReset_Interpol(NEP nep)
453: {
455:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

458:   PEPReset(ctx->pep);
459:   return(0);
460: }

462: PetscErrorCode NEPDestroy_Interpol(NEP nep)
463: {
465:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

468:   PEPDestroy(&ctx->pep);
469:   PetscFree(nep->data);
470:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NULL);
471:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NULL);
472:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
473:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
474:   return(0);
475: }

477: PETSC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
478: {
480:   NEP_INTERPOL   *ctx;

483:   PetscNewLog(nep,&ctx);
484:   nep->data   = (void*)ctx;
485:   ctx->maxdeg = 5;
486:   ctx->tol    = PETSC_DEFAULT;

488:   nep->ops->solve          = NEPSolve_Interpol;
489:   nep->ops->setup          = NEPSetUp_Interpol;
490:   nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
491:   nep->ops->reset          = NEPReset_Interpol;
492:   nep->ops->destroy        = NEPDestroy_Interpol;
493:   nep->ops->view           = NEPView_Interpol;

495:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NEPInterpolSetInterpolation_Interpol);
496:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NEPInterpolGetInterpolation_Interpol);
497:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
498:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
499:   return(0);
500: }