Actual source code: interpol.c

slepc-3.6.1 2015-09-03
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-2015, 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>

 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:   PetscBool      flg,istrivial;
 55:   PetscInt       its;

 58:   if (nep->ncv) { /* ncv set */
 59:     if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
 60:   } else if (nep->mpd) { /* mpd set */
 61:     nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
 62:   } else { /* neither set: defaults depend on nev being small or large */
 63:     if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
 64:     else {
 65:       nep->mpd = 500;
 66:       nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
 67:     }
 68:   }
 69:   if (!nep->mpd) nep->mpd = nep->ncv;
 70:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 71:   if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
 72:   if (!nep->max_funcs) nep->max_funcs = nep->max_it;
 73:   if (!nep->split) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL only available for split operator");

 75:   /* transfer PEP options */
 76:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
 77:   PEPSetBV(ctx->pep,nep->V);
 78:   PEPSetBasis(ctx->pep,PEP_BASIS_CHEBYSHEV1);
 79:   PEPSetWhichEigenpairs(ctx->pep,PEP_TARGET_MAGNITUDE);
 80:   PEPSetTarget(ctx->pep,0.0);
 81:   PEPGetST(ctx->pep,&st);
 82:   STSetType(st,STSINVERT);
 83:   PEPSetDimensions(ctx->pep,nep->nev,nep->ncv?nep->ncv:PETSC_DEFAULT,nep->mpd?nep->mpd:PETSC_DEFAULT);
 84:   tol=ctx->pep->tol;
 85:   if (tol==PETSC_DEFAULT) tol = (nep->rtol==PETSC_DEFAULT)?SLEPC_DEFAULT_TOL/10.0:nep->rtol/10.0;
 86:   its=ctx->pep->max_it;
 87:   if (!its) its = nep->max_it?nep->max_it:PETSC_DEFAULT;
 88:   PEPSetTolerances(ctx->pep,tol,its);

 90:   /* transfer region options */
 91:   RGIsTrivial(nep->rg,&istrivial);
 92:   if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL requires a nontrivial region");
 93:   PetscObjectTypeCompare((PetscObject)nep->rg,RGINTERVAL,&flg);
 94:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for interval regions");
 95:   RGIntervalGetEndpoints(nep->rg,&a,&b,&c,&d);
 96:   if (a<=-PETSC_MAX_REAL || b>=PETSC_MAX_REAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for bounded intervals");
 97:   PEPGetRG(ctx->pep,&rg);
 98:   RGIsTrivial(rg,&istrivial);
 99:   if (istrivial) {   /* user did not give region options */
100:     RGSetType(rg,RGINTERVAL);
101:     s = 2.0/(b-a);
102:     c = (c==0)? -PETSC_MAX_REAL: c*s;
103:     d = (d==0)? PETSC_MAX_REAL: d*s;
104:     RGIntervalSetEndpoints(rg,-1.0,1.0,c,d);
105:   }

107:   NEPAllocateSolution(nep,0);
108:   return(0);
109: }

113: /*
114:   Input: 
115:     d, number of nodes to compute
116:     a,b, interval extrems
117:   Output:
118:     *x, array containing the d Chebyshev nodes of the interval [a,b]
119:     *dct2, coefficients to compute a discrete cosine transformation (DCT-II)
120: */
121: static PetscErrorCode ChebyshevNodes(PetscInt d,PetscReal a,PetscReal b,PetscScalar *x,PetscReal *dct2)
122: {
123:   PetscInt  j,i;
124:   PetscReal t;

127:   for (j=0;j<d+1;j++) {
128:     t = ((2*j+1)*PETSC_PI)/(2*(d+1));
129:     x[j] = (a+b)/2.0+((b-a)/2.0)*PetscCosReal(t);
130:     for (i=0;i<d+1;i++) dct2[j*(d+1)+i] = PetscCosReal(i*t);
131:   }
132:   return(0);
133: }

