Actual source code: interpol.c

slepc-3.18.3 2023-03-24
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>

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

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

 45:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
 47:   if (nep->max_it==PETSC_DEFAULT) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 48:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
 50:   NEPCheckUnsupported(nep,NEP_FEATURE_CALLBACK | NEP_FEATURE_STOPPING | NEP_FEATURE_TWOSIDED);

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

 70:   /* transfer region options */
 71:   RGIsTrivial(nep->rg,&istrivial);
 73:   PetscObjectTypeCompare((PetscObject)nep->rg,RGINTERVAL,&flg);
 75:   RGIntervalGetEndpoints(nep->rg,&a,&b,&c,&d);
 77:   PEPGetRG(ctx->pep,&rg);
 78:   RGSetType(rg,RGINTERVAL);
 80:   s = 2.0/(b-a);
 81:   c = c*s;
 82:   d = d*s;
 83:   RGIntervalSetEndpoints(rg,-1.0,1.0,c,d);
 84:   RGCheckInside(nep->rg,1,&nep->target,&zero,&in);
 86:   PEPSetTarget(ctx->pep,(nep->target-(a+b)/2)*s);

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

 92: /*
 93:   Input:
 94:     d, number of nodes to compute
 95:     a,b, interval extremes
 96:   Output:
 97:     *x, array containing the d Chebyshev nodes of the interval [a,b]
 98:     *dct2, coefficients to compute a discrete cosine transformation (DCT-II)
 99: */
