Actual source code: filter.c
slepc-3.13.0 2020-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: Filter spectral transformation, to encapsulate polynomial filters
12: */
14: #include <slepc/private/stimpl.h> /*I "slepcst.h" I*/
15: #include "filter.h"
17: /*
18: Operator (filter):
19: Op P M
20: if nmat=1: p(A) NULL p(A)
21: */
22: PetscErrorCode STComputeOperator_Filter(ST st)
23: {
25: ST_FILTER *ctx = (ST_FILTER*)st->data;
28: if (st->nmat>1) SETERRQ(PetscObjectComm((PetscObject)st),1,"Only implemented for standard eigenvalue problem");
29: if (ctx->intb >= PETSC_MAX_REAL && ctx->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)st),1,"Must pass an interval with STFilterSetInterval()");
30: if (ctx->right == 0.0 && ctx->left == 0.0) SETERRQ(PetscObjectComm((PetscObject)st),1,"Must pass an approximate numerical range with STFilterSetRange()");
31: if (ctx->left > ctx->inta || ctx->right < ctx->intb) SETERRQ4(PetscObjectComm((PetscObject)st),1,"The requested interval [%g,%g] must be contained in the numerical range [%g,%g]",(double)ctx->inta,(double)ctx->intb,(double)ctx->left,(double)ctx->right);
32: if (!ctx->polyDegree) ctx->polyDegree = 100;
33: ctx->frame[0] = ctx->left;
34: ctx->frame[1] = ctx->inta;
35: ctx->frame[2] = ctx->intb;
36: ctx->frame[3] = ctx->right;
37: STFilter_FILTLAN_setFilter(st,&st->T[0]);
38: st->M = st->T[0];
39: MatDestroy(&st->P);
40: return(0);
41: }
43: PetscErrorCode STSetUp_Filter(ST st)
44: {
48: STSetWorkVecs(st,4);
49: return(0);
50: }
52: PetscErrorCode STSetFromOptions_Filter(PetscOptionItems *PetscOptionsObject,ST st)
53: {
55: PetscReal array[2]={0,0};
56: PetscInt k;
57: PetscBool flg;
60: PetscOptionsHead(PetscOptionsObject,"ST Filter Options");
62: k = 2;
63: PetscOptionsRealArray("-st_filter_interval","Interval containing the desired eigenvalues (two real values separated with a comma without spaces)","STFilterSetInterval",array,&k,&flg);
64: if (flg) {
65: if (k<2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_SIZ,"Must pass two values in -st_filter_interval (comma-separated without spaces)");
66: STFilterSetInterval(st,array[0],array[1]);
67: }
68: k = 2;
69: PetscOptionsRealArray("-st_filter_range","Interval containing all eigenvalues (two real values separated with a comma without spaces)","STFilterSetRange",array,&k,&flg);
70: if (flg) {
71: if (k<2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_SIZ,"Must pass two values in -st_filter_range (comma-separated without spaces)");
72: STFilterSetRange(st,array[0],array[1]);
73: }
74: PetscOptionsInt("-st_filter_degree","Degree of filter polynomial","STFilterSetDegree",100,&k,&flg);
75: if (flg) { STFilterSetDegree(st,k); }
77: PetscOptionsTail();
78: return(0);
79: }
81: static PetscErrorCode STFilterSetInterval_Filter(ST st,PetscReal inta,PetscReal intb)
82: {
83: ST_FILTER *ctx = (ST_FILTER*)st->data;
86: if (inta>intb) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be inta<intb");
87: if (ctx->inta != inta || ctx->intb != intb) {
88: ctx->inta = inta;
89: ctx->intb = intb;
90: st->state = ST_STATE_INITIAL;
91: st->opready = PETSC_FALSE;
92: }
93: return(0);
94: }
96: /*@
97: STFilterSetInterval - Defines the interval containing the desired eigenvalues.
99: Logically Collective on st
101: Input Parameters:
102: + st - the spectral transformation context
103: . inta - left end of the interval
104: - intb - right end of the interval
106: Options Database Key:
107: . -st_filter_interval <a,b> - set [a,b] as the interval of interest
109: Level: intermediate
111: Notes:
112: The filter will be configured to emphasize eigenvalues contained in the given
113: interval, and damp out eigenvalues outside it. If the interval is open, then
114: the filter is low- or high-pass, otherwise it is mid-pass.
116: Common usage is to set the interval in EPS with EPSSetInterval().
118: The interval must be contained within the numerical range of the matrix, see
119: STFilterSetRange().
121: .seealso: STFilterGetInterval(), STFilterSetRange(), EPSSetInterval()
122: @*/
123: PetscErrorCode STFilterSetInterval(ST st,PetscReal inta,PetscReal intb)
124: {
131: PetscTryMethod(st,"STFilterSetInterval_C",(ST,PetscReal,PetscReal),(st,inta,intb));
132: return(0);
133: }
135: static PetscErrorCode STFilterGetInterval_Filter(ST st,PetscReal *inta,PetscReal *intb)
136: {
137: ST_FILTER *ctx = (ST_FILTER*)st->data;
140: if (inta) *inta = ctx->inta;
141: if (intb) *intb = ctx->intb;
142: return(0);
143: }
145: /*@
146: STFilterGetInterval - Gets the interval containing the desired eigenvalues.
148: Not Collective
150: Input Parameter:
151: . st - the spectral transformation context
153: Output Parameter:
154: + inta - left end of the interval
155: - intb - right end of the interval
157: Level: intermediate
159: .seealso: STFilterSetInterval()
160: @*/
161: PetscErrorCode STFilterGetInterval(ST st,PetscReal *inta,PetscReal *intb)
162: {
167: PetscUseMethod(st,"STFilterGetInterval_C",(ST,PetscReal*,PetscReal*),(st,inta,intb));
168: return(0);
169: }
171: static PetscErrorCode STFilterSetRange_Filter(ST st,PetscReal left,PetscReal right)
172: {
173: ST_FILTER *ctx = (ST_FILTER*)st->data;
176: if (left>right) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONG,"Badly defined interval, must be left<right");
177: if (ctx->left != left || ctx->right != right) {
178: ctx->left = left;
179: ctx->right = right;
180: st->state = ST_STATE_INITIAL;
181: st->opready = PETSC_FALSE;
182: }
183: return(0);
184: }
186: /*@
187: STFilterSetRange - Defines the numerical range (or field of values) of the matrix, that is,
188: the interval containing all eigenvalues.
190: Logically Collective on st
192: Input Parameters:
193: + st - the spectral transformation context
194: . left - left end of the interval
195: - right - right end of the interval
197: Options Database Key:
198: . -st_filter_range <a,b> - set [a,b] as the numerical range
200: Level: intermediate
202: Notes:
203: The filter will be most effective if the numerical range is tight, that is, left and right
204: are good approximations to the leftmost and rightmost eigenvalues, respectively.
206: .seealso: STFilterGetRange(), STFilterSetInterval()
207: @*/
208: PetscErrorCode STFilterSetRange(ST st,PetscReal left,PetscReal right)
209: {
216: PetscTryMethod(st,"STFilterSetRange_C",(ST,PetscReal,PetscReal),(st,left,right));
217: return(0);
218: }
220: static PetscErrorCode STFilterGetRange_Filter(ST st,PetscReal *left,PetscReal *right)
221: {
222: ST_FILTER *ctx = (ST_FILTER*)st->data;
225: if (left) *left = ctx->left;
226: if (right) *right = ctx->right;
227: return(0);
228: }
230: /*@
231: STFilterGetRange - Gets the interval containing all eigenvalues.
233: Not Collective
235: Input Parameter:
236: . st - the spectral transformation context
238: Output Parameter:
239: + left - left end of the interval
240: - right - right end of the interval
242: Level: intermediate
244: .seealso: STFilterSetRange()
245: @*/
246: PetscErrorCode STFilterGetRange(ST st,PetscReal *left,PetscReal *right)
247: {
252: PetscUseMethod(st,"STFilterGetRange_C",(ST,PetscReal*,PetscReal*),(st,left,right));
253: return(0);
254: }
256: static PetscErrorCode STFilterSetDegree_Filter(ST st,PetscInt deg)
257: {
258: ST_FILTER *ctx = (ST_FILTER*)st->data;
261: if (deg == PETSC_DEFAULT || deg == PETSC_DECIDE) {
262: ctx->polyDegree = 0;
263: st->state = ST_STATE_INITIAL;
264: st->opready = PETSC_FALSE;
265: } else {
266: if (deg<=0) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of degree. Must be > 0");
267: if (ctx->polyDegree != deg) {
268: ctx->polyDegree = deg;
269: st->state = ST_STATE_INITIAL;
270: st->opready = PETSC_FALSE;
271: }
272: }
273: return(0);
274: }
276: /*@
277: STFilterSetDegree - Sets the degree of the filter polynomial.
279: Logically Collective on st
281: Input Parameters:
282: + st - the spectral transformation context
283: - deg - polynomial degree
285: Options Database Key:
286: . -st_filter_degree <deg> - sets the degree of the filter polynomial
288: Level: intermediate
290: .seealso: STFilterGetDegree()
291: @*/
292: PetscErrorCode STFilterSetDegree(ST st,PetscInt deg)
293: {
299: PetscTryMethod(st,"STFilterSetDegree_C",(ST,PetscInt),(st,deg));
300: return(0);
301: }
303: static PetscErrorCode STFilterGetDegree_Filter(ST st,PetscInt *deg)
304: {
305: ST_FILTER *ctx = (ST_FILTER*)st->data;
308: *deg = ctx->polyDegree;
309: return(0);
310: }
312: /*@
313: STFilterGetDegree - Gets the degree of the filter polynomial.
315: Not Collective
317: Input Parameter:
318: . st - the spectral transformation context
320: Output Parameter:
321: . deg - polynomial degree
323: Level: intermediate
325: .seealso: STFilterSetDegree()
326: @*/
327: PetscErrorCode STFilterGetDegree(ST st,PetscInt *deg)
328: {
334: PetscUseMethod(st,"STFilterGetDegree_C",(ST,PetscInt*),(st,deg));
335: return(0);
336: }
338: static PetscErrorCode STFilterGetThreshold_Filter(ST st,PetscReal *gamma)
339: {
340: ST_FILTER *ctx = (ST_FILTER*)st->data;
343: *gamma = ctx->filterInfo->yLimit;
344: return(0);
345: }
347: /*@
348: STFilterGetThreshold - Gets the threshold value gamma such that rho(lambda)>=gamma for
349: eigenvalues lambda inside the wanted interval and rho(lambda)<gamma for those outside.
351: Not Collective
353: Input Parameter:
354: . st - the spectral transformation context
356: Output Parameter:
357: . gamma - the threshold value
359: Level: developer
360: @*/
361: PetscErrorCode STFilterGetThreshold(ST st,PetscReal *gamma)
362: {
368: PetscUseMethod(st,"STFilterGetThreshold_C",(ST,PetscReal*),(st,gamma));
369: return(0);
370: }
372: PetscErrorCode STReset_Filter(ST st)
373: {
375: ST_FILTER *ctx = (ST_FILTER*)st->data;
378: ctx->left = 0.0;
379: ctx->right = 0.0;
380: MatDestroy(&ctx->T);
381: return(0);
382: }
384: PetscErrorCode STView_Filter(ST st,PetscViewer viewer)
385: {
387: ST_FILTER *ctx = (ST_FILTER*)st->data;
388: PetscBool isascii;
391: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
392: if (isascii) {
393: PetscViewerASCIIPrintf(viewer," Filter: interval of desired eigenvalues is [%g,%g]\n",(double)ctx->inta,(double)ctx->intb);
394: PetscViewerASCIIPrintf(viewer," Filter: numerical range is [%g,%g]\n",(double)ctx->left,(double)ctx->right);
395: PetscViewerASCIIPrintf(viewer," Filter: degree of filter polynomial is %D\n",ctx->polyDegree);
396: if (st->state>=ST_STATE_SETUP) {
397: PetscViewerASCIIPrintf(viewer," Filter: limit to accept eigenvalues: theta=%g\n",ctx->filterInfo->yLimit);
398: }
399: }
400: return(0);
401: }
403: PetscErrorCode STDestroy_Filter(ST st)
404: {
406: ST_FILTER *ctx = (ST_FILTER*)st->data;
409: PetscFree(ctx->opts);
410: PetscFree(ctx->filterInfo);
411: PetscFree(ctx->baseFilter);
412: PetscFree(st->data);
413: PetscObjectComposeFunction((PetscObject)st,"STFilterSetInterval_C",NULL);
414: PetscObjectComposeFunction((PetscObject)st,"STFilterGetInterval_C",NULL);
415: PetscObjectComposeFunction((PetscObject)st,"STFilterSetRange_C",NULL);
416: PetscObjectComposeFunction((PetscObject)st,"STFilterGetRange_C",NULL);
417: PetscObjectComposeFunction((PetscObject)st,"STFilterSetDegree_C",NULL);
418: PetscObjectComposeFunction((PetscObject)st,"STFilterGetDegree_C",NULL);
419: PetscObjectComposeFunction((PetscObject)st,"STFilterGetThreshold_C",NULL);
420: return(0);
421: }
423: SLEPC_EXTERN PetscErrorCode STCreate_Filter(ST st)
424: {
426: ST_FILTER *ctx;
427: FILTLAN_IOP iop;
428: FILTLAN_PFI pfi;
430: PetscNewLog(st,&ctx);
431: st->data = (void*)ctx;
433: st->usesksp = PETSC_FALSE;
435: ctx->inta = PETSC_MIN_REAL;
436: ctx->intb = PETSC_MAX_REAL;
437: ctx->left = 0.0;
438: ctx->right = 0.0;
439: ctx->polyDegree = 0;
440: ctx->baseDegree = 10;
442: PetscNewLog(st,&iop);
443: ctx->opts = iop;
444: iop->intervalWeights[0] = 100.0;
445: iop->intervalWeights[1] = 1.0;
446: iop->intervalWeights[2] = 1.0;
447: iop->intervalWeights[3] = 1.0;
448: iop->intervalWeights[4] = 100.0;
449: iop->transIntervalRatio = 0.6;
450: iop->reverseInterval = PETSC_FALSE;
451: iop->initialPlateau = 0.1;
452: iop->plateauShrinkRate = 1.5;
453: iop->initialShiftStep = 0.01;
454: iop->shiftStepExpanRate = 1.5;
455: iop->maxInnerIter = 30;
456: iop->yLimitTol = 0.0001;
457: iop->maxOuterIter = 50;
458: iop->numGridPoints = 200;
459: iop->yBottomLine = 0.001;
460: iop->yRippleLimit = 0.8;
462: PetscNewLog(st,&pfi);
463: ctx->filterInfo = pfi;
465: st->ops->apply = STApply_Generic;
466: st->ops->setup = STSetUp_Filter;
467: st->ops->computeoperator = STComputeOperator_Filter;
468: st->ops->setfromoptions = STSetFromOptions_Filter;
469: st->ops->destroy = STDestroy_Filter;
470: st->ops->reset = STReset_Filter;
471: st->ops->view = STView_Filter;
473: PetscObjectComposeFunction((PetscObject)st,"STFilterSetInterval_C",STFilterSetInterval_Filter);
474: PetscObjectComposeFunction((PetscObject)st,"STFilterGetInterval_C",STFilterGetInterval_Filter);
475: PetscObjectComposeFunction((PetscObject)st,"STFilterSetRange_C",STFilterSetRange_Filter);
476: PetscObjectComposeFunction((PetscObject)st,"STFilterGetRange_C",STFilterGetRange_Filter);
477: PetscObjectComposeFunction((PetscObject)st,"STFilterSetDegree_C",STFilterSetDegree_Filter);
478: PetscObjectComposeFunction((PetscObject)st,"STFilterGetDegree_C",STFilterGetDegree_Filter);
479: PetscObjectComposeFunction((PetscObject)st,"STFilterGetThreshold_C",STFilterGetThreshold_Filter);
480: return(0);
481: }