Support routines

Abstract Classes

class sails.modelfit.AbstractLinearModel[source]

This is a base class defining data storage and generic methods for LinearModel classes. New LinearModels should inherit from AbstractLinearModel and overload the fit_model method

companion

Parameter matrix in companion form

compute_diagnostics(data, compute_pc=False)[source]

Compute several diagnostic metrics from a fitted MVAR model and a dataset.

Parameters:
  • data (ndarrat) – The data to compute the model diagnostics with, typically the same data used during model fit.
  • compute_pc (bool) – Flag indicating whether the compute the percent consistency, this is typically the the longest computation of all the metrics here (Default value = False)
Returns:

object containing LL, AIC, BIC, stability ratio, durbin-watson, R squared and percent consistency measures

Return type:

sails.ModelDiagnostics

get_residuals(data, mode='valid')[source]

Returns the prediction error from a fitted model. This is a wrapper function for get_residuals()

Parameters:data (ndarray) – Data to compute the residuals of, of size [nchannels x nsamples x nrealisations]
Returns:Residual data
Return type:ndarray
nsignals

Number of signals in fitted model

order

Order of fitted model

class sails.mvar_metrics.AbstractMVARMetrics[source]

This is a base class for computing various connectivity metrics from fitted autoregressive models

PSD

Power Spectral Density matrix

S

Spectral Matrix of shape [nsignals, nsignals, order, epochs]

coherency

Complex-valued coherency

d_directed_transfer_function

Direct directed transfer function

directed_transfer_function

Directed transfer function

ff_directed_transfer_function

Full-frequency directed transfer function

get_node_weight(metric='S')[source]

Compute the node weight for a given connectivity metric.

Parameters:metric (str) – Name of the connectivity metric to plot (Default value = ‘S’)
Returns:Node-weight matrix of size [nchannels x nchannels]
Return type:ndarray
get_spectral_generalisation(metric='S')[source]

Compute the spectral generalisation matrix for a given metric. Returns the spatial similarity of different frequencies.

Parameters:metric (str) – Name of the connectivity metric to plot (Default value = ‘S’)
Returns:Matrix of spectral generalisation [nfrequencies x nfrequencies]
Return type:ndarray
geweke_granger_causality

Geweke granger causality

imaginary_coherence

Imaginary part of complex coherency

inv_S

Inverse Power Spectrl Density matrix

isolated_effective_coherence

Isolated effective coherence

magnitude_squared_coherence

Magnitude squared coherence

partial_coherence

Partialled magnitude-squared coherence

partial_directed_coherence

Partial directed coherence

phase_coherence

Phase angle of complex coherency

plot_diags(metric='S', inds=None, F=None, ax=None)[source]

Helper function for plotting the diagonal part of a connectivty metric. Calls sails.plotting.plot_diagonal.

Parameters:
  • metric (str) – Name of the connectivity metric to plot (Default value = ‘S’)
  • F (matplotlib figure handle) – Handle for matplotlib figure to plot in (Default value = None)
  • ax (matplotlib axes handle) – Handle for matplotlib axes to plot in (Default value = None)
plot_summary(metric='S', ind=1)[source]

Helper function for plotting a summary of a connectivty metric. Calls sails.plotting.plot_summary.

Parameters:
  • metric (str) – Name of the connectivity metric to plot (Default value = ‘S’)
  • ind (int) – Index into frequency dimension to plot connectivity matrix (Default value = 1)

Statistical Functions

sails.stats.find_cov_multitrial(a, b, ddof=1)[source]

Estimate the average covariance between two datasets across many trials

Parameters:
  • a (ndarray) – a multivariate time series of shape nsignals x nsamples x ntrials
  • b (ndarray) – a multivariate time series of shape nsignals x nsamples x ntrials
  • ddof (int) – if None, no normalisation will be performed. If an integer, covariance will be normalised by (nsamples - ddof) (Default value = 1)
Returns:

A covariance matrix between a and b of shape nsignals x nsignals

Return type:

type

sails.stats.find_cov(a, b, ddof=1)[source]

Estimate the covariance between two datasets.

Parameters:
  • a (ndarray) – a multivariate time series of shape nsignals x nsamples
  • b (ndarray) – a multivariate time series of shape nsignals x nsamples
  • ddof (int) – if None, no normalisation will be performed. If an integer, covariance will be normalised by (nsamples - ddof) (Default value = 1)