100: static PetscErrorCode ChebyshevNodes(PetscInt d,PetscReal a,PetscReal b,PetscScalar *x,PetscReal *dct2)
101: {
102:   PetscInt  j,i;
103:   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: {
115:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
116:   Mat            *A,*P;
117:   PetscScalar    *x,*fx,t;
118:   PetscReal      *cs,a,b,s,aprox,aprox0=1.0,*matnorm;
119:   PetscInt       i,j,k,deg=ctx->maxdeg;
120:   PetscBool      hasmnorm=PETSC_FALSE;
121:   Vec            vr,vi=NULL;
122:   ST             st;

124:   PetscMalloc4((deg+1)*(deg+1),&cs,deg+1,&x,(deg+1)*nep->nt,&fx,nep->nt,&matnorm);
125:   for  (j=0;j<nep->nt;j++) {
126:     MatHasOperation(nep->A[j],MATOP_NORM,&hasmnorm);
127:     if (!hasmnorm) break;
128:     MatNorm(nep->A[j],NORM_INFINITY,matnorm+j);
129:   }
130:   if (!hasmnorm) for (j=0;j<nep->nt;j++) matnorm[j] = 1.0;
131:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
132:   ChebyshevNodes(deg,a,b,x,cs);
133:   for (j=0;j<nep->nt;j++) {
134:     for (i=0;i<=deg;i++) FNEvaluateFunction(nep->f[j],x[i],&fx[i+j*(deg+1)]);
135:   }
136:   /* Polynomial coefficients */
137:   PetscMalloc1(deg+1,&A);
138:   if (nep->P) PetscMalloc1(deg+1,&P);
139:   ctx->deg = deg;
140:   for (k=0;k<=deg;k++) {
141:     MatDuplicate(nep->A[0],MAT_COPY_VALUES,&A[k]);
142:     if (nep->P) MatDuplicate(nep->P[0],MAT_COPY_VALUES,&P[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:     if (nep->P) MatScale(P[k],t);
150:     for (j=1;j<nep->nt;j++) {
151:       t = 0.0;
152:       for (i=0;i<deg+1;i++) t += fx[i+j*(deg+1)]*cs[i*(deg+1)+k];
153:       t *= 2.0/(deg+1);
154:       if (k==0) t /= 2.0;
155:       aprox += matnorm[j]*PetscAbsScalar(t);
156:       MatAXPY(A[k],t,nep->A[j],nep->mstr);
157:       if (nep->P) MatAXPY(P[k],t,nep->P[j],nep->mstrp);
158:     }
159:     if (k==0) aprox0 = aprox;
160:     if (k>1 && aprox/aprox0<ctx->tol) { ctx->deg = k; deg = k; break; }
161:   }
162:   PEPSetOperators(ctx->pep,deg+1,A);
163:   MatDestroyMatrices(deg+1,&A);
164:   if (nep->P) {
165:     PEPGetST(ctx->pep,&st);
166:     STSetSplitPreconditioner(st,deg+1,P,nep->mstrp);
167:     MatDestroyMatrices(deg+1,&P);
168:   }
169:   PetscFree4(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) BVInsertVec(nep->V,++i,vi);
190: #endif
191:   }
192:   VecDestroy(&vr);
193:   VecDestroy(&vi);

195:   nep->state = NEP_STATE_EIGENVECTORS;
196:   return 0;
197: }

199: static PetscErrorCode PEPMonitor_Interpol(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
200: {
201:   PetscInt       i,n;
202:   NEP            nep = (NEP)ctx;
203:   PetscReal      a,b,s;
204:   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(NEP nep,PetscOptionItems *PetscOptionsObject)
226: {
227:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
228:   PetscInt       i;
229:   PetscBool      flg1,flg2;
230:   PetscReal      r;

232:   PetscOptionsHeadBegin(PetscOptionsObject,"NEP Interpol Options");

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

240:   PetscOptionsHeadEnd();

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

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

251:   if (tol == PETSC_DEFAULT) {
252:     ctx->tol   = PETSC_DEFAULT;
253:     nep->state = NEP_STATE_INITIAL;
254:   } else {
256:     ctx->tol = tol;
257:   }
258:   if (degree == PETSC_DEFAULT || degree == PETSC_DECIDE) {
259:     ctx->maxdeg = 0;
260:     if (nep->state) NEPReset(nep);
261:     nep->state = NEP_STATE_INITIAL;
262:   } else {
264:     if (ctx->maxdeg != degree) {
265:       ctx->maxdeg = degree;
266:       if (nep->state) NEPReset(nep);
267:       nep->state = NEP_STATE_INITIAL;
268:     }
269:   }
270:   return 0;
271: }

273: /*@
274:    NEPInterpolSetInterpolation - Sets the tolerance and maximum degree when building
275:    the interpolation polynomial.

277:    Collective on nep

279:    Input Parameters:
280: +  nep - nonlinear eigenvalue solver
281: .  tol - tolerance to stop computing polynomial coefficients
282: -  deg - maximum degree of interpolation

284:    Options Database Key:
285: +  -nep_interpol_interpolation_tol <tol> - Sets the tolerance to stop computing polynomial coefficients
286: -  -nep_interpol_interpolation_degree <degree> - Sets the maximum degree of interpolation

288:    Notes:
289:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

291:    Level: advanced

293: .seealso: NEPInterpolGetInterpolation()
294: @*/
295: PetscErrorCode NEPInterpolSetInterpolation(NEP nep,PetscReal tol,PetscInt deg)
296: {
300:   PetscTryMethod(nep,"NEPInterpolSetInterpolation_C",(NEP,PetscReal,PetscInt),(nep,tol,deg));
301:   return 0;
302: }

304: static PetscErrorCode NEPInterpolGetInterpolation_Interpol(NEP nep,PetscReal *tol,PetscInt *deg)
305: {
306:   NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;

308:   if (tol) *tol = ctx->tol;
309:   if (deg) *deg = ctx->maxdeg;
310:   return 0;
311: }

313: /*@
314:    NEPInterpolGetInterpolation - Gets the tolerance and maximum degree when building
315:    the interpolation polynomial.

317:    Not Collective

319:    Input Parameter:
320: .  nep - nonlinear eigenvalue solver

322:    Output Parameters:
323: +  tol - tolerance to stop computing polynomial coefficients
324: -  deg - maximum degree of interpolation

326:    Level: advanced

328: .seealso: NEPInterpolSetInterpolation()
329: @*/
330: PetscErrorCode NEPInterpolGetInterpolation(NEP nep,PetscReal *tol,PetscInt *deg)
331: {
333:   PetscUseMethod(nep,"NEPInterpolGetInterpolation_C",(NEP,PetscReal*,PetscInt*),(nep,tol,deg));
334:   return 0;
335: }

337: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
338: {
339:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

341:   PetscObjectReference((PetscObject)pep);
342:   PEPDestroy(&ctx->pep);
343:   ctx->pep = pep;
344:   nep->state = NEP_STATE_INITIAL;
345:   return 0;
346: }

348: /*@
349:    NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
350:    nonlinear eigenvalue solver.

352:    Collective on nep

354:    Input Parameters:
355: +  nep - nonlinear eigenvalue solver
356: -  pep - the polynomial eigensolver object

358:    Level: advanced

360: .seealso: NEPInterpolGetPEP()
361: @*/
362: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
363: {
367:   PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
368:   return 0;
369: }

371: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
372: {
373:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

375:   if (!ctx->pep) {
376:     PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
377:     PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
378:     PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
379:     PEPAppendOptionsPrefix(ctx->pep,"nep_interpol_");
380:     PetscObjectSetOptions((PetscObject)ctx->pep,((PetscObject)nep)->options);
381:     PEPMonitorSet(ctx->pep,PEPMonitor_Interpol,nep,NULL);
382:   }
383:   *pep = ctx->pep;
384:   return 0;
385: }

387: /*@
388:    NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
389:    associated with the nonlinear eigenvalue solver.

391:    Not Collective

393:    Input Parameter:
394: .  nep - nonlinear eigenvalue solver

396:    Output Parameter:
397: .  pep - the polynomial eigensolver object

399:    Level: advanced

401: .seealso: NEPInterpolSetPEP()
402: @*/
403: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
404: {
407:   PetscUseMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
408:   return 0;
409: }

411: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
412: {
413:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
414:   PetscBool      isascii;

416:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
417:   if (isascii) {
418:     if (!ctx->pep) NEPInterpolGetPEP(nep,&ctx->pep);
419:     PetscViewerASCIIPrintf(viewer,"  polynomial degree %" PetscInt_FMT ", max=%" PetscInt_FMT "\n",ctx->deg,ctx->maxdeg);
420:     PetscViewerASCIIPrintf(viewer,"  tolerance for norm of polynomial coefficients %g\n",(double)ctx->tol);
421:     PetscViewerASCIIPushTab(viewer);
422:     PEPView(ctx->pep,viewer);
423:     PetscViewerASCIIPopTab(viewer);
424:   }
425:   return 0;
426: }

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

432:   PEPReset(ctx->pep);
433:   return 0;
434: }

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

440:   PEPDestroy(&ctx->pep);
441:   PetscFree(nep->data);
442:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NULL);
443:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NULL);
444:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
445:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
446:   return 0;
447: }

449: SLEPC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
450: {
451:   NEP_INTERPOL   *ctx;

453:   PetscNew(&ctx);
454:   nep->data   = (void*)ctx;
455:   ctx->maxdeg = 5;
456:   ctx->tol    = PETSC_DEFAULT;

458:   nep->ops->solve          = NEPSolve_Interpol;
459:   nep->ops->setup          = NEPSetUp_Interpol;
460:   nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
461:   nep->ops->reset          = NEPReset_Interpol;
462:   nep->ops->destroy        = NEPDestroy_Interpol;
463:   nep->ops->view           = NEPView_Interpol;

465:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetInterpolation_C",NEPInterpolSetInterpolation_Interpol);
466:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetInterpolation_C",NEPInterpolGetInterpolation_Interpol);
467:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
468:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
469:   return 0;
470: }