137: PetscErrorCode NEPSolve_Interpol(NEP nep)
138: {
140:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
141:   Mat            *A;   /*T=nep->function,Tp=nep->jacobian;*/
142:   PetscScalar    *x,*fx,t;
143:   PetscReal      *cs,a,b,s;
144:   PetscInt       i,j,k,deg=ctx->deg;

147:   PetscMalloc4(deg+1,&A,(deg+1)*(deg+1),&cs,deg+1,&x,(deg+1)*nep->nt,&fx);
148:   RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
149:   ChebyshevNodes(deg,a,b,x,cs);
150:   for (j=0;j<nep->nt;j++) {
151:     for (i=0;i<=deg;i++) {
152:       FNEvaluateFunction(nep->f[j],x[i],&fx[i+j*(deg+1)]);
153:     }
154:   }

156:   /* Polynomial coefficients */
157:   for (k=0;k<=deg;k++) {
158:     MatDuplicate(nep->A[0],MAT_COPY_VALUES,&A[k]);
159:     t = 0.0;
160:     for (i=0;i<deg+1;i++) t += fx[i]*cs[i*(deg+1)+k];
161:     t *= 2.0/(deg+1); 
162:     if (k==0) t /= 2.0;
163:     MatScale(A[k],t);
164:     for (j=1;j<nep->nt;j++) {
165:       t = 0.0;
166:       for (i=0;i<deg+1;i++) t += fx[i+j*(deg+1)]*cs[i*(deg+1)+k];
167:       t *= 2.0/(deg+1); 
168:       if (k==0) t /= 2.0;
169:       MatAXPY(A[k],t,nep->A[j],SUBSET_NONZERO_PATTERN);
170:     }
171:   }

173:   PEPSetOperators(ctx->pep,deg+1,A);
174:   for (k=0;k<=deg;k++) {
175:     MatDestroy(&A[k]);
176:   }
177:   PetscFree4(A,cs,x,fx);

179:   /* Solve polynomial eigenproblem */
180:   PEPSolve(ctx->pep);
181:   PEPGetConverged(ctx->pep,&nep->nconv);
182:   PEPGetIterationNumber(ctx->pep,&nep->its);
183:   PEPGetConvergedReason(ctx->pep,(PEPConvergedReason*)&nep->reason);
184:   s = 2.0/(b-a);
185:   for (i=0;i<nep->nconv;i++) {
186:     PEPGetEigenpair(ctx->pep,i,&nep->eigr[i],&nep->eigi[i],NULL,NULL);
187:     nep->eigr[i] /= s;
188:     nep->eigr[i] += (a+b)/2.0;
189:     nep->eigi[i] /= s;
190:   }
191:   nep->state = NEP_STATE_EIGENVECTORS;
192:   return(0);
193: }

197: PetscErrorCode NEPSetFromOptions_Interpol(PetscOptions *PetscOptionsObject,NEP nep)
198: {
200:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

203:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
204:   PEPSetFromOptions(ctx->pep);
205:   PetscOptionsHead(PetscOptionsObject,"NEP Interpol Options");
206:   PetscOptionsInt("-nep_interpol_degree","Degree of interpolation polynomial","NEPInterpolSetDegree",ctx->deg,&ctx->deg,NULL);
207:   PetscOptionsTail();
208:   return(0);
209: }

213: static PetscErrorCode NEPInterpolSetDegree_Interpol(NEP nep,PetscInt deg)
214: {
215:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

218:   ctx->deg = deg;
219:   return(0);
220: }

224: /*@
225:    NEPInterpolSetDegree - Sets the degree of the interpolation polynomial.

227:    Collective on NEP

229:    Input Parameters:
230: +  nep - nonlinear eigenvalue solver
231: -  deg - polynomial degree

233:    Level: advanced

235: .seealso: NEPInterpolGetDegree()
236: @*/
237: PetscErrorCode NEPInterpolSetDegree(NEP nep,PetscInt deg)
238: {

244:   PetscTryMethod(nep,"NEPInterpolSetDegree_C",(NEP,PetscInt),(nep,deg));
245:   return(0);
246: }