Returns:

A covariance matrix between a and b of shape nsignals x nsignals

Return type:

type

sails.stats.mutual_information(A, B, r=None, nbins=21, base=<ufunc 'log2'>)[source]

Calculate the mutual information between two signals.

Parameters:
  • A (ndarray) – First signal, of shape [channels x samples]
  • B (ndarray) – Second signal, of shape [channels x samples]
  • r (tuple (min,max)) – Range of the random variables used for bin calculations. If not set, takes the minimum and maximum values from the combination of A and B (Default value = None)
  • nbins (int) – Number of bins to compute the frequencies within (Default value = 21)
  • base (function) – Base to use for MI calculation - defaults to np.log2
Returns:

the mutual information within each signal passed

Return type:

type

sails.stats.durbin_watson(residuals, step=1)[source]

Test for autocorrelation in the residuals of a model (see [Durbin1950] and [Durbin1951])

The result is a value between 0 and 4. A value of 0 indicates a positive autocorrelaiton, 2 indicates no correlation and a value of 4 indicates a negative autocorrelation.

If many trials are found the statistic is estimated for each trial separately.

Parameters:
  • residuals (ndarray) – Residuals of shape [channels, samples, ntrials] If only two dimensions, will be treated as [channels, samples] If only one dimension, will be treated as [samples]
  • step (int) – step-size to use through residuals when calculating DW. Defaults to 1.
Returns:

durbin-watson index

Return type:

type

sails.stats.percent_consistency(data, residuals)[source]

Estimate how well a model retains the correlational structure of a dataset.

This uses the percent consistency measure outlined in [Ding2000]. All auto and cross correlations within a multivariate data set are computed and compared to the auto and cross correlations of the data as predicted by a model (in this case the data - model residuals, as these are easily obtained).

Parameters:
  • data (ndarray) – An array containing the observed data, of shape [nsignals x nsamples x ntrials]
  • residuals (ndarray) – An array containing the remaining variance of data after model fitting, of shape [nsignals x nsamples x ntrials]
Returns:

A value between 0 and 100 indicating the percentage of the autocorrelation in the data that is retained in the model.

Return type:

type

sails.stats.rsquared_from_residuals(data, residuals, per_signal=False)[source]

Estimate the variance explained by a model from the original data and the model residuals

If data and residuals are 2D, they will be treated as [channels, nsamples]

If data and residuals are 1D, they will be treated as [nsamples]

Parameters:
  • data (ndarray) – Array containing the original data, of shape [channels, nsamples, ntrials]
  • residuals (ndarray) – Array containing the model residuals, of shape [channels, nsamples, ntrials]
  • per_signal (bool) – Boolean indicating whether the variance explained should be computed separately per channel (Default value = False)
Returns:

Value between 0 and 1 indicating the proportion of variance explained by the model.

Return type:

type

sails.stats.mahalanobis(X, Y=None, sigma=None)[source]

Estimate the Mahalanobis distance of each point in an array across the samples dimension.

Parameters:
  • X (ndarray) – Array containing data of points: shape [channels x samples]
  • Y (ndarray (Optional)) – Optional Array containing data of second points: shape [channels x samples]. If not provided, the origin will be used as the point for comparison. (Default value = None)
  • sigma (ndarray (Optional)) – A precomputed covariance matrix to use in the calculation, of shape [channels x channels]. If not provided, sigma will be calculated from X (Default value = None)
Returns:

The Mahalanobis distance for each sample, of shape [samples]

Return type:

type

sails.stats.profile_likelihood(X)[source]

Find the ‘elbow’ point in a curve.

This function uses the automatic dimensionality selection method from [Zhu2006] to estimate the elbow or inflection point in a given curve.

Parameters:X (ndarray) – 1d array containing the curve to be analysed, of shape [samples]
Returns:The log likelihood that each point on the curve separates the distributions to the left and right, of shape [samples]
Return type:type
sails.stats.stability_index(A)[source]

Compute the stability index as defined in [Lutkephol2006] pages 15-16.

This is a proxy for the stationarity assumption of MVAR models, if the magnitude of the largest principal component of the model parameters is less than 1, the model is stable and thus, stationary.

