Actual source code: interpol.c
slepc-3.5.2 2014-10-10
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-2014, 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*/
39: typedef struct {
40: PEP pep;
41: PetscInt deg;
42: } NEP_INTERPOL;
46: PetscErrorCode NEPSetUp_Interpol(NEP nep)
47: {
49: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
50: ST st;
51: RG rg;
52: PetscReal a,b,c,d,s;
53: PetscBool flg,istrivial;
56: if (nep->ncv) { /* ncv set */
57: if (nep->ncv<nep->nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
58: } else if (nep->mpd) { /* mpd set */
59: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
60: } else { /* neither set: defaults depend on nev being small or large */
61: if (nep->nev<500) nep->ncv = PetscMin(nep->n,PetscMax(2*nep->nev,nep->nev+15));
62: else {
63: nep->mpd = 500;
64: nep->ncv = PetscMin(nep->n,nep->nev+nep->mpd);
65: }
66: }
67: if (!nep->mpd) nep->mpd = nep->ncv;
68: if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
69: if (!nep->max_it) nep->max_it = PetscMax(5000,2*nep->n/nep->ncv);
70: if (!nep->max_funcs) nep->max_funcs = nep->max_it;
71: if (!nep->split) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL only available for split operator");
73: /* transfer PEP options */
74: if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
75: PEPSetBV(ctx->pep,nep->V);
76: PEPSetBasis(ctx->pep,PEP_BASIS_CHEBYSHEV1);
77: PEPSetWhichEigenpairs(ctx->pep,PEP_TARGET_MAGNITUDE);
78: PEPSetTarget(ctx->pep,0.0);
79: PEPGetST(ctx->pep,&st);
80: STSetType(st,STSINVERT);
81: PEPSetDimensions(ctx->pep,nep->nev,nep->ncv?nep->ncv:PETSC_DEFAULT,nep->mpd?nep->mpd:PETSC_DEFAULT);
82: PEPSetTolerances(ctx->pep,nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:nep->rtol/10.0,nep->max_it?nep->max_it:PETSC_DEFAULT);
84: /* transfer region options */
85: RGIsTrivial(nep->rg,&istrivial);
86: if (istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEPINTERPOL requires a nontrivial region");
87: PetscObjectTypeCompare((PetscObject)nep->rg,RGINTERVAL,&flg);
88: if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for interval regions");
89: RGIntervalGetEndpoints(nep->rg,&a,&b,&c,&d);
90: if (a<=-PETSC_MAX_REAL || b>=PETSC_MAX_REAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Only implemented for bounded intervals");
91: PEPGetRG(ctx->pep,&rg);
92: RGIsTrivial(rg,&istrivial);
93: if (istrivial) { /* user did not give region options */
94: RGSetType(rg,RGINTERVAL);
95: s = 2.0/(b-a);
96: c = (c==0)? -PETSC_MAX_REAL: c*s;
97: d = (d==0)? PETSC_MAX_REAL: d*s;
98: RGIntervalSetEndpoints(rg,-1.0,1.0,c,d);
99: }
101: NEPAllocateSolution(nep,0);
102: return(0);
103: }
107: /*
108: Input:
109: d, number of nodes to compute
110: a,b, interval extrems
111: Output:
112: *x, array containing the d Chebyshev nodes of the interval [a,b]
113: *dct2, coefficients to compute a discrete cosine transformation (DCT-II)
114: */
115: static PetscErrorCode ChebyshevNodes(PetscInt d,PetscReal a,PetscReal b,PetscScalar *x,PetscReal *dct2)
116: {
117: PetscInt j,i;
118: PetscReal t;
121: for (j=0;j<d+1;j++) {
122: t = ((2*j+1)*PETSC_PI)/(2*(d+1));
123: x[j] = (a+b)/2.0+((b-a)/2.0)*PetscCosReal(t);
124: for (i=0;i<d+1;i++) dct2[j*(d+1)+i] = PetscCosReal(i*t);
125: }
126: return(0);
127: }
131: PetscErrorCode NEPSolve_Interpol(NEP nep)
132: {
134: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
135: Mat *A; /*T=nep->function,Tp=nep->jacobian;*/
136: PetscScalar *x,*fx,t;
137: PetscReal *cs,a,b,s;
138: PetscInt i,j,k,deg=ctx->deg;
141: PetscMalloc4(deg+1,&A,(deg+1)*(deg+1),&cs,deg+1,&x,(deg+1)*nep->nt,&fx);
142: RGIntervalGetEndpoints(nep->rg,&a,&b,NULL,NULL);
143: ChebyshevNodes(deg,a,b,x,cs);
144: for (j=0;j<nep->nt;j++) {
145: for (i=0;i<=deg;i++) {
146: FNEvaluateFunction(nep->f[j],x[i],&fx[i+j*(deg+1)]);
147: }
148: }
150: /* Polynomial coefficients */
151: for (k=0;k<=deg;k++) {
152: MatDuplicate(nep->A[0],MAT_COPY_VALUES,&A[k]);
153: t = 0.0;
154: for (i=0;i<deg+1;i++) t += fx[i]*cs[i*(deg+1)+k];
155: t *= 2.0/(deg+1);
156: if (k==0) t /= 2.0;
157: MatScale(A[k],t);
158: for (j=1;j<nep->nt;j++) {
159: t = 0.0;
160: for (i=0;i<deg+1;i++) t += fx[i+j*(deg+1)]*cs[i*(deg+1)+k];
161: t *= 2.0/(deg+1);
162: if (k==0) t /= 2.0;
163: MatAXPY(A[k],t,nep->A[j],SUBSET_NONZERO_PATTERN);
164: }
165: }
167: PEPSetOperators(ctx->pep,deg+1,A);
168: for (k=0;k<=deg;k++) {
169: MatDestroy(&A[k]);
170: }
171: PetscFree4(A,cs,x,fx);
173: /* Solve polynomial eigenproblem */
174: PEPSolve(ctx->pep);
175: PEPGetConverged(ctx->pep,&nep->nconv);
176: PEPGetIterationNumber(ctx->pep,&nep->its);
177: PEPGetConvergedReason(ctx->pep,(PEPConvergedReason*)&nep->reason);
178: s = 2.0/(b-a);
179: for (i=0;i<nep->nconv;i++) {
180: PEPGetEigenpair(ctx->pep,i,&nep->eigr[i],&nep->eigi[i],NULL,NULL);
181: nep->eigr[i] /= s;
182: nep->eigr[i] += (a+b)/2.0;
183: nep->eigi[i] /= s;
184: }
185: nep->state = NEP_STATE_EIGENVECTORS;
186: return(0);
187: }
191: PetscErrorCode NEPSetFromOptions_Interpol(NEP nep)
192: {
194: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
197: if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
198: PEPSetFromOptions(ctx->pep);
199: PetscOptionsHead("NEP Interpol Options");
200: PetscOptionsInt("-nep_interpol_degree","Degree of interpolation polynomial","NEPInterpolSetDegree",ctx->deg,&ctx->deg,NULL);
201: return(0);
202: }
206: static PetscErrorCode NEPInterpolSetDegree_Interpol(NEP nep,PetscInt deg)
207: {
208: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
211: ctx->deg = deg;
212: return(0);
213: }
217: /*@
218: NEPInterpolSetDegree - Sets the degree of the interpolation polynomial.
220: Collective on NEP
222: Input Parameters:
223: + nep - nonlinear eigenvalue solver
224: - deg - polynomial degree
226: Level: advanced
228: .seealso: NEPInterpolGetDegree()
229: @*/
230: PetscErrorCode NEPInterpolSetDegree(NEP nep,PetscInt deg)
231: {
237: PetscTryMethod(nep,"NEPInterpolSetDegree_C",(NEP,PetscInt),(nep,deg));
238: return(0);
239: }
243: static PetscErrorCode NEPInterpolGetDegree_Interpol(NEP nep,PetscInt *deg)
244: {
245: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
248: *deg = ctx->deg;
249: return(0);
250: }
254: /*@
255: NEPInterpolGetDegree - Gets the degree of the interpolation polynomial.
257: Not Collective
259: Input Parameter:
260: . nep - nonlinear eigenvalue solver
262: Output Parameter:
263: . pep - the polynomial degree
265: Level: advanced
267: .seealso: NEPInterpolSetDegree()
268: @*/
269: PetscErrorCode NEPInterpolGetDegree(NEP nep,PetscInt *deg)
270: {
276: PetscTryMethod(nep,"NEPInterpolGetDegree_C",(NEP,PetscInt*),(nep,deg));
277: return(0);
278: }
282: static PetscErrorCode NEPInterpolSetPEP_Interpol(NEP nep,PEP pep)
283: {
285: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
288: PetscObjectReference((PetscObject)pep);
289: PEPDestroy(&ctx->pep);
290: ctx->pep = pep;
291: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
292: nep->state = NEP_STATE_INITIAL;
293: return(0);
294: }
298: /*@
299: NEPInterpolSetPEP - Associate a polynomial eigensolver object (PEP) to the
300: nonlinear eigenvalue solver.
302: Collective on NEP
304: Input Parameters:
305: + nep - nonlinear eigenvalue solver
306: - pep - the polynomial eigensolver object
308: Level: advanced
310: .seealso: NEPInterpolGetPEP()
311: @*/
312: PetscErrorCode NEPInterpolSetPEP(NEP nep,PEP pep)
313: {
320: PetscTryMethod(nep,"NEPInterpolSetPEP_C",(NEP,PEP),(nep,pep));
321: return(0);
322: }
326: static PetscErrorCode NEPInterpolGetPEP_Interpol(NEP nep,PEP *pep)
327: {
329: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
330: ST st;
333: if (!ctx->pep) {
334: PEPCreate(PetscObjectComm((PetscObject)nep),&ctx->pep);
335: PEPSetOptionsPrefix(ctx->pep,((PetscObject)nep)->prefix);
336: PEPAppendOptionsPrefix(ctx->pep,"nep_");
337: PEPGetST(ctx->pep,&st);
338: STSetOptionsPrefix(st,((PetscObject)ctx->pep)->prefix);
339: PetscObjectIncrementTabLevel((PetscObject)ctx->pep,(PetscObject)nep,1);
340: PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->pep);
341: }
342: *pep = ctx->pep;
343: return(0);
344: }
348: /*@
349: NEPInterpolGetPEP - Retrieve the polynomial eigensolver object (PEP)
350: associated with the nonlinear eigenvalue solver.
352: Not Collective
354: Input Parameter:
355: . nep - nonlinear eigenvalue solver
357: Output Parameter:
358: . pep - the polynomial eigensolver object
360: Level: advanced
362: .seealso: NEPInterpolSetPEP()
363: @*/
364: PetscErrorCode NEPInterpolGetPEP(NEP nep,PEP *pep)
365: {
371: PetscTryMethod(nep,"NEPInterpolGetPEP_C",(NEP,PEP*),(nep,pep));
372: return(0);
373: }
377: PetscErrorCode NEPView_Interpol(NEP nep,PetscViewer viewer)
378: {
380: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
383: if (!ctx->pep) { NEPInterpolGetPEP(nep,&ctx->pep); }
384: PetscViewerASCIIPrintf(viewer," Interpol: polynomial degree %D\n",ctx->deg);
385: PetscViewerASCIIPushTab(viewer);
386: PEPView(ctx->pep,viewer);
387: PetscViewerASCIIPopTab(viewer);
388: return(0);
389: }
393: PetscErrorCode NEPReset_Interpol(NEP nep)
394: {
396: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
399: if (!ctx->pep) { PEPReset(ctx->pep); }
400: return(0);
401: }
405: PetscErrorCode NEPDestroy_Interpol(NEP nep)
406: {
408: NEP_INTERPOL *ctx = (NEP_INTERPOL*)nep->data;
411: PEPDestroy(&ctx->pep);
412: PetscFree(nep->data);
413: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NULL);
414: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NULL);
415: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NULL);
416: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NULL);
417: return(0);
418: }
422: PETSC_EXTERN PetscErrorCode NEPCreate_Interpol(NEP nep)
423: {
425: NEP_INTERPOL *ctx;
428: PetscNewLog(nep,&ctx);
429: ctx->deg = 5;
430: nep->data = (void*)ctx;
432: nep->ops->solve = NEPSolve_Interpol;
433: nep->ops->setup = NEPSetUp_Interpol;
434: nep->ops->setfromoptions = NEPSetFromOptions_Interpol;
435: nep->ops->reset = NEPReset_Interpol;
436: nep->ops->destroy = NEPDestroy_Interpol;
437: nep->ops->view = NEPView_Interpol;
438: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetDegree_C",NEPInterpolSetDegree_Interpol);
439: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetDegree_C",NEPInterpolGetDegree_Interpol);
440: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolSetPEP_C",NEPInterpolSetPEP_Interpol);
441: PetscObjectComposeFunction((PetscObject)nep,"NEPInterpolGetPEP_C",NEPInterpolGetPEP_Interpol);
442: return(0);
443: }