analysis.spectral

Module: analysis.spectral

Inheritance diagram for nitime.analysis.spectral:

Inheritance diagram of nitime.analysis.spectral

Classes

FilterAnalyzer

class nitime.analysis.spectral.FilterAnalyzer(time_series, lb=0, ub=None, boxcar_iterations=2, filt_order=64, gpass=1, gstop=60, iir_ftype='ellip', fir_win='hamming')

Bases: nitime.descriptors.ResetMixin

A class for performing filtering operations on time-series and producing the filtered versions of the time-series

Parameters:

time_series: A nitime TimeSeries object. :

lb,ub: float (optional) :

Lower and upper band of a pass-band into which the data will be filtered. Default: 0, Nyquist

boxcar_iterations: int (optional) :

For box-car filtering, how many times to iterate over the data while convolving with a box-car function. Default: 2

gpass: float (optional) :

For iir filtering, the pass-band maximal ripple loss (default: 1)

gstop: float (optional) :

For iir filtering, the stop-band minimal attenuation (default: 60).

filt_order: int (optional) :

For iir/fir filtering, the order of the filter. Note for fir filtering, this needs to be an even number. Default: 64

iir_ftype: str (optional) :

The type of filter to be used in iir filtering (see scipy.signal.iirdesign for details). Default ‘ellip’

fir_win: str :

The window to be used in fir filtering (see scipy.signal.firwin for details). Default: ‘hamming’

__init__(time_series, lb=0, ub=None, boxcar_iterations=2, filt_order=64, gpass=1, gstop=60, iir_ftype='ellip', fir_win='hamming')

Initialize self. See help(type(self)) for accurate signature.

filtered_boxcar()

Filter the time-series by a boxcar filter.

The low pass filter is implemented by convolving with a boxcar function of the right length and amplitude and the high-pass filter is implemented by subtracting a low-pass version (as above) from the signal

filtered_fourier()

Filter the time-series by passing it to the Fourier domain and null out the frequency bands outside of the range [lb,ub]

filtfilt(b, a, in_ts=None)

Zero-phase delay filtering (either iir or fir).

Parameters:

a,b: filter coefficients :

in_ts: time-series object. :

This allows to replace the input. Instead of analyzing this analyzers input data, analyze some other time-series object

fir()

Filter the time-series using an FIR digital filter. Filtering is done back and forth (using scipy.signal.filtfilt) to achieve zero phase delay

iir()

Filter the time-series using an IIR filter. Filtering is done back and forth (using scipy.signal.filtfilt) to achieve zero phase delay

HilbertAnalyzer

class nitime.analysis.spectral.HilbertAnalyzer(input=None)

Bases: nitime.analysis.base.BaseAnalyzer

Analyzer class for extracting the Hilbert transform

__init__(input=None)

Constructor function for the Hilbert analyzer class.

Parameters:input: TimeSeries :
amplitude()
analytic()

The natural output for this analyzer is the analytic signal

imag()
phase()
real()

MorletWaveletAnalyzer

class nitime.analysis.spectral.MorletWaveletAnalyzer(input=None, freqs=None, sd_rel=0.2, sd=None, f_min=None, f_max=None, nfreqs=None, log_spacing=False, log_morlet=False)

Bases: nitime.analysis.base.BaseAnalyzer

Analyzer class for extracting the (complex) Morlet wavelet transform

__init__(input=None, freqs=None, sd_rel=0.2, sd=None, f_min=None, f_max=None, nfreqs=None, log_spacing=False, log_morlet=False)

Constructor function for the Wavelet analyzer class.

Parameters:

freqs: list or float :

List of center frequencies for the wavelet transform, or a scalar for a single band-passed signal.

sd: list or float :

List of filter bandwidths, given as standard-deviation of center frequencies. Alternatively sd_rel can be specified.

sd_rel: float :

Filter bandwidth, given as a fraction of the center frequencies.

f_min: float :

Minimal frequency.

f_max: float :

Maximal frequency.

nfreqs: int :

Number of frequencies.

log_spacing: bool :

If true, frequencies will be evenly spaced on a log-scale. Default: False

log_morlet: bool :

If True, a log-Morlet wavelet is used, if False, a regular Morlet wavelet is used. Default: False

amplitude()
analytic()

The natural output for this analyzer is the analytic signal

imag()
phase()
real()

SpectralAnalyzer

class nitime.analysis.spectral.SpectralAnalyzer(input=None, method=None, BW=None, adaptive=False, low_bias=False)

Bases: nitime.analysis.base.BaseAnalyzer

Analyzer object for spectral analysis

__init__(input=None, method=None, BW=None, adaptive=False, low_bias=False)

The initialization of the

Parameters:

input: time-series objects :

method: dict (optional), :

The method spec used in calculating ‘psd’ see algorithms.get_spectra() for details.

BW: float (optional), :

In ‘spectrum_multi_taper’ The bandwidth of the windowing function will determine the number tapers to use. This parameters represents trade-off between frequency resolution (lower main lobe BW for the taper) and variance reduction (higher BW and number of averaged estimates).

adaptive : {True/False}

In ‘spectrum_multi_taper’, use an adaptive weighting routine to combine the PSD estimates of different tapers.

low_bias: {True/False} :

In spectrum_multi_taper, use bias correction

Examples

>>> np.set_printoptions(precision=4)  # for doctesting
>>> t1 = ts.TimeSeries(data = np.arange(0,1024,1).reshape(2,512),
... sampling_rate=np.pi)
>>> s1 = SpectralAnalyzer(t1)
>>> s1.method['this_method']
'welch'
>>> s1.method['Fs'] # doctest: +ELLIPSIS
3.1415926535... Hz
>>> f,s = s1.psd
>>> f
array([ 0.    ,  0.0491,  0.0982,  0.1473,  0.1963,  0.2454,  0.2945,
        0.3436,  0.3927,  0.4418,  0.4909,  0.54  ,  0.589 ,  0.6381,
        0.6872,  0.7363,  0.7854,  0.8345,  0.8836,  0.9327,  0.9817,
        1.0308,  1.0799,  1.129 ,  1.1781,  1.2272,  1.2763,  1.3254,
        1.3744,  1.4235,  1.4726,  1.5217,  1.5708])
>>> s[0,0]   # doctest: +ELLIPSIS
1128276.92538360...
cpsd()

This outputs both the PSD and the CSD calculated using algorithms.get_spectra().

Returns:

(f,s): tuple :

f: Frequency bands over which the psd/csd are calculated and s: the n by n by len(f) matrix of PSD (on the main diagonal) and CSD (off diagonal)

periodogram()

This is the spectrum estimated as the FFT of the time-series

Returns:

(f,spectrum): f is an array with the frequencies and spectrum is the :

complex-valued FFT. :

psd()

The standard output for this analyzer is a tuple f,s, where: f is the frequency bands associated with the discrete spectral components and s is the PSD calculated using mlab.psd().

spectrum_fourier()

This is the spectrum estimated as the FFT of the time-series

Returns:

(f,spectrum): f is an array with the frequencies and spectrum is the :

complex-valued FFT. :

spectrum_multi_taper()

The spectrum and cross-spectra, computed using :func:`multi_taper_csd’