NIPY logo

Site Navigation

NIPY Community

Table Of Contents

Previous topic

algorithms

Next topic

descriptors

analysis

Module: analysis

Inheritance diagram for nitime.analysis:

System Message: WARNING/2 (/tmp/buildd/nitime-0.2/doc/api/generated/nitime.analysis.rst, line None)

Could not execute ‘dot’. Are you sure you have ‘graphviz’ installed?

Nitime analysis

This module implements an analysis interface between between time-series objects implemented in the timeseries module and the algorithms provided in the algorithms library and other algorithms.

The general pattern of use of Analyzer objects is that they an object is initialized with a TimeSeries object as input. Depending on the analysis methods implemented in the particular analysis object, additional inputs may also be required.

The methods of the object are then implemented as instances of OneTimeProperty, which means that they are only calculated when they are needed and then cached for further use.

Analyzer objects are generally implemented inheriting the desc.ResetMixin(), which means that they have a reset() method. This method resets the object to its initialized state, in which none of the OneTimeProperty methods have been calculated. This allows to change parameter settings of the object and recalculating the quantities in these methods with the new parameter setting.

Classes

BaseAnalyzer

class nitime.analysis.BaseAnalyzer(input=None)

Bases: nitime.descriptors.ResetMixin

Analyzer that implements the default data flow.

All analyzers inherit from this class at least have to * implement a __init__ function to set parameters * define the ‘output’ property

__init__(input=None)
static parameterlist()
parameters
set_input(input)

Set the input of the analyzer, if you want to reuse the analyzer with a different input than the original

CoherenceAnalyzer

class nitime.analysis.CoherenceAnalyzer(input=None, method=None, unwrap_phases=False)

Bases: nitime.analysis.BaseAnalyzer

Analyzer object for coherence/coherency analysis

__init__(input=None, method=None, unwrap_phases=False)
Parameters :

input: TimeSeries object :

Containing the data to analyze.

method: dict, optional, :

This is the method used for spectral analysis of the signal for the coherence caclulation. See algorithms.get_spectra() documentation for details.

unwrap_phases: bool, optional :

Whether to unwrap the phases. This should be True if you assume that the time-delay is the same for all the frequency bands. See _[Sun2005] for details. Default : False

Examples

>>> t1 = ts.TimeSeries(data = np.arange(0,1024,1).reshape(2,512),sampling_rate=np.pi)
>>> c1 = CoherenceAnalyzer(t1)
>>> c1.method['Fs']
3.14159265359 Hz
>>> c1.method['this_method']
'welch'
>>> c1.coherency[0,1]
array([ 0.9499+0.j    ,  0.9495-0.0332j,  0.8696-0.4571j,  0.8918-0.3848j,
        0.9099-0.3191j,  0.9217-0.2679j,  0.9294-0.2285j,  0.9346-0.1977j,
        0.9382-0.1732j,  0.9407-0.1533j,  0.9426-0.1367j,  0.9441-0.1226j,
        0.9452-0.1106j,  0.9460-0.1001j,  0.9467-0.0908j,  0.9473-0.0826j,
        0.9477-0.0752j,  0.9481-0.0684j,  0.9484-0.0623j,  0.9487-0.0566j,
        0.9489-0.0513j,  0.9491-0.0465j,  0.9492-0.0419j,  0.9494-0.0376j,
        0.9495-0.0336j,  0.9496-0.0298j,  0.9497-0.0263j,  0.9497-0.0232j,
        0.9498-0.0205j,  0.9498-0.0187j,  0.9498-0.0188j,  0.9497-0.0266j,
        1.0000+0.j    ])
>>> c1.phase[0,1]
array([ 0.    , -0.035 , -0.4839, -0.4073, -0.3373, -0.2828, -0.241 ,
       -0.2085, -0.1826, -0.1615, -0.144 , -0.1292, -0.1164, -0.1054,
       -0.0956, -0.0869, -0.0791, -0.072 , -0.0656, -0.0596, -0.0541,
       -0.0489, -0.0441, -0.0396, -0.0353, -0.0314, -0.0277, -0.0244,
       -0.0216, -0.0197, -0.0198, -0.028 ,  0.    ])
static coherence()

The coherence between the different channels in the input TimeSeries object

static coherence_partial()

The partial coherence between data[i] and data[j], given data[k], as a function of frequency band

static coherency()

The standard output for this kind of analyzer is the coherency

static delay()

The delay in seconds between the two time series

static frequencies()

The central frequencies in the bands

