In file FilterDesign/FilterDesign.hh:

class FilterDesign : public FilterParse

Filter design

Inheritance:


Public Methods

explicit FilterDesign (double fsample = 1.0, const char* name = "filter")
Constructor.
explicit FilterDesign (const char* spec, double fsample = 1.0, const char* name = "filter")
Constructor.
FilterDesign (const FilterDesign& design)
Copy constructor.
virtual ~FilterDesign ()
Destructor.
FilterDesign& operator= (const FilterDesign& design)
Assignment operator.
void init (double fsample = 1.0)
Initialization.
bool isUnityGain () const
Check for unity gain filter.
void setName (const char* name)
Set name.
const char* getName () const
Get name.
void setFSample (double fsample)
Set sampling frequency.
double getFSample () const
Get sampling frequency.
double getFOut () const
Get output sampling frequency.
bool getHeterodyne () const
Get heterodyne flag.
void setPrewarp (bool set = true)
Set prewarping flag
bool getPrewarp () const
Get prewarping flag
const char* getFilterSpec () const
Get filter specification.
ligogui::TLGMultiPad* getPad () const
Get plot pad.
Pipe& get ()
Get the filter
const Pipe& get () const
Get the filter (const)
Pipe& operator) ()
Get the filter
const Pipe& operator) () const
Get the filter (const)
Pipe* copy () const
Get a copy of the filter
Pipe* release ()
Get the current filter
void reset ()
Reset the current filter
void set (const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
Set the filter
bool add (const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
Add a filter
bool filter (const char* formula)
Add a filter
bool gain (double g, const char* format = "scalar")
Add a gain
bool pole (double f, double gain, const char* plane = "s")
Add pole
bool zero (double f, double gain, const char* plane = "s")
Add zero
bool pole2 (double f, double Q, double gain, const char* plane = "s")
Add a complex pole pair
bool zero2 (double f, double Q, double gain, const char* plane = "s")
Add a complex zero pair
bool zpk (int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain, const char* plane = "s")
Add zero-pole-gain (zpk)
bool rpoly (int nnumer, const double* numer, int ndenom, const double* denom, double gain)
Add a rational polynomial (rpoly)
bool biquad (double b0, double b1, double b2, double a1, double a2)
Add a second order section
bool sos (int nba, const double* ba, const char* format = "s")
Add a list of second order sections
bool zroots (int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain = 1.0)
Add a filter from z-plane roots
bool direct (int nb, const double* b, int na, const double* a)
Add a filter from the direct form
bool ellip ( Filter_Type type, int order, double rp, double as, double f1, double f2 = 0.0 )
Add an elliptic filter
bool cheby1 ( Filter_Type type, int order, double rp, double f1, double f2 = 0.0 )
Add a chebyshev filter
bool cheby2 ( Filter_Type type, int order, double as, double f1, double f2 = 0.0 )
Add a chebyshev filter
bool butter ( Filter_Type type, int order, double f1, double f2 = 0.0 )
Add a butterworth filter
bool notch ( double f0, double Q, double depth = 0.0 )
Add a notch filter
bool resgain ( double f0, double Q, double height = 30.0 )
Add a resonant gain filter
bool comb ( double f0, double Q, double amp = 0.0, int N = 0 )
Add a comb filter
bool remez (int N, int nBand, const double* Bands, const double* Func, const double* Weight = 0)
FIR filter design using the McClellan-Parks algorithm.
bool 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.
bool difference ()
Add a difference filter
bool decimateBy2 (int N, int FilterID = 1)
Add a decimation by 2^N.
bool linefilter (double f, double T = 0.0, int fid = 1, int nT = 1)
Add a line removal filter
bool limiter (const char* type, double l1, double l2 = 0, double l3 = 0)
Add a limiter
bool multirate (const char* type, double m1, double m2 = 1E-3, double atten = 80)
Add a resampling filter
bool mixer (double fmix, double phase = 0)
Add a mixer
bool frontend (const char* file, const char* module, int section = -1)
Add a filter from a configuaryion file
bool setgain (double f, double gain)
Set filter gain
bool closeloop (double k = 1.0)
Closed loop response
bool Xfer (fComplex& coeff, double f) const
Get a transfer coefficent of a Filter.
bool Xfer (fComplex* tf, const float* freqs, int points) const
Get the transfer function of a Filter.
bool Xfer (FSeries& fseries, double Fmin = 0.0, double Fmax = 0.0, double dF = 1.0) const
Get the transfer function of a Filter.
bool 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.
bool Xfer (float* freqs, fComplex* tf, const SweptSine& param) const
Get the transfer function of a Filter.
bool response (TSeries& output, const char* type = "step", const Interval& duration = 1.0) const
Compute filter response.
bool response (TSeries& output, const Chirp& func, const Interval& duration = 1.0) const
Compute filter response.
bool response (TSeries& output, const TSeries& input) const
Compute filter response.
bool bode (double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const
Bode plot.
bool bode (const float* freqs, int points = 101) const
Bode plot.
bool bode (const SweptSine& param) const
Bode plot.
bool resp (const char* type = "step", const Interval& duration = 1.0) const
Plot filter response.
bool resp (const Chirp& func, const Interval& duration = 1.0) const
Plot filter response.
bool resp (const TSeries& input) const
Plot filter response.
bool wizard ()
Filter wizard.

Documentation

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 filter
    

    Example 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 response
    

    Example 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 noise
    

The 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 configuration
    
where
    f             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 section
    
For a more detailed description of the designer functions see the documentation of the corresponding method below.

explicit FilterDesign(double fsample = 1.0, const char* name = "filter")
Constructs a filter design class. The filter design class keeps track of the sampling frequency and the current filter.
Parameters:
fsample - Sampling rate
name - Name of filter

explicit FilterDesign(const char* spec, double fsample = 1.0, const char* name = "filter")
Constructs a filter design class. The current filter is initialized with the specified string---see also AddFilter(const char*).
Parameters:
spec - Filter design specification
fsample - Sampling rate
name - Name of filter

FilterDesign(const FilterDesign& design)
Copy constructor.

virtual ~FilterDesign()
Destructs the filter design class.

FilterDesign& operator= (const FilterDesign& design)
Assignment operator.

void init(double fsample = 1.0)
Initializs the filter design class with a new sampling rate. The output rate is set equal to the input rate and the current filter is reset to NULL. The prewarp flag is left unchanged.
Returns:
void
Parameters:
fsample - Sampling rate

bool isUnityGain() const
Check if we have a trivial unity gain filter.
Returns:
true if unity gain

void setName(const char* name)
Set the filter name.
Returns:
void
Parameters:
name - Name of filter

const char* getName() const
Get the filter name.
Returns:
filter name

void setFSample(double fsample)
Set the sampling frequency. The current filter gets reset.
Returns:
void
Parameters:
fsample - Sampling rate

double getFSample() const
Get the (intput) sampling frequency.
Returns:
Sampling rate

double getFOut() const
Get the output sampling frequency.
Returns:
Sampling rate

bool getHeterodyne() const
Get the heterodyne flag. This flag is set to false initially and will be set to true if one of the filters is a mixer. This is then an indication that the output time series will be complex for a real input time series.
Returns:
Sampling rate

void setPrewarp(bool set = true)
Set the prewarping flag
Returns:
void
Parameters:
prewarping - flag

bool getPrewarp() const
Is prewarping enabled?
Returns:
true if prewarping is enabled

const char* getFilterSpec() const
Get the filter specification.
Returns:
Filter specification

ligogui::TLGMultiPad* getPad() const
Get pad returned by the last call to a plotting function.
Returns:
plot pad

Pipe& get()
Get the filter. Returns the identity filter, if no filter has been added.
Returns:
current filter

const Pipe& get() const
Get the filter. Returns the identity filter, if no filter has been added.
Returns:
current filter

Pipe& operator) ()
Get the filter. Returns the identity filter, if no filter has been added.
Returns:
current filter

const Pipe& operator) () const
Get the filter. Returns the identity filter, if no filter has been added.
Returns:
current filter

Pipe* copy() const
Get a copy of the filter. User must delete it! Returns the identity filter, if no filter has been added.
Returns:
copy of current filter

Pipe* release()
Release the current filter and return it. Returns the identity filter, if no filter has been added. After this call the current filter is set back to an identity filter.
Returns:
current filter

void reset()
Resets the current filter to identity.
Returns:
void

void set(const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
Set the filter
Parameters:
filter - New current filter
resampling - Resampling factor
heterodyne - Filter is heterodyning

bool add(const Pipe& filter, double resampling = 1.0, bool heterodyne = false)
Add a filter to the current one.
Returns:
pointer to new current filter (owned by design class)
Parameters:
filter - Filter to be added
heterodyne - Filter is heterodyning

bool filter(const char* formula)
Make a filter from a descriptive string and add it to the current one.
Returns:
true if successful
Parameters:
spec - New filter specification

bool gain(double g, const char* format = "scalar")
Multiply the filter by the specified gain. The gain can be specified as a scalar or in dB. The default is scalar.
Returns:
true if successful
Parameters:
g - gain.
format - scalar or dB

bool pole(double f, double gain, const char* plane = "s")
Make a simple pole. Poles and zero can be specifed in the s-plane using units of rad/s (set plane to "s") or units of Hz ("f"). In both cases the frequency is expected to be POSITIVE. Alternativly, one can also specify normalized poles and zeros ("n"). Normalized poles and zeros are expected to have POSITIVE frequencies and their respective low and high frequency gains are set to 1---unless they are located a 0Hz. The default loaction is the s-plane. For the s-plane the formula is ["s"] or ["f"]. For a normalized pole ["n"] the formula is ; and ; .
Returns:
true if successful
Parameters:
f - Pole frequency.
g - gain.
plane - location where poles/zeros are specified
See Also:
IIRFilter.hh

bool zero(double f, double gain, const char* plane = "s")
Make a simple zero. Poles and zero can be specifed in the s-plane using units of rad/s (set plane to "s") or units of Hz ("f"). In both cases the frequency is expected to be POSITIVE. Alternativly, one can also specify normalized poles and zeros ("n"). Normalized poles and zeros are expected to have POSITIVE frequencies and their respective low and high frequency gains are set to 1---unless they are located a 0Hz. The default loaction is the s-plane. For the s-plane the formula is zero(f0) = (s + w0)zero(f0) = (s + 2 pi f0)zero(f0) = (1 + i f/f0)f0 > 0zero(f0) = (i f)f0 = 0$.
Returns:
true if successful
Parameters:
f - Zero frequency.
g - gain.
plane - location where poles/zeros are specified
See Also:
IIRFilter.hh

bool pole2(double f, double Q, double gain, const char* plane = "s")
Make a complex pole pair. Poles and zero can be specifed in the s-plane using units of rad/s (set plane to "s") or units of Hz ("f"). In both cases the frequency is expected to be POSITIVE. Alternativly, one can also specify normalized poles and zeros ("n"). Normalized poles and zeros are expected to have POSITIVE frequencies and their respective low and high frequency gains are set to 1---unless they are located a 0Hz. The default loaction is the s-plane. For the s-plane the formula is ["s"] or ["f"]. For a normalized complex pole pair ["n"] the formula is with . The quality factor Q must be greater than 0.5.
Returns:
true if successful
Parameters:
f - Pole frequency.
Q - Quality factor.
g - gain.
plane - location where poles/zeros are specified
See Also:
IIRFilter.hh

bool zero2(double f, double Q, double gain, const char* plane = "s")
Make a complex zero pair. Poles and zero can be specifed in the s-plane using units of rad/s (set plane to "s") or units of Hz ("f"). In both cases the frequency is expected to be POSITIVE. Alternativly, one can also specify normalized poles and zeros ("n"). Normalized poles and zeros are expected to have POSITIVE frequencies and their respective low and high frequency gains are set to 1---unless they are located a 0Hz. The default loaction is the s-plane. For the s-plane the formula is ["s"] or ["f"]. For a normalized complex pole pair ["n"] the formula is with . The quality factor Q must be greater than 0.5.
Returns:
true if successful
Parameters:
f - Pole frequency.
Q - Quality factor.
g - gain.
plane - location where poles/zeros are specified
See Also:
IIRFilter.hh

bool zpk(int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain, const char* plane = "s")
Make an filter from zeros and poles. Poles and zero can be specifed in the s-plane using units of rad/s (set plane to "s") or units of Hz ("f"). In both cases the real part of the roots are expected to be NEGATIVE. Alternativly, one can also specify normalized poles and zeros ("n"). Normalized poles and zeros are expected to have POSITIVE real parts and their respective low and high frequency gains are set to 1---unless they are located a 0Hz. The default loaction is the s-plane.

For the s-plane the formula is

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.


Returns:
true if successful
Parameters:
nzeros - Number of zeros
zero - Array of zeros
npoles - Number of poles
pole - Array of poles
gain - Gain
plane - location where poles/zeros are specified
See Also:
IIRFilter.hh

bool rpoly(int nnumer, const double* numer, int ndenom, const double* denom, double gain)
Make a filter from a rational polynomial in s. A rational polynomial in s is specified by the polynomial coefficients in the numerator and the denominator in descending order of s.

The formula is

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.


Returns:
true if successful
Parameters:
nnumer - Number of coefficients in numerator
numer - Array of numerator coefficents
ndenom - Number of coeffcients in denominator
denom - Array of denominator coeffcients
gain - Gain
See Also:
IIRFilter.hh

bool biquad(double b0, double b1, double b2, double a1, double a2)
Make an filter from the biquad coefficients and add it to the current one.
Returns:
true if successful
Parameters:
b0 - b0
b1 - b1
b2 - b2
a1 - a1
a2 - a2
See Also:
IIRdesign.hh

bool sos(int nba, const double* ba, const char* format = "s")
Make a IIR filter from a list of second order sections. If the format is 's' (standard), the coefficents are ordered like:
 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.

Returns:
true if successful
Parameters:
nba - Number of coefficients
ba - List of coefficients
format - Coefficient format
See Also:
IIRdesign.hh

bool zroots(int nzeros, const dComplex* zero, int npoles, const dComplex* pole, double gain = 1.0)
Make a IIR filter from the list of poles and zeros in the z-plane. To be stable the z-plane poles must lie within the unity cirlce.
Returns:
true if successful
Parameters:
nzeros - Number of zeros
zero - Array of zeros
npoles - Number of poles
pole - Array of poles
gain - Gain
See Also:
IIRdesign.hh

bool direct(int nb, const double* b, int na, const double* a)
Make a filter from the direct form. The direct form can be written as
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.


Returns:
true if successful
Parameters:
nb - Number of coefficients in numerator
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 )
Make an elliptic filter and add it to the current one.
Returns:
true if successful
Parameters:
type - Filter type
order - Filter order
rp - Pass band ripple (dB)
as - Stop band attenuation (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
See Also:
IIRdesign.hh

bool cheby1( Filter_Type type, int order, double rp, double f1, double f2 = 0.0 )
Make an chebyshev filter of type 1 and add it to the current one.
Returns:
true if successful
Parameters:
type - Filter type
order - Filter order
rp - Pass band ripple (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
See Also:
IIRdesign.hh

bool cheby2( Filter_Type type, int order, double as, double f1, double f2 = 0.0 )
Make an chebyshev filter of type 2 and add it to the current one.
Returns:
true if successful
Parameters:
type - Filter type
order - Filter order
as - Stop band attenuation (dB)
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
See Also:
IIRdesign.hh

bool butter( Filter_Type type, int order, double f1, double f2 = 0.0 )
Make an butterworth filter and add it to the current one.
Returns:
true if successful
Parameters:
type - Filter type
order - Filter order
f1 - Pass band edge (Hz)
f2 - Another pass band edge (Hz)
See Also:
IIRdesign.hh

bool notch( double f0, double Q, double depth = 0.0 )
Make an notch filter and add it to the current one.
Returns:
true if successful
Parameters:
f0 - Center frequency.
Q - Quality factor ( Q = (Center freq)/(Width) ).
depth - Depth of the notch (dB).
See Also:
IIRdesign.hh

bool resgain( double f0, double Q, double height = 30.0 )
Make a resonant gain filter and add it to the current one.
Returns:
true if successful
Parameters:
f0 - Center frequency.
Q - Quality factor ( Q = (Center freq)/(Width) ).
height - Height of the peak (dB).
See Also:
IIRdesign.hh

bool comb( double f0, double Q, double amp = 0.0, int N = 0 )
Make an comb filter and add it to the current one.
Returns:
true if successful
Parameters:
f0 - Fundamental frequency.
Q - Quality factor ( Q = (Center freq)/(Width) ).
amp - Depth/height of notches/peaks (dB).
N - Number of harmonics.
See Also:
IIRdesign.hh

bool remez(int N, int nBand, const double* Bands, const double* Func, const double* Weight = 0)
Design an FIR filter using the McClellan-Parks algorithm and add it to the current one. Design a optimal equiripple linear phase FIR filter from a list of bands and the desired amplitudes in each band. nBand defines the number of bands in which the response will be specified. 'Bands' contains a minimum and a maximum frequency for each band. 'Func' contains the desired amplitude for each band and 'Weight' contains an optional weighting factor for each band. The resulting filter coefficients are stored in the Filter object. The filter history is left unmodified.
Returns:
true if successful
Parameters:
N - Number of filter coefficients.
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.
See Also:
FIRdesign.hh

bool firw(int N, Filter_Type type, const char* window, double Flow, double Fhigh = 0, double Ripple = 0, double dF = 0)
Design an FIR filter using the window method and add it to the current one. Design a linear phase FIR filter using the window function technique. The windowing functions that may be used are "Rectangle", "Triangle", "Hamming", "Hanning", "Kaiser" or "Cheby". Filter types that may be designed are "LowPass", "HighPass", "BandPass" or "BandStop". N or dF will be inferred from the other parameters if omitted or set to zero in Kaiser and Chebyshev designs. The same is true of Ripple if it is omitted or set to zero in Chebyshev designs. The designed Filter coefficients replace the current coefficients. The Filter history is not affected.
Returns:
true if successful
Parameters:
N - Number of filter coefficients.
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.
See Also:
FIRdesign.hh

bool difference()
Make a difference filter and add it to the current one. Differential filter takes the difference between each successive pair of elements.
Returns:
true if successful
See Also:
Difference.hh

bool decimateBy2(int N, int FilterID = 1)
Make a decimate-by-2 filter and add it to the current one. The DecimateBy2 class decimates a TSeries by 2^N. The TSeries data are filtered before decimation to remove aliasing. This class is essentially a wrapper of the C code written for LIGO GDS by Peter Fritschel.

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:
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)
The antialiasing filters add an effective delay of , 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.

Returns:
true if successful
Parameters:
N - number of stages
FilterID - filter type
See Also:
DecimateBy2.hh

bool linefilter(double f, double T = 0.0, int fid = 1, int nT = 1)
Make a line removal filter and add it to the current one. The LineFilter class containes methods to track and remove quasi-monochromatic lines. A filter id of 0 corresponds to a pure "comb" filter, a filter if of 1 corresponds to an optimal Wiener filter with 1-st method of noise estimation (default), and a filter id of -1
Returns:
true if successful
Parameters:
f - Approximate line frequency
T - time interval to remove lines.
id - Select filter ID (0/1), 1 by default
nT - number of subdivisions of time interval w
See Also:
LineFilter.hh

bool limiter(const char* type, double l1, double l2 = 0, double l3 = 0)
Make a limiter filter and add it to the current one. The type can be "" for no limits, "val" for high/low limits on the sampled values, "sym" for a symmetric limits on the values, "slew" for a slew rate limit and "val/slew" or "sym/slew" for both. If the limiter type is "val", l1 is the lower bound and l2 is the upper bound, for "sym" l1 is the upper bound and -l1 is the lower bound, for "slew" l1 is the slew rate limit in counts/s, for "val/slew" l1, l2 and l3 are the lower bound, the upper bound and the slew rate limit, respectively, and for "sym/slew" l1 and l2 are the symmetric limit and the slew rate limit, respectively.
Returns:
true if successful
Parameters:
type - Limiter type
l1 - First limit
l2 - Second limit
l3 - Third limit
See Also:
Limiter.hh

bool multirate(const char* type, double m1, double m2 = 1E-3, double atten = 80)
Make a resampling filter and add it to the current one. The filter type is either "abs" or "rel" depending if the required resampling rate is specified realtive or absolute. For absolute rate specifications the first parameter is the desired sampling rate and the second parameter is the maximally allowed error. For a relative rate specifiaction the first parameter is the interpolation factor and the second parameter is the decimation factor (both are integer parameters). The third parameter is in both cases the required stop band attenuation.

Returns:
true if successful
Parameters:
type - Multirate type
m1 - First parameter (desired fS or interpolation factor)
m2 - Second parameter (error or decimation factor)
attn - Stopband attenuation
See Also:
MultiRate.hh

bool mixer(double fmix, double phase = 0)
Make a mixer and add it to the current filter. This will transform a real time series into a complex one. The mixer will heterodyne (multiply) a time series by exp(2*pi*I*fc*t+phase).
Returns:
true if successful
Parameters:
fmix - Mixer frequency
phase - Initial phase
See Also:
Mixer.hh

bool frontend(const char* file, const char* module, int section = -1)
Read filter coefficients from front-end file, make filter and add it to the current one.
Returns:
true if successful
Parameters:
file - Front-end configuration file
module - Module name or number
section - Filter section or -1 for all

bool setgain(double f, double gain)
Set the filter gain at specified frequency.
Returns:
true if successful
Parameters:
f - Frequency
gain - Set point

bool closeloop(double k = 1.0)
Form the closed loop response of the current filter.
Returns:
true if successful
Parameters:
k - Additional gain: 1/(1+k*G(f))

bool Xfer(fComplex& coeff, double f) const
The transfer coefficient of the filter at the specified frequency is calculated and returned as a complex number.
Returns:
true if successful
Parameters:
coeff - a complex number representing the Filter response at the specified frequency (return)
F - Frequency at which to sample the transfer function.

bool Xfer(fComplex* tf, const float* freqs, int points) const
The transfer function of the filter at the specified frequency points is calculated and returned as a complex array. The frequency points are user supplied. The return array must be at least of length points.
Returns:
true if successful
Parameters:
tf - Transfer function (return)
freqs - Frequency points
points - Number of points.

bool Xfer(FSeries& fseries, double Fmin = 0.0, double Fmax = 0.0, double dF = 1.0) const
The transfer function of the filter in the specified frequency interval is calculated and returned as a complex frequency series.
Returns:
true if successful
Parameters:
fseries - a complex FSeries containing the Filter response at each frequency step (return)
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
The transfer function of the filter at the specified frequency points is calculated and returned as a complex array. The 'sweep type' can be "log" for logarithmically spaced frequency points or "linear" for linear spacing. The specifed arrays must be at least of length points.
Returns:
true if successful
Parameters:
freqs - Frequency points (return)
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
The transfer function of the filter in the specified frequency interval is calculated using a swept sine. This method works with all filters but may take significantly longer than the direct method implemented for FIR and IIR filters.
Returns:
true if successful
Parameters:
freqs - Frequency points (return)
tf - Transfer function (return)
param - Swept sine parameters

bool response(TSeries& output, const char* type = "step", const Interval& duration = 1.0) const
Compute the response of the filter. The reqested reponse type must be either "step" for a step response, "impulse" for an impulse response, or "ramp" for a ramp response.
Returns:
time series of response
Parameters:
output - Time series (return)
type - Response type
duration - Length of returned time series.

bool response(TSeries& output, const Chirp& func, const Interval& duration = 1.0) const
Compute the response of the filter.
Returns:
time series of response
Parameters:
output - Time series (return)
func - Signal function
duration - Length of returned time series.

bool response(TSeries& output, const TSeries& input) const
Compute the response of the filter.
Returns:
time series of response
Parameters:
output - Time series (return)
input - Input time series

bool bode(double fstart, double fstop = 0.0, int points = 101, const char* type = "log") const
Bode plot of filter. Make a bode plot of the filter.
Returns:
true if successful
Parameters:
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 bode(const float* freqs, int points = 101) const
Bode plot of filter. Make a bode plot of the filter.
Returns:
true if successful
Parameters:
freqs. -
points - Number of points.
type - Type of frequency spacing.

bool bode(const SweptSine& param) const
Bode plot of filter. Make a bode plot of the filter using a swept sine method. This method works with all filters but may take significantly longer than the direct method implemented for FIR and IIR filters.
Returns:
true if successful
Parameters:
param - Swept sine parameters
freqs - Frequency points (return)
tf - Transfer function (return)

bool resp(const char* type = "step", const Interval& duration = 1.0) const
Plot the response of the filter. The reqested reponse type must be either "step" for a step response, "impulse" for an impulse response, or "ramp" for a ramp response.
Returns:
true if successful
Parameters:
type - Response type
duration - Length of returned time series.

bool resp(const Chirp& func, const Interval& duration = 1.0) const
Plot the response of the filter.
Returns:
true if successful
Parameters:
func - Signal function
duration - Length of returned time series.

bool resp(const TSeries& input) const
Plot the response of the filter.
Returns:
true if successful
Parameters:
input - Input time series

bool wizard()
Interactive design of a filter. Not yet implemented.
Returns:
true if a new filter was designed


This class has no child classes.
Author:
Written July 2002 by Ed Daw, Masahiro Ito, Daniel Sigg
Version:
1.0

alphabetic index hierarchy of classes


Please send questions and comments to zweizig_j@ligo.caltech.edu


generated by doc++