Functions

src/filter.c File Reference

#include "filter.h"
#include "eeg.h"
#include "distances.h"
#include "fidlib/fidlib.h"

Go to the source code of this file.

Functions

EEGeeg_filter_fidlib (EEG *eeg, const char *spec, bool alloc)
EEGeeg_filter_running_median (EEG *eeg, int win, bool alloc)
EEGeeg_filter_weighted_running_median (EEG *eeg, int win, bool alloc)
double * moving_average (double *s, int n, int win)
double * running_median (double *d, int n, int win)
double * weighted_running_median (double *d, int n, int win)

Function Documentation

EEG* eeg_filter_fidlib ( EEG eeg,
const char *  spec,
bool  alloc 
)

Filters EEG data using Fidlib from Jim Peters (http://uazu.net/fidlib/). Filtering is done in place or on a copy of the struct (depending on alloc).

Parameters:
eeg 
spec filter-specification as given in Fidlib
alloc allocate or not allocate, that's the question

     The spec consists of a series of letters usually followed by the order
     of the filter and then by any other parameters required, preceded by
     slashes.  For example:
     
     LpBu4/20.4    Lowpass butterworth, 4th order, -3.01dB at 20.4Hz
     BpBu2/3-4     Bandpass butterworth, 2nd order, from 3 to 4Hz
     BpBu2/=3-4    Same filter, but adjusted exactly to the range given
     BsRe/1000/10  Bandstop resonator, Q=1000, frequency 10Hz
     
     Note that filter frequencies should be specified in the same units as
     the sampling rate, i.e. normally in Hz.  However, all filters work
     internally with the ratio of freq/sampling_rate, so you can use a
     sampling rate of 1.0 and frequency values of 0 to 0.5 (Nyquist) if you
     prefer.

     The auto-adjust option, selected by prefixing the frequency or frequency
     range with '=', automatically adjusts the -3.01dB points to the given
     frequencies with a kind of binary search.  
     This might be useful with some lower-order filters where the
     -3.01dB points don't come out quite right otherwise.

Available Filters:

     BpRe/<value>/<freq>
    Bandpass resonator, Q=<value> (0 means Inf), frequency <freq>
     BsRe/<value>/<freq>
    Bandstop resonator, Q=<value> (0 means Inf), frequency <freq>
     ApRe/<value>/<freq>
    Allpass resonator, Q=<value> (0 means Inf), frequency <freq>
     Pi/<freq>
    Proportional-integral filter, frequency <freq>
     PiZ/<freq>
    Proportional-integral filter, matched z-transform, frequency <freq>
     LpBe<order>/<freq>
    Lowpass Bessel filter, order <order>, -3.01dB frequency <freq>
     HpBe<order>/<freq>
    Highpass Bessel filter, order <order>, -3.01dB frequency <freq>
     BpBe<order>/<range>
    Bandpass Bessel filter, order <order>, -3.01dB frequencies <range>
     BsBe<order>/<range>
    Bandstop Bessel filter, order <order>, -3.01dB frequencies <range>
     LpBu<order>/<freq>
    Lowpass Butterworth filter, order <order>, -3.01dB frequency <freq>
     HpBu<order>/<freq>
    Highpass Butterworth filter, order <order>, -3.01dB frequency <freq>
     BpBu<order>/<range>
    Bandpass Butterworth filter, order <order>, -3.01dB frequencies <range>
     BsBu<order>/<range>
    Bandstop Butterworth filter, order <order>, -3.01dB frequencies <range>
     LpCh<order>/<value>/<freq>
    Lowpass Chebyshev filter, order <order>, passband ripple <value>dB,
     -3.01dB frequency <freq>
     HpCh<order>/<value>/<freq>
    Highpass Chebyshev filter, order <order>, passband ripple <value>dB,
     -3.01dB frequency <freq>
     BpCh<order>/<value>/<range>
    Bandpass Chebyshev filter, order <order>, passband ripple <value>dB,
     -3.01dB frequencies <range>
     BsCh<order>/<value>/<range>
    Bandstop Chebyshev filter, order <order>, passband ripple <value>dB,
     -3.01dB frequencies <range>
     LpBeZ<order>/<freq>
    Lowpass Bessel filter, matched z-transform, order <order>, -3.01dB
     frequency <freq>
     HpBeZ<order>/<freq>
    Highpass Bessel filter, matched z-transform, order <order>, -3.01dB
     frequency <freq>
     BpBeZ<order>/<range>
    Bandpass Bessel filter, matched z-transform, order <order>, -3.01dB
     frequencies <range>
     BsBeZ<order>/<range>
    Bandstop Bessel filter, matched z-transform, order <order>, -3.01dB
     frequencies <range>
     LpBuZ<order>/<freq>
    Lowpass Butterworth filter, matched z-transform, order <order>, -3.01dB
     frequency <freq>
     HpBuZ<order>/<freq>
    Highpass Butterworth filter, matched z-transform, order <order>, -3.01dB
     frequency <freq>
     BpBuZ<order>/<range>
    Bandpass Butterworth filter, matched z-transform, order <order>, -3.01dB
     frequencies <range>
     BsBuZ<order>/<range>
    Bandstop Butterworth filter, matched z-transform, order <order>, -3.01dB
     frequencies <range>
     LpChZ<order>/<value>/<freq>
    Lowpass Chebyshev filter, matched z-transform, order <order>, passband
     ripple <value>dB, -3.01dB frequency <freq>
     HpChZ<order>/<value>/<freq>
    Highpass Chebyshev filter, matched z-transform, order <order>, passband
     ripple <value>dB, -3.01dB frequency <freq>
     BpChZ<order>/<value>/<range>
    Bandpass Chebyshev filter, matched z-transform, order <order>, passband
     ripple <value>dB, -3.01dB frequencies <range>
     BsChZ<order>/<value>/<range>
    Bandstop Chebyshev filter, matched z-transform, order <order>, passband
     ripple <value>dB, -3.01dB frequencies <range>
     LpBuBe<order>/<value>/<freq>
    Lowpass Butterworth-Bessel <value>% cross, order <order>, -3.01dB
     frequency <freq>
     LpBq<optional-order>/<value>/<freq>
    Lowpass biquad filter, order <order>, Q=<value>, -3.01dB frequency <freq>
     HpBq<optional-order>/<value>/<freq>
    Highpass biquad filter, order <order>, Q=<value>, -3.01dB frequency
     <freq>
     BpBq<optional-order>/<value>/<freq>
    Bandpass biquad filter, order <order>, Q=<value>, centre frequency <freq>
     BsBq<optional-order>/<value>/<freq>
    Bandstop biquad filter, order <order>, Q=<value>, centre frequency <freq>
     ApBq<optional-order>/<value>/<freq>
    Allpass biquad filter, order <order>, Q=<value>, centre frequency <freq>
     PkBq<optional-order>/<value>/<value>/<freq>
    Peaking biquad filter, order <order>, Q=<value>, dBgain=<value>,
     frequency <freq>
     LsBq<optional-order>/<value>/<value>/<freq>
    Lowpass shelving biquad filter, S=<value>, dBgain=<value>, frequency
     <freq>
     HsBq<optional-order>/<value>/<value>/<freq>
    Highpass shelving biquad filter, S=<value>, dBgain=<value>, frequency
     <freq>
     LpBl/<freq>
    Lowpass Blackman window, -3.01dB frequency <freq>
     LpHm/<freq>
    Lowpass Hamming window, -3.01dB frequency <freq>
     LpHn/<freq>
    Lowpass Hann window, -3.01dB frequency <freq>
     LpBa/<freq>
    Lowpass Bartlet (triangular) window, -3.01dB frequency <freq>

Definition at line 312 of file filter.c.

EEG* eeg_filter_running_median ( EEG eeg,
int  win,
bool  alloc 
)

Simple Running Median filter. Is applied to each channel and trial separately.

Parameters:
s input eeg
win the window size in which to calculate the median
alloc TRUE: a new struct is returned; FALSE: in-place modification
Returns:
filtered EEG

Definition at line 36 of file filter.c.

EEG* eeg_filter_weighted_running_median ( EEG eeg,
int  win,
bool  alloc 
)

Weighted Running Median filter. Is applied to each channel and trial separately.

Parameters:
s input eeg
win the window size in which to calculate the median
alloc TRUE: a new struct is returned; FALSE: in-place modification
Returns:
filtered EEG

Definition at line 121 of file filter.c.

double* moving_average ( double *  s,
int  n,
int  win 
)

Simple Moving Average filter in the time-domain. s[i] = mean{ s[i-win], ..., s[i+win] }

Definition at line 94 of file filter.c.

double* running_median ( double *  d,
int  n,
int  win 
)

Simple Running Median filter.

 d[i] = median{ d[i-win], ..., d[i+win] } 

-- choice of time-window 'win' is crucial use 2*l+1 where l is the size of to-be-preserved details in sampling units Ref: Schelter et al, 2005, chap. 6 (6.2.1)

Definition at line 60 of file filter.c.

double* weighted_running_median ( double *  d,
int  n,
int  win 
)

Weighted Running Median filter. Ref: Schelter et al, 2005, chap. 6 (6.2.1)

Definition at line 143 of file filter.c.