|
FILTLAN
1.0a
|
Classes | |
| struct | IntervalOptions |
| An instance of this class, taken by GetIntervals() and PolynomialFilterInterface::setFilter(), is a collection of options to determine the intervals which decides the base filter. More... | |
| struct | PolynomialFilterInfo |
| The routine GetIntervals() returns an instance of this class, which gives the information of a polynomial filter. More... | |
Functions | |
| Matrix | HermiteBaseFilterInChebyshevBasis (const Vector &intv, const int *HiLowFlags, mkIndex baseDeg) |
| Compute a base filter P(z) which is a continuous, piecewise polynomial P(z) expanded in a basis of ‘translated’ (i.e. scale-and-shift) Chebyshev polynomials in each interval. More... | |
| PolynomialFilterInfo | GetIntervals (Vector &intv, const Vector &frame, mkIndex polyDeg, mkIndex baseDeg, IntervalOptions &opts) |
| This routine determines the intervals (including the transition one(s)) by an interative process. More... | |
| PolynomialFilterInfo | GetIntervals (Vector &intervals, const Vector &frame, mkIndex polyDeg, mkIndex baseDeg) |
| Same as GetIntervals(), with default IntervalOptions. | |
| PolynomialFilterInfo GetIntervals | ( | Vector & | intv, |
| const Vector & | frame, | ||
| mkIndex | polyDeg, | ||
| mkIndex | baseDeg, | ||
| IntervalOptions & | opts | ||
| ) |
This routine determines the intervals (including the transition one(s)) by an interative process.
| frame | is a vector consisting of 4 ordered elements:
|
| baseDeg | is the left-and-right degree of the base filter for each interval. |
| polyDeg | is the (maximum possible) degree of s(z), with z*s(z) the polynomial filter. |
| intv | is the output; the jth interval is [intv(j),intv(j+1)). |
| opts | is a collection of interval options. |
The base filter P(z) is a piecewise polynomial from Hermite interpolation with degree baseDeg at each end point of intervals.
The polynomial filter Q(z) is in the form z*s(z), i.e. Q(0)==0, such that ||(1-P(z))-z*s(z)||w is minimized with s(z) a polynomial of degree up to polyDeg.
The resulting polynomial filter Q(z) satisfies Q(x)>=Q(y) for x in [frame(2),frame(3)], and y in [frame(1),frame(4)] but not in [frame(2),frame(3)].
| Matrix HermiteBaseFilterInChebyshevBasis | ( | const Vector & | intv, |
| const int * | HiLowFlags, | ||
| mkIndex | baseDeg | ||
| ) |
Compute a base filter P(z) which is a continuous, piecewise polynomial P(z) expanded in a basis of ‘translated’ (i.e. scale-and-shift) Chebyshev polynomials in each interval.
The base filter P(z) equals Pj(z) for z in the jth interval [intv(j), intv(j+1)), where Pj(z) a Hermite interpolating polynomial.
| intv | is a vector which defines the intervals. The jth interval is [intv(j), intv(j+1)). |
| HiLowFlags | determines the shape of the base filter P(z). Consider the jth interval.
|
| baseDeg | defines the smoothness of the piecewise polynomial P(z). To be precise, To be precise, the ith derivative of P(z) is zero, i.e. diP(z)/dzi==0, at all interval end points z=intv(j) for i=1,...,baseDeg. |
The base filter P(z) expanded in a basis of ‘translated’ (scale-and-shift) Chebyshev polynomials, presented by a matrix of Chebyshev coefficients pp.
For z in the jth interval [intv(j), intv(j+1)), P(z) equals Pj(z) = pp(1,j)*S0(z) + pp(2,j)*S1(z) + ... + pp(n,j)*Sn-1(z).
Here Si(z) is the ‘translated’ Chebyshev polynomial in the jth interval, with Si((z-c)/h) = Ti(z), c = (intv(j) + intv(j+1))/2, h = (intv(j+1) - intv(j))/2,
where Ti(z) is the Chebyshev polynomial of the first kind, T0(z) = 1, T1(z) = z, and Ti(z) = 2*z*Ti-1(z) - Ti-2(Z) for i>=2.
1.8.14