static phase()

The frequency-dependent phase relationship between all the pairwise combinations of time-series in the data

static spectrum()

The spectra of each of the channels and cross-spectra between different channles in the input TimeSeries object

CorrelationAnalyzer

class nitime.analysis.CorrelationAnalyzer(input=None)

Bases: nitime.analysis.BaseAnalyzer

Analyzer object for correlation analysis. Has the same API as the CoherenceAnalyzer

__init__(input=None)
Parameters :

input: TimeSeries object :

Containing the data to analyze.

Examples

>>> t1 = ts.TimeSeries(data = np.sin(np.arange(0,10*np.pi,10*np.pi/100)).reshape(2,50),sampling_rate=np.pi)
>>> c1 = CorrelationAnalyzer(t1)
>>> c1 = CorrelationAnalyzer(t1)
>>> c1.corrcoef
array([[ 1., -1.],
       [-1.,  1.]])
>>> c1.xcorr.sampling_rate
3.1415926536 Hz
>>> c1.xcorr.t0
-15.915494309150001 s
static corrcoef()

The correlation coefficient between every pairwise combination of time-series contained in the object

static xcorr()

The cross-correlation between every pairwise combination time-series in the object. Uses np.correlation(‘full’).

Returns :

TimeSeries: the time-dependent cross-correlation, with zero-lag :

at time=0 :

static xcorr_norm()

The cross-correlation between every pairwise combination time-series in the object, where the zero lag correlation is normalized to be equal to the correlation coefficient between the time-series

Returns :

TimeSeries: A TimeSeries object :

the time-dependent cross-correlation, with zero-lag at time=0

EventRelatedAnalyzer

class nitime.analysis.EventRelatedAnalyzer(time_series, events, len_et, zscore=False, correct_baseline=False, offset=0)

Bases: nitime.descriptors.ResetMixin

Analyzer object for reverse-correlation/event-related analysis.

Note: right now, this class assumes the input time series is only two-dimensional. If your input data is something like (nchannels,nsubjects, ...) with more dimensions, things are likely to break in hard to understand ways.

__init__(time_series, events, len_et, zscore=False, correct_baseline=False, offset=0)
Parameters :

time_series: a time-series object :

A time-series with data on which the event-related analysis proceeds

events_time_series: a TimeSeries object or an Events object :

The events which occured in tandem with the time-series in the :

EventRelatedAnalyzer. This object’s data has to have the same :

dimensions as the data in the EventRelatedAnalyzer object. In each :

sample in the time-series, there is an integer, which denotes the kind :

of event which occured at that time. In time-bins in which :

no event occured, a 0 should be entered. The data in this time series :

object needs to have the same dimensionality as the data in the data :

time-series :

len_et: int :

The expected length of the event-triggered quantity (in the same :

time-units as the events are represented (presumably number of TRs, for :

fMRI data). For example, the size of the block dedicated in the :

fir_matrix to each type of event :

zscore: a flag to return the result in zscore (where relevant) :

correct_baseline: a flag to correct the baseline according to the first :

point in the event-triggered average (where possible) :

offset: the offset of the beginning of the event-related time-series, :

relative to the event occurence :

static FIR()

Calculate the FIR event-related estimated of the HRFs for different kinds of events

Returns :

A time-series object, shape[:-2] are dimensions corresponding to the to :

shape[:-2] of the EventRelatedAnalyzer data, shape[-2] corresponds to :

the different kinds of events used (ordered according to the sorted :

order of the unique components in the events time-series). shape[-1] :

corresponds to time, and has length = len_et :

XXX code needs to be changed to use flattening (see ‘eta’ below) :

static FIR_estimate()

Calculate back the LTI estimate of the time-series, from FIR

static et_data()

The event-triggered data (all occurences).

This gets the time-series corresponding to the inidividual event occurences. Returns a list of lists of time-series. The first dimension is the different channels in the original time-series data and the second dimension is each type of event in the event time series

The time-series itself has the first diemnsion of the data being the specific occurence, with time 0 locked to the that occurence of the event and the last dimension is time.e

This complicated structure is so that it can deal with situations where each channel has different events and different events have different # of occurences

static eta()

The event-triggered average activity.

static ets()

The event-triggered standard error of the mean

static xcorr_eta()

Compute the normalized cross-correlation estimate of the HRFs for different kinds of events

Returns :

A time-series object, shape[:-2] are dimensions corresponding to the to :

shape[:-2] of the EventRelatedAnalyzer data, shape[-2] corresponds to :

the different kinds of events used (ordered according to the sorted :

order of the unique components in the events time-series). shape[-1] :

corresponds to time, and has length = len_et (xcorr looks both back :

and forward for half of this length) :

FilterAnalyzer

class nitime.analysis.FilterAnalyzer(time_series, lb=0, ub=None, boxcar_iterations=2)

Bases: nitime.descriptors.ResetMixin

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

__init__(time_series, lb=0, ub=None, boxcar_iterations=2)
static filtered_boxcar()

Filte 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

static 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]