Parameters:A (ndarray) – The MVAR parameter matrix of size [nchannels x nchannels x model_order]
Returns:The stability index
Return type:type
sails.stats.stability_ratio(A)[source]

Compute the stability ratio as a measure of stationarity

A stronger test for stationarity than the stability index. Computes the ratio between the largest eigenvalue of the parameters and all the others. If the largest parameter is larger than all other parameters (indicated by a value < 1) we have a super-stable system.

Parameters:A (ndarray) – The MVAR parameter matrix of size [nchannels x nchannels x model_order]
Returns:The stability ratio
Return type:type
sails.modelvalidation.approx_model_selection_criteria(A, EF, method='ss')[source]

Estimate Akaike’s Information Criterion and Bayes Information Criterion from a model fit using an approximation of the loglikelihood.

Two approximations are available:
‘ss’: The sum square of the residuals. This approximation is outlined in [Burnham2002] and [Chatfield2003]. ‘det’: The determinant of the residual covariance matrix.

The sum squares method is used as default.

Parameters:
  • A (ndarray) – Array containing the parameters of a fitted model, of shape [channels x channels x maxorder]
  • EF (ndarray) – Array containing the residuals of the fitted model, of shape [channels x samples x trials]
  • method ({'ss','det'}) – string indicating which method should be used to approximate the loglikelihood (‘ss’ or ‘det’) (Default value = ‘ss’)
Returns:

tuple of (aLL, AIC, BIC) aLL: Approximation of the loglikelihood (determined by the method defined above) AIC: Estimate of Akaike’s Information Criterion BIC: Estimate of the Bayesian Information Criterion

Return type:

tuple

Data Processing Functions

sails.modal.adjust_phase(x)[source]

Method for normalising eigenvector matrices to meet the conditions in eqn 21 of [Neumaier2001].

Parameters:x (ndarray) – Vector of complex numbers to adjust the phase of
Returns:Vector of numbers with adjusted phase
Return type:type
sails.orthogonalise.symmetric_orthonormal(A, maintain_mag=False)[source]

Implement the [Colclough2015] algorithm for multivariate leakage correction in MEG. Also see the OHBA toolbox (https://github.com/OHBA-analysis/MEG-ROI-nets) and http://empslocal.ex.ac.uk/people/staff/reverson/uploads/Site/procrustes.pdf

Parameters:
  • A (ndarray) – Data matrix to orthogonalise, of size [channels, samples]
  • verbose (bool) – Flag indicating whether to show detailed output (Default value = False)
  • maintain_mag (bool) – Flag indicating whether to maintain magnitudes (Default value = False)
Returns:

  • data (numpy.ndarray) – [channels, samples] of orthogonalised data
  • W (numpy.ndarray)

sails.orthogonalise.innovations(A, order, mvar_class=None, verbose=False)[source]

Implementation of the innovations-orthogonalisation defined in [Pascual-Marqui2017] . This creates the mixing matrix on the input data after variance explained by an MVAR model is removed. The orthogonalisation is defined on the residuals but applied to the raw data.

Parameters:
  • A (ndarray) – Data matrix to orthogonalise, of size [channels, samples]
  • order (int) – Model order used during MVAR model fit
  • mvar_class (sails LinearModel class) – Model class to do model fitting (Default value = None)
  • verbose (bool) – Flag indicating whether to show detailed output (Default value = False)
Returns:

Orthogonalised data matrix, of size [channels, samples]

Return type:

ndarray

Tutorial helper functions

sails.tutorial_utils.generate_pink_roots(power, order=24)[source]

Generates the poles of a system with a pink spectrum according to eqn 116 in [Kasdin1995]

Parameters:
  • power (float) – Power of pole must be between -2 and 2.
  • order – Order of system to generate (Default value = 24)
Returns:

Roots of the system (length of order)

Return type:

ndarray

sails.tutorial_utils.model_from_roots(roots)[source]

Sets up a univariate model based on the roots of a polynomial

Parameters:roots (ndarray) – Vector of polynomial roots.
Returns:Model fit class containing the parameters defined by the input roots
Return type:AbstractLinearModel
sails.tutorial_utils.find_support_path()[source]

Returns the path to the support files directory

sails.tutorial_utils.find_example_path()[source]

Returns the path to the tutorial example data directory (if available)

sails.tutorial_utils.set_sail()[source]