Actual source code: interpol.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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>
 27: #include <slepc/private/pepimpl.h>

 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==PETSC_DEFAULT) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 51:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
 52:   if (nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target magnitude eigenvalues");
 53:   NEPCheckUnsupported(nep,NEP_FEATURE_CALLBACK | NEP_FEATURE_STOPPING | NEP_FEATURE_TWOSIDED);

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

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

 91:   NEPAllocateSolution(nep,0);
 92:   return(0);
 93: }

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

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

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

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

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

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:   BVSetActiveColumns(nep->V,0,nep->nconv);
177:   BVCreateVec(nep->V,&vr);
178: #if !defined(PETSC_USE_COMPLEX)
179:   VecDuplicate(vr,&vi);
180: #endif
181:   s = 2.0/(b-a);
182:   for (i=0;i<nep->nconv;i++) {
183:     PEPGetEigenpair(ctx->pep,i,&nep->eigr[i],&nep->eigi[i],vr,vi);
184:     nep->eigr[i] /= s;
185:     nep->eigr[i] += (a+b)/2.0;
186:     nep->eigi[i] /= s;
187:     BVInsertVec(nep->V,i,vr);
188: #if !defined(PETSC_USE_COMPLEX)
189:     if (nep->eigi[i]!=0.0) {
190:       BVInsertVec(nep->V,++i,vi);
191:     }
192: #endif
193:   }
194:   VecDestroy(&vr);
195:   VecDestroy(&vi);

197:   nep->state = NEP_STATE_EIGENVECTORS;
198:   return(0);
199: }

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

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

229: PetscErrorCode NEPSetFromOptions_Interpol(PetscOptionItems *PetscOptionsObject,NEP nep)
230: {
232:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
233:   PetscInt       i;
234:   PetscBool      flg1,flg2;
235:   PetscReal      r;

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

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

246:   PetscOptionsTail();

248:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
249:   PEPSetFromOptions(ctx->pep);
250:   return(0);
251: }

253: static PetscErrorCode NEPInterpolSetInterpolation_Interpol(NEP nep,PetscReal tol,PetscInt degree)
254: {
256:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

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

281: /*@
282:    NEPInterpolSetInterpolation - Sets the tolerance and maximum degree when building
283:    the interpolation polynomial.

285:    Collective on nep

287:    Input Parameters:
288: +  nep - nonlinear eigenvalue solver
289: .  tol - tolerance to stop computing polynomial coefficients
290: -  deg - maximum degree of interpolation

292:    Options Database Key:
293: +  -nep_interpol_interpolation_tol <tol> - Sets the tolerance to stop computing polynomial coefficients
294: -  -nep_interpol_interpolation_degree <degree> - Sets the maximum degree of interpolation

296:    Notes:
297:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

299:    Level: advanced

301: .seealso: NEPInterpolGetInterpolation()
302: @*/
303: PetscErrorCode NEPInterpolSetInterpolation(NEP nep,PetscReal tol,PetscInt deg)
304: {

311:   PetscTryMethod(nep,"NEPInterpolSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,deg));
312:   return(0);
313: }

315: static PetscErrorCode NEPInterpolGetInterpolation_Interpol(NEP nep,PetscReal *tol,PetscInt *deg)
316: {
317:   NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;

320:   if (tol) *tol = ctx->tol;
321:   if (deg) *deg = ctx->maxdeg;
322:   return(0);
323: }

325: /*@
326:    NEPInterpolGetInterpolation - Gets the tolerance and maximum degree when building
327:    the interpolation polynomial.

329:    Not Collective

331:    Input Parameter:
332: .  nep - nonlinear eigenvalue solver

334:    Output Parameter:
335: +  tol - tolerance to stop computing polynomial coefficients
336: -  deg - maximum degree of interpolation

338:    Level: advanced

340: .seealso: NEPInterpolSetInterpolation()
341: @*/
342: PetscErrorCode NEPInterpolGetInterpolation(NEP nep,PetscReal *tol,PetscInt *deg)
343: {

348:   PetscUseMethod(nep,"NEPInterpolGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,deg));
349:   return(0);
350: }

352: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
353: {
355:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

358:   PetscObjectReference((PetscObject)pep);
359:   PEPDestroy(&ctx->pep);
360:   ctx->pep = pep;
361:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
362:   nep->state = NEP_STATE_INITIAL;
363:   return(0);
364: }

366: /*@
367:    NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
368:    nonlinear eigenvalue solver.

370:    Collective on nep

372:    Input Parameters:
373: +  nep - nonlinear eigenvalue solver
374: -  pep - the polynomial eigensolver object

376:    Level: advanced

378: .seealso: NEPInterpolGetPEP()
379: @*/
380: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
381: {

388:   PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
389:   return(0);
390: }

392: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
393: {
395:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

398:   if (!ctx->pep) {
399:     PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
400:     PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
401:     PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
402:     PEPAppendOptionsPrefix(ctx->pep,"nep_interpol_");
403:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
404:     PetscObjectSetOptions((PetscObject)ctx->pep,((PetscObject)nep)->options);
405:     PEPMonitorSet(ctx->pep,PEPMonitor_Interpol,nep,NULL);
406:   }
407:   *pep = ctx->pep;
408:   return(0);
409: }

411: /*@
412:    NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
413:    associated with the nonlinear eigenvalue solver.

415:    Not Collective

417:    Input Parameter:
418: .  nep - nonlinear eigenvalue solver

420:    Output Parameter:
421: .  pep - the polynomial eigensolver object

423:    Level: advanced

425: .seealso: NEPInterpolSetPEP()
426: @*/
427: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
428: {

434:   PetscUseMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
435:   return(0);
436: }

438: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
439: {
441:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
442:   PetscBool      isascii;

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

457: PetscErrorCode NEPReset_Interpol(NEP nep)
458: {
460:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

463:   PEPReset(ctx->pep);
464:   return(0);
465: }

467: PetscErrorCode NEPDestroy_Interpol(NEP nep)
468: {
470:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

473:   PEPDestroy(&ctx->pep);
474:   PetscFree(nep->data);
475:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NULL);
476:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NULL);
477:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
478:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
479:   return(0);
480: }

482: SLEPC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
483: {
485:   NEP_INTERPOL   *ctx;

488:   PetscNewLog(nep,&ctx);
489:   nep->data   = (void*)ctx;
490:   ctx->maxdeg = 5;
491:   ctx->tol    = PETSC_DEFAULT;

493:   nep->ops->solve          = NEPSolve_Interpol;
494:   nep->ops->setup          = NEPSetUp_Interpol;
495:   nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
496:   nep->ops->reset          = NEPReset_Interpol;
497:   nep->ops->destroy        = NEPDestroy_Interpol;
498:   nep->ops->view           = NEPView_Interpol;

500:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NEPInterpolSetInterpolation_Interpol);
501:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NEPInterpolGetInterpolation_Interpol);
502:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
503:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
504:   return(0);
505: }