250: static PetscErrorCode NEPInterpolGetDegree_Interpol(NEP nep,PetscInt *deg)
251: {
252:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

255:   *deg = ctx->deg;
256:   return(0);
257: }

261: /*@
262:    NEPInterpolGetDegree - Gets the degree of the interpolation polynomial.

264:    Not Collective

266:    Input Parameter:
267: .  nep - nonlinear eigenvalue solver

269:    Output Parameter:
270: .  pep - the polynomial degree

272:    Level: advanced

274: .seealso: NEPInterpolSetDegree()
275: @*/
276: PetscErrorCode NEPInterpolGetDegree(NEP nep,PetscInt *deg)
277: {

283:   PetscTryMethod(nep,"NEPInterpolGetDegree_C",(NEP,PetscInt*),(nep,deg));
284:   return(0);
285: }

289: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
290: {
292:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

295:   PetscObjectReference((PetscObject)pep);
296:   PEPDestroy(&ctx->pep);
297:   ctx->pep = pep;
298:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
299:   nep->state = NEP_STATE_INITIAL;
300:   return(0);
301: }

305: /*@
306:    NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
307:    nonlinear eigenvalue solver.

309:    Collective on NEP

311:    Input Parameters:
312: +  nep - nonlinear eigenvalue solver
313: -  pep - the polynomial eigensolver object

315:    Level: advanced

317: .seealso: NEPInterpolGetPEP()
318: @*/
319: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
320: {

327:   PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
328:   return(0);
329: }

333: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
334: {
336:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;
337:   ST             st;

340:   if (!ctx->pep) {
341:     PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
342:     PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
343:     PEPAppendOptionsPrefix(ctx->pep,"nep_");
344:     PEPGetST(ctx->pep,&st);
345:     STSetOptionsPrefix(st,((PetscObject)ctx->pep)->prefix);
346:     PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
347:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
348:   }
349:   *pep = ctx->pep;
350:   return(0);
351: }

355: /*@
356:    NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
357:    associated with the nonlinear eigenvalue solver.

359:    Not Collective

361:    Input Parameter:
362: .  nep - nonlinear eigenvalue solver

364:    Output Parameter:
365: .  pep - the polynomial eigensolver object

367:    Level: advanced

369: .seealso: NEPInterpolSetPEP()
370: @*/
371: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
372: {

378:   PetscTryMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
379:   return(0);
380: }

384: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
385: {
387:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

390:   if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
391:   PetscViewerASCIIPrintf(viewer,"  Interpol: polynomial degree %D\n",ctx->deg);
392:   PetscViewerASCIIPushTab(viewer);
393:   PEPView(ctx->pep,viewer);
394:   PetscViewerASCIIPopTab(viewer);
395:   return(0);
396: }

400: PetscErrorCode NEPReset_Interpol(NEP nep)
401: {
403:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

406:   if (!ctx->pep) { PEPReset(ctx->pep); }
407:   return(0);
408: }

412: PetscErrorCode NEPDestroy_Interpol(NEP nep)
413: {
415:   NEP_INTERPOL   *ctx = (NEP_INTERPOL*)nep->data;

418:   PEPDestroy(&ctx->pep);
419:   PetscFree(nep->data);
420:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NULL);
421:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NULL);
422:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
423:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
424:   return(0);
425: }

429: PETSC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
430: {
432:   NEP_INTERPOL   *ctx;

435:   PetscNewLog(nep,&ctx);
436:   ctx->deg  = 5;
437:   nep->data = (void*)ctx;

439:   nep->ops->solve          = NEPSolve_Interpol;
440:   nep->ops->setup          = NEPSetUp_Interpol;
441:   nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
442:   nep->ops->reset          = NEPReset_Interpol;
443:   nep->ops->destroy        = NEPDestroy_Interpol;
444:   nep->ops->view           = NEPView_Interpol;
445:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NEPInterpolSetDegree_Interpol);
446:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NEPInterpolGetDegree_Interpol);
447:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
448:   PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
449:   return(0);
450: }