connectivityanalysis#

syncopy.connectivityanalysis(data, method='coh', keeptrials=False, output='abs', foi=None, foilim=None, pad='maxperlen', polyremoval=0, tapsmofrq=None, nTaper=None, taper='hann', taper_opt=None, parallel=None, select=None, **kwargs)[source]#

Perform connectivity analysis of Syncopy SpectralData OR directly AnalogData objects

In case the input is an AnalogData object, a (multi-)tapered Fourier analysis is performed implicitly to arrive at the cross spectral densities needed for coherence and Granger causality estimates. Relevant parameters are the same as for freqanalysis() with method='mtmfft':

(‘foi’, ‘foilim’, ‘pad’, ‘tapsmofrq’, ‘nTaper’, ‘taper’, ‘taper_opt’, ‘polyremoval’)

If the input is already in the spectral domain, so data is of class SpectralData, no additional modification of the spectra is performed, and all parameters above to control the spectral estimation have no effect.

Usage Summary

List of available analysis methods and respective distinct options:

“coh”(Multi-) tapered coherency estimate

Compute the normalized cross spectral densities between all channel combinations

  • output : one of (‘abs’, ‘pow’, ‘complex’, ‘angle’, ‘imag’ or ‘real’)

Spectral analysis (input is AnalogData):

  • taper : one of availableTapers

  • tapsmofrq : spectral smoothing box for slepian tapers (in Hz)

  • nTaper : (optional) number of orthogonal tapers for slepian tapers

  • pad: either pad to an absolute length in seconds or set to ‘nextpow2’

“csd”(‘Multi-) tapered cross spectral density estimate

Computes the cross spectral estimates between all channel combinations

output : complex spectrum

Spectral analysis (input is AnalogData):

  • taper : one of availableTapers

  • tapsmofrq : spectral smoothing box for slepian tapers (in Hz)

  • nTaper : (optional) number of orthogonal tapers for slepian tapers

  • pad: either pad to an absolute length in seconds or set to ‘nextpow2’

“corr”Cross-correlations

Computes the one sided (positive lags) cross-correlations between all channel combinations of AnalogData. The maximal lag is half the trial lengths.

  • keeptrials : set to True for single trial cross-correlations

“granger”Spectral Granger-Geweke causality

Computes linear causality estimates between all channel combinations.

WARNING: When inputting SpectralData directly, it is very important that the previous spy.freqanalysis was done without foi/foilim specification as Granger causality needs all attainable frequencies (0, f_Nyquist)!

Spectral analysis (input is AnalogData):

  • taper : one of availableTapers

  • tapsmofrq : spectral smoothing box for slepian tapers (in Hz)

  • nTaper : (optional, not recommended) number of slepian tapers

  • pad: either pad to an absolute length in seconds or set to ‘nextpow2’

After the computation, information about the convergence and potential regularization of the cross-spectral densities can be obtained by inspecting the .info property of the resulting :class:~`syncopy.CrossSpectralData` object. Keys of that info-dict are: `{'converged', 'max rel. err', 'reg. factor', 'initial cond. num'}`.

The parameters listed below can be provided as is or a via a cfg configuration ‘structure’, see Notes for details.

Parameters
  • data (~syncopy.AnalogData) – A non-empty Syncopy SpectralData or AnalogData object

  • method (str) – Connectivity estimation method, one of 'coh'`, 'corr', 'granger', 'csd'

  • output (str) – Relevant for cross-spectral density estimation (method='coh') Use 'pow' for absolute squared coherence, 'abs' for absolute value of coherence , 'complex' for the complex valued coherency or 'angle', 'imag' or 'real' to extract the phase difference, imaginary or real part of the coherency respectively.

  • keeptrials (bool) – Relevant for cross-correlations (method='corr'). If True single-trial cross-correlations are returned.

  • foi (array-like or None) – Frequencies of interest (Hz) for output. If desired frequencies cannot be matched exactly, the closest possible frequencies are used. If foi is None or foi = "all", all attainable frequencies (i.e., zero to Nyquist / 2) are selected.

  • foilim (array-like [fmin, fmax] or None or "all") – Frequency-window [fmin, fmax] (in Hz) of interest. The foi array will be constructed in 1Hz steps from fmin to fmax (inclusive).

  • pad ('maxperlen', float or 'nextpow2' -) – For the default 'maxperlen', no padding is performed in case of equal length trials, while trials of varying lengths are padded to match the longest trial. If pad is a number all trials are padded so that pad indicates the absolute length of all trials after padding (in seconds). For instance pad = 2 pads all trials to an absolute length of 2000 samples, if and only if the longest trial contains at maximum 2000 samples and the samplerate is 1kHz. If pad is 'nextpow2' all trials are padded to the nearest power of two (in samples) of the longest trial.

  • tapsmofrq (float or None) – Only valid if method is 'coh' or 'granger'. Enables multi-tapering and sets the amount of spectral smoothing with slepian tapers in Hz.

  • nTaper (int or None) – Only valid if method is 'coh' or 'granger' and tapsmofrq is set. Number of orthogonal tapers to use for multi-tapering. It is not recommended to set the number of tapers manually! Leave at None for the optimal number to be set automatically.

  • taper (str or None, optional) – Only valid if method is 'coh' or 'granger'. Windowing function, one of availableTapers For multi-tapering with slepian tapers use tapsmofrq directly.

  • taper_opt (dict or None) – Dictionary with keys for additional taper parameters. For example kaiser() has the additional parameter ‘beta’. For multi-tapering set tapsmofrq directly.

  • parallel (None or bool) – If None (recommended), processing is automatically performed in parallel (i.e., concurrently across trials/channel-groups), provided a dask parallel processing client is running and available. Parallel processing can be manually disabled by setting parallel to False. If parallel is True but no parallel processing client is running, computing will be performed sequentially.

  • select (dict or StructDict or str) – In-place selection of subset of input data for processing. Please refer to syncopy.selectdata() for further usage details.

Returns

out – The analyis result with dims [‘time’, ‘freq’, ‘channel_i’, channel_j’]. The out may contain additional metadata, based on the method used to compute it:

  • For method=’granger’, the out.info property contains the keys listed and explained in metadata_keys.

Return type

~syncopy.CrossSpectralData

Examples

In the following adata is an instance of AnalogData

Calculate the coherence between all channels with 2Hz spectral smoothing, and plot the results for two combinations between 30Hz and 90Hz:

>>> coh = spy.connectivityanalysis(adata, method='coh', tapsmofrq=2)
>>> coh.singlepanelplot(channel_i=0, channel_j=1, foilim=[30,90])
>>> coh.singlepanelplot(channel_i=1, channel_j=2, foilim=[30,90])

Compute the cross-correlation between channel 8 and 12 and plot the results for the first 200ms:

>>> cfg = spy.StructDict()
>>> cfg.method = 'corr'
>>> cfg.select = {'channel': ['channel8', 'channel12']}
>>> corr = spy.connectivityanalysis(adata, cfg)
>>> corr.singlepanelplot(channel_i='channel8', channel_j='channel12', toilim=[0, 0.2])

Estimate Granger causality between the same channels (re-using the cfg from above):

>>> cfg.method = 'granger'
>>> granger = spy.connectivityanalysis(adata, cfg)

Plot the results between 15Hz and 30Hz:

>>> granger.singlepanelplot(channel_i='channel8', channel_j='channel12', foilim=[15, 25])