HilbertAnalyzer

class nitime.analysis.HilbertAnalyzer(input=None)

Bases: nitime.analysis.BaseAnalyzer

Analyzer class for extracting the Hilbert transform

__init__(input=None)

Constructor function for the Hilbert analyzer class.

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

The natural output for this analyzer is the analytic signal

static imag()
static phase()
static real()

MTCoherenceAnalyzer

class nitime.analysis.MTCoherenceAnalyzer(input=None, bandwidth=None, alpha=0.050000000000000003, adaptive=True)

Bases: nitime.analysis.BaseAnalyzer

Analyzer for multi-taper coherence analysis, including jack-knife estimate of confidence interval

__init__(input=None, bandwidth=None, alpha=0.050000000000000003, adaptive=True)

Initializer function for the MTCoherenceAnalyzer

Parameters :

input: TimeSeries object :

bandwidth: float, :

The bandwidth of the windowing function will determine the number tapers to use. This parameters represents trade-off between frequency resolution (lower main lobe bandwidth for the taper) and variance reduction (higher bandwidth and number of averaged estimates). Per default will be set to 4 times the fundamental frequency, such that NW=4

alpha: float, default =0.05 :

This is the alpha used to construct a confidence interval around the multi-taper csd estimate, based on a jack-knife estimate of the variance [Thompson2007].

adaptive: bool, default to True :

Whether to set the weights for the tapered spectra according to the adaptive algorithm [Thompson2007].

[Thompson2007]

Thompson, DJ Jackknifing multitaper spectrum

estimates. IEEE Signal Processing Magazing. 24: 20-30

static coherence()
static confidence_interval()

The size of the 1-alpha confidence interval

static df()
static eigs()
static frequencies()
static spectra()
static tapers()
static weights()

MorletWaveletAnalyzer

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

Bases: nitime.analysis.BaseAnalyzer

Analyzer class for extracting the (complex) Morlet wavelet transform

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

Constructor function for the Hilbert 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

static amplitude()
static analytic()

The natural output for this analyzer is the analytic signal

static imag()
static phase()
static real()

NormalizationAnalyzer

class nitime.analysis.NormalizationAnalyzer(input=None)

Bases: nitime.analysis.BaseAnalyzer

A class for performing normalization operations on time-series and producing the renormalized versions of the time-series

__init__(input=None)

Constructor function for the Normalization analyzer class.

Parameters :input: TimeSeries object :
static percent_change()
static z_score()

SNRAnalyzer

class nitime.analysis.SNRAnalyzer(input=None, bandwidth=None, adaptive=False, low_bias=False)

Bases: nitime.analysis.BaseAnalyzer

Calculate SNR for a response to repetitions of the same stimulus, according to [Borst1999] (Figure 2) and [Hsu2004].

[Hsu2004]Hsu A, Borst A and Theunissen, FE (2004) Quantifying

variability in neural responses ans its application for the validation of model predictions. Network: Comput Neural Syst 15:91-109

[Borst1999]Borst A and Theunissen FE (1999) Information theory and

neural coding. Nat Neurosci 2:947-957

__init__(input=None, bandwidth=None, adaptive=False, low_bias=False)

Initializer for the multi_taper_SNR object

Parameters :

input: TimeSeries object :

bandwidth: float, :

The bandwidth of the windowing function will determine the number tapers to use. This parameters represents trade-off between frequency resolution (lower main lobe bandwidth for the taper) and variance reduction (higher bandwidth and number of averaged estimates). Per default will be set to 4 times the fundamental frequency, such that NW=4

adaptive: bool, default to False :

Whether to set the weights for the tapered spectra according to the adaptive algorithm [Thompson2007].

low_bias : bool, default to False

Rather than use 2NW tapers, only use the tapers that have better than 90% spectral concentration within the bandwidth (still using a maximum of 2NW tapers)

[Thompson2007]

Thompson, DJ Jackknifing multitaper spectrum

estimates. IEEE Signal Processing Magazing. 24: 20-30

