class FilterDesign : public FilterParse Filter design
| | FilterDesign (double fsample = 1.0, const char* name = "filter") Constructor. |
| | FilterDesign (const char* spec, double fsample = 1.0, const char* name = "filter") Constructor. |
| | FilterDesign (const FilterDesign& design) Copy constructor. |
| | ~FilterDesign () Destructor. |
| | operator= (const FilterDesign& design) Assignment operator. |
| | init (double fsample = 1.0) Initialization. |
| | isUnityGain () const Check for unity gain filter. |
| | setName (const char* name) Set name. |
| | getName () const Get name. |
| | setFSample (double fsample) Set sampling frequency. |
| | getFSample () const Get sampling frequency. |
| | getFOut () const Get output sampling frequency. |
| | getHeterodyne () const Get heterodyne flag. |
| | setPrewarp (bool set = true) Set prewarping flag |
| | getPrewarp () const Get prewarping flag |
| | getFilterSpec () const Get filter specification. |
| | getPad () const Get plot pad. |
| | get () Get the filter |
| | get () const Get the filter (const) |
| | operator) () Get the filter |
| | operator) () const Get the filter (const) |
| | copy () const Get a copy of the filter |
| | release () Get the current filter |
| | reset () Reset the current filter |
| | set (const Pipe& filter, double resampling = 1.0, bool heterodyne = false) Set the filter |
| | add (const Pipe& filter, double resampling = 1.0, bool heterodyne = false) Add a filter |
| | filter (const char* formula) Add a filter |
| | gain (double g, const char* format = "scalar") Add a gain |
| | pole (double f, double gain, const char* plane = "s") Add pole |
| | zero (double f, double gain, const char* plane = "s") Add zero |
| | pole2 (double f, double Q, double gain, const char* plane = "s") Add a complex pole pair |
| | zero2 (double f, double Q, double gain, const char* plane = "s") Add a complex zero pair |
| | zpk (int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain, const char* plane = "s") Add zero-pole-gain (zpk) |
| | rpoly (int nnumer, const double* numer, int ndenom, const double* denom, double gain) Add a rational polynomial (rpoly) |
| | biquad (double b0, double b1, double b2, double a1, double a2) Add a second order section |
| | sos (int nba, const double* ba, const char* format = "s") Add a list of second order sections |
| | zroots (int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain = 1.0) Add a filter from z-plane roots |
| | direct (int nb, const double* b, int na, const double* a) Add a filter from the direct form |
| | ellip ( Filter_Type type, int order, double rp, double as, double f1, double f2 = 0.0 ) Add an elliptic filter |
| | cheby1 ( Filter_Type type, int order, double rp, double f1, double f2 = 0.0 ) Add a chebyshev filter |
| | cheby2 ( Filter_Type type, int order, double as, double f1, double f2 = 0.0 ) Add a chebyshev filter |
| | butter ( Filter_Type type, int order, double f1, double f2 = 0.0 ) Add a butterworth filter |
| | notch ( double f0, double Q, double depth = 0.0 ) Add a notch filter |
| | resgain ( double f0, double Q, double height = 30.0 ) Add a resonant gain filter |
| | comb ( double f0, double Q, double amp = 0.0, int N = 0 ) Add a comb filter |
| | remez (int N, int nBand, const double* Bands, const double* Func, const double* Weight = 0) FIR filter design using the McClellan-Parks algorithm. |
| | firw (int N, Filter_Type type, const char* window, double Flow, double Fhigh = 0, double Ripple = 0, double dF = 0) Filter design from window function. |
| | difference () Add a difference filter |
| | decimateBy2 (int N, int FilterID = 1) Add a decimation by 2^N. |
| | linefilter (double f, double T = 0.0, int fid = 1, int nT = 1) Add a line removal filter |
| | limiter (const char* type, double l1, double l2 = 0, double l3 = 0) Add a limiter |
| | multirate (const char* type, double m1, double m2 = 1E-3, double atten = 80) Add a resampling filter |
| | mixer (double fmix, double phase = 0) Add a mixer |
| | frontend (const char* file, const char* module, int section = -1) Add a filter from a configuaryion file |
| | setgain (double f, double gain) Set filter gain |
| | closeloop (double k = 1.0) Closed loop response |
| | Xfer (fComplex& coeff, double f) const Get a transfer coefficent of a Filter. |
| | Xfer (fComplex* tf, const float* freqs, int points) const Get the transfer function of a Filter. |
| | Xfer (FSeries& fseries, double Fmin = 0.0, double Fmax = 0.0, double dF = 1.0) const Get the transfer function of a Filter. |
| | Xfer (float* freqs, fComplex* tf, double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const Get the transfer function of a Filter. |
| | Xfer (float* freqs, fComplex* tf, const SweptSine& param) const Get the transfer function of a Filter. |
| | response (TSeries& output, const char* type = "step", const Interval& duration = 1.0) const Compute filter response. |
| | response (TSeries& output, const Chirp& func, const Interval& duration = 1.0) const Compute filter response. |
| | response (TSeries& output, const TSeries& input) const Compute filter response. |
| | bode (double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const Bode plot. |
| | bode (const float* freqs, int points = 101) const Bode plot. |
| | bode (const SweptSine& param) const Bode plot. |
| | resp (const char* type = "step", const Interval& duration = 1.0) const Plot filter response. |
| | resp (const Chirp& func, const Interval& duration = 1.0) const Plot filter response. |
| | resp (const TSeries& input) const Plot filter response. |
| | wizard () Filter wizard. |
Design a filter. The filter design class allows to design complex filters. It is not a filter class by itself but rather a "factory" class that produces filters depending on a user specifiaction. It supports:1. IIR: elliptic, chebyshev, butterworth, zero/pole/gain, notches, resonant gain, pedndulum and biquad coefficients. 2. FIR: window and McClellan-Parks algorithms. 3. Others: line tracker, decimation and resampling.Any number of filters can be added in series. Filters can be designed by calling their corresponding creation methods (such as ellip), by parsing a string formula or by reading the coefficients of a configuration file used by the front-end code.The filter design class also supports a configuration line parser which follows the matlab syntax. This way a filter can be written by a product of designer functions. For example, "zpk([100;100],[1;1],1) /. setgain(0,1)" represents a filter consisting of two poles at 100Hz, two zeros at 1Hz and a dc gain of 1. For a list of supported designer functions see below.
For filters that support it such as IIR, FIR and combinations of them the filter design class can also calculate the transfer function and generate a Bode plot. All filters can tested for step, impulse and ramp response.
By default the filter design class will return a unity gain filter (identity operator), if no filter has been added.
Example 1: Make an elliptic high pass filter FilterDesign fdesign (16384); // specify sampling frequency fdesign.ellip (kLowPass,8,1,60,100); // make elliptic low pass // 8th order, 1dB ripple, 60dB attenuation, 100Hz corner Pipe* pipe = fdesign.release(); // get filterExample 2: Make a whitening filter with a resonant gain, generate a bode plot and investigate its step response FilterDesign fdesign ("zpk([1;1],[100;100],1)*resgain(12,30)", 16384); // specify filter string and sampling frequency fdesign.bode (0.1, 1000, 1001); // make a bode plot fdesign.response ("step", 10.0) // plot step responseExample 3: Read DARM filter from online system, add a 60Hz comb and test it on random noise FilterDesign fdesign (16384); fdesign.frontend ("H2LSC.txt", "DARM"); fdesign.comb (60, 100); // add Q=100 comb TSeries test (Time (0), Interval (1./fdesign.getFSample()), 10*16384, GaussNoise (fdesign.getFSample())); // make some random noise fdesign.response (test) // plot filter response to random noiseThe filter types are: kLowPass for low pass, kHighPass for high pass, kBandPass for band pass and kBandStop for a band stop filter, respectively. When using a string specification, the corresponding filter types are "LowPass", "HighPass", "BandPass" and "BandStop".
FILTER FORMULA SPECIFICATION
Filter formulas are written as a product of real numbers and designer functions. Optionally one can add a gain condition at the end. Designer functions ar of the form "f(a1,a2,...aN)" where "f" must be a supported function and "a1", .. "aN" must be the corresponding arguments. A complex number argument is written as "real+i*imag" or "real+imag*i". Number arguments must be constants and expressions are not allowed. A vector of numbers is written as "[e1;e2;...eN]". The gain condition is introduce by the "/." symbol and it must be followed by the "setgain(f,g)" where g is the desired gain and f is the frequency at which it should be true. The multiplication operator "*" is optional and can be omitted or replaced by a space.
List of designer functions:
gain(g,gfmt='scalar') multiply by g pole(f,g=1,pl='s') a simple pole pole(f,pl='s') a simple pole zero(f,g=1,pl='s') a simple zero zero(f,pl='s') a simple zero pole2(f,Q,g=1,pl='s') a complex pole pair pole2(f,Q,pl='s') a complex pole pair zero2(f,Q,g=1,pl='s') a complex zero pair zero2(f,Q,pl='s') a complex zero pair zpk([z1;...zN],[p1;...pM],g=1,pl) a set of s-plane poles and zeros rpoly([z1;...zN],[p1;...pM],g=1) a rational polynomial biquad(b0,b1,b2,a1,a2) second order section sos(g,[b1,b2,a1,a2,...],fmt='s') list of second order sections zroots([z1;...zN],[p1;...pM],g=1) a set of z-plane poles and zeros direct([b0;...bN],[a1;...aM]) direct form ellip(type,order,rp,as,f1,f2=0) elliptic filter butter(type,order,f1,f2=0) butterworth filter cheby1(type,order,rp,f1,f2=0) chebyshev type 1 filter cheby2(type,order,as,f1,f2=0) chebyshev type 2 filter notch(f,Q,depth=0) notch filter resgain(f,q,height=0) resonant gain filter comb(f,Q,depth=0,hN=0) comb filter remez(order,[...],[...],[...]) FIR filter according to McClellan-Parks firw(...) FIR filter with window method difference() difference filter limiter(type,l1,l2=0,l3=0) limiter decimateBy2(stages,id=1) decimate by factors of 2 (FIR) multirate(type,m1,m2,attn) resample the time series mixer(fmix,phase) mixer linefilter(...) line removal frontend(file,module,sect=-1) read front-end configurationwheref corner frequency, g gain, gfmt gain format Q quality factor z zero location p pole location pl pole/zero location a1,a2 denominator coefficents b0,b1,b2 numerator coefficents fmt format of second order sections type filter type order filter order rp passband ripple as stopband attenuation f1 corner frequency f2 corner frequency depth depth in dB height depth in dB hN Number of harmonics l1,l2,l3 limits (bounds and/or slew rate) stages number of stages id filter id m1,m2 resampling parameters (abs or rel) attn attenuation in dB fmix heterodyne frequency phase phase file file name module filter module name sect filter sectionFor a more detailed description of the designer functions see the documentation of the corresponding method below.
For the s-plane the formula is
The formula is
DecimateBy2 is based on the DMT FilterBase class and as
such is can be treated as a generic DMT Filter. The
decimation may be applied to the TSeries with either the
apply method or the parentheses operator. These methods
return a complex TSeries if the input series is complex,
and a single precision floating point TSeries for all other
input series data types. Four anti-aliasinglow-pass FIR filters are available,
each of which cuts of at 0.9 F_{Nyquist}. The options are
tabulated below:
explicit FilterDesign(const char* spec, double fsample = 1.0, const char* name = "filter")
fsample - Sampling rate
name - Name of filter
FilterDesign(const FilterDesign& design)
virtual ~FilterDesign()
FilterDesign& operator= (const FilterDesign& design)
void init(double fsample = 1.0)
bool isUnityGain() const
void setName(const char* name)
const char* getName() const
void setFSample(double fsample)
double getFSample() const
double getFOut() const
bool getHeterodyne() const
void setPrewarp(bool set = true)
bool getPrewarp() const
const char* getFilterSpec() const
ligogui::TLGMultiPad* getPad() const
Pipe& get()
const Pipe& get() const
Pipe& operator) ()
const Pipe& operator) () const
Pipe* copy() const
Pipe* release()
void reset()
void set(const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
resampling - Resampling factor
heterodyne - Filter is heterodyning
bool add(const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
heterodyne - Filter is heterodyning
bool filter(const char* formula)
bool gain(double g, const char* format = "scalar")
format - scalar or dB
bool pole(double f, double gain, const char* plane = "s")
["s"]
or
["f"]. For a normalized pole
["n"] the formula is
;
and
;
.
g - gain.
plane - location where poles/zeros are specified
bool zero(double f, double gain, const char* plane = "s")
zero(f0) = (s + 2 pi f0)
zero(f0) = (1 + i f/f0)
f0 > 0
zero(f0) = (i f)
f0 = 0$.
g - gain.
plane - location where poles/zeros are specified
bool pole2(double f, double Q, double gain, const char* plane = "s")
["s"] or
["f"].
For a normalized complex pole pair ["n"] the formula is
with
.
The quality factor Q must be greater than 0.5.
Q - Quality factor.
g - gain.
plane - location where poles/zeros are specified
bool zero2(double f, double Q, double gain, const char* plane = "s")
["s"] or
["f"].
For a normalized complex pole pair ["n"] the formula is
with
.
The quality factor Q must be greater than 0.5.
Q - Quality factor.
g - gain.
plane - location where poles/zeros are specified
bool zpk(int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain, const char* plane = "s")
zpk(s) = k \frac{(s - s_1)(s - s_2) \cdots}{(s - p_1)(s - p_2) \cdots}$$
with $s_1$, $s_2$,... the loactions of the zeros and
$p_1$, $p_2$,... the loaction of the poles. By replacing
s_1->$2 \pi fz_1$and p_1->$2 \pi fp_1$ one obatines the "f"
representation. For normalize poles and zeros ["n"] one uses
poles of the form $pole(f0) = 1/(1 + i f/f0)$; $f0 > 0$ and
$pole(f0) = 1/(i f)$; $f0 = 0$. The zeros are then of the form
$zero(f0) = (1 + i f/f0)$; $f0 > 0$ and $zero(f0) = (i f)$;
$f0 = 0$. Poles and zeros with a non-zero imaginary part must
come in pairs of complex conjugates.
zero - Array of zeros
npoles - Number of poles
pole - Array of poles
gain - Gain
plane - location where poles/zeros are specified
bool rpoly(int nnumer, const double* numer, int ndenom, const double* denom, double gain)
zpk(s) = k \frac{a_n s^{n_z} + a_{n-1} s^{n_z - 1} \cdots}
{b_n s^{n_p} + b_{n-1} s^{n_p - 1} \cdots}$$
where $a_n$, $a_{n-1}$,..., $a_0$ are the coeffciients of the
polynomial in the numerator and $b_n$, $b_{n-1}$,..., $b_0$ are
the coefficients of the polynomial in the denominator.
The polynomial coefficients have to be real.
numer - Array of numerator coefficents
ndenom - Number of coeffcients in denominator
denom - Array of denominator coeffcients
gain - Gain
bool biquad(double b0, double b1, double b2, double a1, double a2)
b1 - b1
b2 - b2
a1 - a1
a2 - a2
bool sos(int nba, const double* ba, const char* format = "s")
gain, b1_1, b2_1, a1_1, a2_1, b1_2, b2_2, a1_2, a2_2,... $$
whereas for the the format 'o' (online) the order is
$$ gain, a1_1, a2_1, b1_1, b2_1, a1_2, a2_2, b1_2, b2_2,... $$
The number of coefficients must be 4 times the number of second
order sections plus one.
ba - List of coefficients
format - Coefficient format
bool zroots(int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain = 1.0)
zero - Array of zeros
npoles - Number of poles
pole - Array of poles
gain - Gain
bool direct(int nb, const double* b, int na, const double* a)
H(z)=\Sum_{k=0}^{nb}b_k z^{-k} / (1-\Sum_{k=1}^{na}a_k z^{-k})$$
Cascaded second order sections are formed by finding the roots
of the direct form. The specified coefficients are $b_0$, $b_1$,...,
$b_{nb}$ for the numerator and $a_1$, $a_2$,..., $a_{na}$ for
the denominator. The argument $nb$ specifies the number of b
coeffcients supplied to the function minus 1, whereas $na$ is
exactly the number of coeffcients since $a_0$ is always 1 and
omitted from the list.
Avoid the direct form since even fairly simple filters will
run into precision problems.
b - Array of b coefficents
na - Number of coeffcients in denominator
a - Array of denominator coeffcients exclusive a0
bool ellip( Filter_Type type, int order, double rp, double as, double f1, double f2 = 0.0 )
order - Filter order
rp - Pass band ripple (dB)
as - Stop band attenuation (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
bool cheby1( Filter_Type type, int order, double rp, double f1, double f2 = 0.0 )
order - Filter order
rp - Pass band ripple (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
bool cheby2( Filter_Type type, int order, double as, double f1, double f2 = 0.0 )
order - Filter order
as - Stop band attenuation (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
bool butter( Filter_Type type, int order, double f1, double f2 = 0.0 )
order - Filter order
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
bool notch( double f0, double Q, double depth = 0.0 )
Q - Quality factor ( Q = (Center freq)/(Width) ).
depth - Depth of the notch (dB).
bool resgain( double f0, double Q, double height = 30.0 )
Q - Quality factor ( Q = (Center freq)/(Width) ).
height - Height of the peak (dB).
bool comb( double f0, double Q, double amp = 0.0, int N = 0 )
Q - Quality factor ( Q = (Center freq)/(Width) ).
amp - Depth/height of notches/peaks (dB).
N - Number of harmonics.
bool remez(int N, int nBand, const double* Bands, const double* Func, const double* Weight = 0)
nBand - Number of bands.
Bands - Band limit frequencies in Hz.
Func - Desired response in each band.
Weight - Weighting factor for each band. If not specified,
Weight is assumed to be 1.0 for each band.
bool firw(int N, Filter_Type type, const char* window, double Flow, double Fhigh = 0, double Ripple = 0, double dF = 0)
window - Type of window on which to base the design.
type - Type of filter.
Flow - Transition frequency of a low-pass or
high-pass filter, or lower edge of a band-pass
or band-stop filter.
Fhigh - Upper edge of a band-pass or band-stop filter.
Ripple - Desired attenuation for Kaiser filters or the
desired pass-band ripple for Chebyshev filters
(in dB).
dF - Transition band width in Chebyshev filters.
bool difference()
bool decimateBy2(int N, int FilterID = 1)
The antialiasing filters add an effective delay of
ID Order Comments
1 42 Least-Squares design, ripple ~.015 (-0.2 dB)
2 42 Equiripple design, ripple ~0.06 dB
3 22 Least-squares design, ripple ~0.1 (-1 dB)
4 82 Least-squares design, ripple ~0.0006 (-0.01 dB)
,
where N is the number of stages (the decimation factor
exponent) and
is the filter order. This delay is NOT
added to the TSeries start time.
FilterID - filter type
bool linefilter(double f, double T = 0.0, int fid = 1, int nT = 1)
T - time interval to remove lines.
id - Select filter ID (0/1), 1 by default
nT - number of subdivisions of time interval w
bool limiter(const char* type, double l1, double l2 = 0, double l3 = 0)
l1 - First limit
l2 - Second limit
l3 - Third limit
bool multirate(const char* type, double m1, double m2 = 1E-3, double atten = 80)
m1 - First parameter (desired fS or interpolation factor)
m2 - Second parameter (error or decimation factor)
attn - Stopband attenuation
bool mixer(double fmix, double phase = 0)
phase - Initial phase
bool frontend(const char* file, const char* module, int section = -1)
module - Module name or number
section - Filter section or -1 for all
bool setgain(double f, double gain)
gain - Set point
bool closeloop(double k = 1.0)
bool Xfer(fComplex& coeff, double f) const
F - Frequency at which to sample the transfer function.
bool Xfer(fComplex* tf, const float* freqs, int points) const
freqs - Frequency points
points - Number of points.
bool Xfer(FSeries& fseries, double Fmin = 0.0, double Fmax = 0.0, double dF = 1.0) const
Fmin - Minimum frequency at which to sample the
transfer function.
Fmax - Maximum frequency at which to sample the
transfer function.
dF - Frequency step.
bool Xfer(float* freqs, fComplex* tf, double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const
tf - Transfer function (return)
fstart - Minimum frequency at which to sample the
transfer function.
fstop - Maximum frequency at which to sample the
transfer function.
points - Number of points.
type - Type of frequency spacing.
bool Xfer(float* freqs, fComplex* tf, const SweptSine& param) const
tf - Transfer function (return)
param - Swept sine parameters
bool response(TSeries& output, const char* type = "step", const Interval& duration = 1.0) const
type - Response type
duration - Length of returned time series.
bool response(TSeries& output, const Chirp& func, const Interval& duration = 1.0) const
func - Signal function
duration - Length of returned time series.
bool response(TSeries& output, const TSeries& input) const
input - Input time series
bool bode(double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const
fstop - Maximum frequency at which to sample the
transfer function.
points - Number of points.
type - Type of frequency spacing.
bool bode(const float* freqs, int points = 101) const
points - Number of points.
type - Type of frequency spacing.
bool bode(const SweptSine& param) const
freqs - Frequency points (return)
tf - Transfer function (return)
bool resp(const char* type = "step", const Interval& duration = 1.0) const
duration - Length of returned time series.
bool resp(const Chirp& func, const Interval& duration = 1.0) const
duration - Length of returned time series.
bool resp(const TSeries& input) const
bool wizard()
alphabetic index hierarchy of classes
Please send questions and comments to zweizig_j@ligo.caltech.edu
generated by doc++