static correlation()

The correlation between all combinations of trials

Returns :

(r,e) : tuple

r is the mean correlation and e is the mean error of the correlation (with df = n_trials - 1)

static mt_coherence()
static mt_frequencies()
static mt_information()
static mt_noise_psd()
static mt_signal_psd()
static mt_snr()

SeedCoherenceAnalyzer

class nitime.analysis.SeedCoherenceAnalyzer(seed_time_series=None, target_time_series=None, method=None, lb=0, ub=None, prefer_speed_over_memory=True, scale_by_freq=True)

Bases: nitime.analysis.BaseAnalyzer

This analyzer takes two time-series. The first is designated as a time-series of seeds. The other is designated as a time-series of targets. The analyzer performs a coherence analysis between each of the channels in the seed time-series and all of the channels in the target time-series.

__init__(seed_time_series=None, target_time_series=None, method=None, lb=0, ub=None, prefer_speed_over_memory=True, scale_by_freq=True)

The constructor for the SeedCoherenceAnalyzer

Parameters :

seed_time_series: a time-series object :

target_time_series: a time-series object :

lb,ub: float,optional, default: lb=0, ub=None (max frequency) :

define a frequency band of interest

prefer_speed_over_memory: Boolean, optional, default=True :

Makes things go a bit faster, if you have enough memory

static coherence()

The coherence between each of the channels of the seed time series and all the channels of the target time-series.

static coherency()
static delay()

The delay in seconds between the two time series

static frequencies()

Get the central frequencies for the frequency bands, given the method of estimating the spectrum

static relative_phases()

The frequency-band dependent relative phase between the two time-series

static target_cache()

SparseCoherenceAnalyzer

class nitime.analysis.SparseCoherenceAnalyzer(time_series=None, ij=(0, 0), method=None, lb=0, ub=None, prefer_speed_over_memory=True, scale_by_freq=True)

Bases: nitime.analysis.BaseAnalyzer

This analyzer is intended for analysis of large sets of data, in which possibly only a subset of combinations of time-series needs to be compared. The constructor for this class receives as input not only a time-series object, but also a list of tuples with index combinations (i,j) for the combinations. Importantly, this class implements only the mlab csd function and cannot use other methods of spectral estimation

__init__(time_series=None, ij=(0, 0), method=None, lb=0, ub=None, prefer_speed_over_memory=True, scale_by_freq=True)

The constructor for the SparseCoherenceAnalyzer

Parameters :

time_series: a time-series object :

ij: a list of tuples, each containing a pair of indices. :

The resulting cache will contain the fft of time-series in the rows indexed by the unique elements of the union of i and j

lb,ub: float,optional, default: lb=0, ub=None (max frequency) :

define a frequency band of interest

prefer_speed_over_memory: Boolean, optional, default=True :

Does exactly what the name implies. If you have enough memory

method: optional, dict :

The method for spectral estimation (see :func:`algorithms.get_spectra`) :

static cache()

Caches the fft windows required by the other methods of the SparseCoherenceAnalyzer. Calculate only once and reuse

static coherence()

The coherence values for the output

static coherency()

The default behavior is to calculate the cache, extract it and then output the coherency

static delay()

The delay in seconds between the two time series

static frequencies()

Get the central frequencies for the frequency bands, given the method of estimating the spectrum

static phases()

The frequency-band dependent phases of the spectra of each of the time -series i,j in the analyzer

static relative_phases()

The frequency-band dependent relative phase between the two time-series

static spectrum()

get the spectrum for the collection of time-series in this analyzer

SpectralAnalyzer

class nitime.analysis.SpectralAnalyzer(input=None, method=None)

Bases: nitime.analysis.BaseAnalyzer

Analyzer object for spectral analysis

__init__(input=None, method=None)

The initialization of the

Parameters :

input: time-series objects :

method: dict optional, see :func:`algorithms.get_spectra` for :

specification of the spectral analysis method :

Examples

>>> 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']
3.14159265359 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]   
1128276.92538360...
static 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)

static 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. :

static 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().

static 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. :

static spectrum_multi_taper()

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

XXX This method needs to be improved to include a clever way of figuring out how many tapers to generate and a way to extract the estimate of error based on the tapers.

Function

nitime.analysis.signal_noise(response)

Signal and noise as defined in Borst and Theunissen 1999, Figure 2

Parameters :

response: nitime TimeSeries object :

The data here are individual responses of a single unit to the same stimulus, with repetitions being the first dimension and time as the last dimension