brainiak.funcalign package

Functional alignment of volumes from different subjects.

Submodules

brainiak.funcalign.rsrm module

Robust Shared Response Model (RSRM)

The implementation is based on the following publications:

[Turek2017](1, 2) “Capturing Shared and Individual Information in fMRI Data”, J. Turek, C. Ellis, L. Skalaban, N. Turk-Browne, T. Willke under review, 2017.
class brainiak.funcalign.rsrm.RSRM(n_iter=10, features=50, gamma=1.0, rand_seed=0)

Bases: sklearn.base.BaseEstimator, sklearn.base.TransformerMixin

Robust Shared Response Model (RSRM)

Given multi-subject data, factorize it as a shared response R among all subjects, an orthogonal transform W per subject, and an individual (outlying) sparse component S per subject:

\[X_i \approx W_i R + S_i, \forall i=1 \dots N\]

This unsupervised model allows to learn idiosyncratic information for subjects and simultaneously improve the shared response estimation. The model has similar properties to the Shared Response Model (SRM) with the addition of the individual components.

The model is estimated solving the following optimization problem:

\[\min_{W_i, S_i, R}\sum_i \frac{1}{2}\|X_i - W_i R - S_i\|_F^2\]
\[+ \gamma\|S_i\|_1\]
\[s.t. \qquad W_i^TW_i = I \quad \forall i=1 \dots N\]

The solution to this problem is obtained by applying a Block-Coordinate Descent procedure. More details can be found in [Turek2017].

Parameters:
  • n_iter (int, default: 10) – Number of iterations to run the algorithm.
  • features (int, default: 50) – Number of features to compute.
  • gamma (float, default: 1.0) – Regularization parameter for the sparseness of the individual components. Higher values yield sparser individual components.
  • rand_seed (int, default: 0) – Seed for initializing the random number generator.
w_

list of array, element i has shape=[voxels_i, features] – The orthogonal transforms (mappings) for each subject.

r_

array, shape=[features, timepoints] – The shared response.

s_

list of array, element i has shape=[voxels_i, timepoints] – The individual components for each subject.

random_state_

RandomState – Random number generator initialized using rand_seed

Note

The number of voxels may be different between subjects. However, the number of timepoints for the alignment data must be the same across subjects.

The Robust Shared Response Model is approximated using the Block-Coordinate Descent (BCD) algorithm proposed in [Turek2017].

This is a single node version.

fit(X)

Compute the Robust Shared Response Model

Parameters:X (list of 2D arrays, element i has shape=[voxels_i, timepoints]) – Each element in the list contains the fMRI data of one subject.
transform(X)

Use the model to transform new data to Shared Response space

Parameters:X (list of 2D arrays, element i has shape=[voxels_i, timepoints_i]) – Each element in the list contains the fMRI data of one subject.
Returns:
  • r (list of 2D arrays, element i has shape=[features_i, timepoints_i]) – Shared responses from input data (X)
  • s (list of 2D arrays, element i has shape=[voxels_i, timepoints_i]) – Individual data obtained from fitting model to input data (X)
transform_subject(X)

Transform a new subject using the existing model

Parameters:X (2D array, shape=[voxels, timepoints]) – The fMRI data of the new subject.
Returns:
  • w (2D array, shape=[voxels, features]) – Orthogonal mapping W_{new} for new subject
  • s (2D array, shape=[voxels, timepoints]) – Individual term S_{new} for new subject

brainiak.funcalign.srm module

Shared Response Model (SRM)

The implementations are based on the following publications:

[Chen2015](1, 2) “A Reduced-Dimension fMRI Shared Response Model”, P.-H. Chen, J. Chen, Y. Yeshurun-Dishon, U. Hasson, J. Haxby, P. Ramadge Advances in Neural Information Processing Systems (NIPS), 2015. http://papers.nips.cc/paper/5855-a-reduced-dimension-fmri-shared-response-model
[Anderson2016]“Enabling Factor Analysis on Thousand-Subject Neuroimaging Datasets”, Michael J. Anderson, Mihai Capotă, Javier S. Turek, Xia Zhu, Theodore L. Willke, Yida Wang, Po-Hsuan Chen, Jeremy R. Manning, Peter J. Ramadge, Kenneth A. Norman, IEEE International Conference on Big Data, 2016. https://doi.org/10.1109/BigData.2016.7840719
class brainiak.funcalign.srm.SRM(n_iter=10, features=50, rand_seed=0, comm=<mpi4py.MPI.Intracomm object>)

Bases: sklearn.base.BaseEstimator, sklearn.base.TransformerMixin

Probabilistic Shared Response Model (SRM)

Given multi-subject data, factorize it as a shared response S among all subjects and an orthogonal transform W per subject:

\[X_i \approx W_i S, \forall i=1 \dots N\]
Parameters:
  • n_iter (int, default: 10) – Number of iterations to run the algorithm.
  • features (int, default: 50) – Number of features to compute.
  • rand_seed (int, default: 0) – Seed for initializing the random number generator.
  • comm (mpi4py.MPI.Intracomm) – The MPI communicator containing the data
w_

list of array, element i has shape=[voxels_i, features] – The orthogonal transforms (mappings) for each subject.

s_

array, shape=[features, samples] – The shared response.

sigma_s_

array, shape=[features, features] – The covariance of the shared response Normal distribution.

mu_

list of array, element i has shape=[voxels_i] – The voxel means over the samples for each subject.

rho2_

array, shape=[subjects] – The estimated noise variance \(\rho_i^2\) for each subject

comm

mpi4py.MPI.Intracomm – The MPI communicator containing the data

random_state_

RandomState – Random number generator initialized using rand_seed

Note

The number of voxels may be different between subjects. However, the number of samples must be the same across subjects.

The probabilistic Shared Response Model is approximated using the Expectation Maximization (EM) algorithm proposed in [Chen2015]. The implementation follows the optimizations published in [Anderson2016].

This is a single node version.

The run-time complexity is \(O(I (V T K + V K^2 + K^3))\) and the memory complexity is \(O(V T)\) with I - the number of iterations, V - the sum of voxels from all subjects, T - the number of samples, and K - the number of features (typically, \(V \gg T \gg K\)).

fit(X, y=None)

Compute the probabilistic Shared Response Model

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, samples]) – Each element in the list contains the fMRI data of one subject.
  • y (not used) –
transform(X, y=None)

Use the model to transform matrix to Shared Response space

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, samples_i]) – Each element in the list contains the fMRI data of one subject note that number of voxels and samples can vary across subjects
  • y (not used (as it is unsupervised learning)) –
Returns:

s – Shared responses from input data (X)

Return type:

list of 2D arrays, element i has shape=[features_i, samples_i]

transform_subject(X)

Transform a new subject using the existing model. The subject is assumed to have recieved equivalent stimulation

Parameters:X (2D array, shape=[voxels, timepoints]) – The fMRI data of the new subject.
Returns:w – Orthogonal mapping W_{new} for new subject
Return type:2D array, shape=[voxels, features]
class brainiak.funcalign.srm.DetSRM(n_iter=10, features=50, rand_seed=0)

Bases: sklearn.base.BaseEstimator, sklearn.base.TransformerMixin

Deterministic Shared Response Model (DetSRM)

Given multi-subject data, factorize it as a shared response S among all subjects and an orthogonal transform W per subject:

\[X_i \approx W_i S, \forall i=1 \dots N\]
Parameters:
  • n_iter (int, default: 10) – Number of iterations to run the algorithm.
  • features (int, default: 50) – Number of features to compute.
  • rand_seed (int, default: 0) – Seed for initializing the random number generator.
w_

list of array, element i has shape=[voxels_i, features] – The orthogonal transforms (mappings) for each subject.

s_

array, shape=[features, samples] – The shared response.

random_state_

RandomState – Random number generator initialized using rand_seed

Note

The number of voxels may be different between subjects. However, the number of samples must be the same across subjects.

The Deterministic Shared Response Model is approximated using the Block Coordinate Descent (BCD) algorithm proposed in [Chen2015].

This is a single node version.

The run-time complexity is \(O(I (V T K + V K^2))\) and the memory complexity is \(O(V T)\) with I - the number of iterations, V - the sum of voxels from all subjects, T - the number of samples, K - the number of features (typically, \(V \gg T \gg K\)), and N - the number of subjects.

fit(X, y=None)

Compute the Deterministic Shared Response Model

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, samples]) – Each element in the list contains the fMRI data of one subject.
  • y (not used) –
transform(X, y=None)

Use the model to transform data to the Shared Response subspace

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, samples_i]) – Each element in the list contains the fMRI data of one subject.
  • y (not used) –
Returns:

s – Shared responses from input data (X)

Return type:

list of 2D arrays, element i has shape=[features_i, samples_i]

transform_subject(X)

Transform a new subject using the existing model. The subject is assumed to have recieved equivalent stimulation

Parameters:X (2D array, shape=[voxels, timepoints]) – The fMRI data of the new subject.
Returns:w – Orthogonal mapping W_{new} for new subject
Return type:2D array, shape=[voxels, features]

brainiak.funcalign.sssrm module

Semi-Supervised Shared Response Model (SS-SRM)

The implementations are based on the following publications:

[Turek2016](1, 2) “A Semi-Supervised Method for Multi-Subject fMRI Functional Alignment”, J. S. Turek, T. L. Willke, P.-H. Chen, P. J. Ramadge IEEE International Conference on Acoustics, Speech and Signal Processing (ICASSP), 2017, pp. 1098-1102. https://doi.org/10.1109/ICASSP.2017.7952326
class brainiak.funcalign.sssrm.SSSRM(n_iter=10, features=50, gamma=1.0, alpha=0.5, rand_seed=0)

Bases: sklearn.base.BaseEstimator, sklearn.base.ClassifierMixin, sklearn.base.TransformerMixin

Semi-Supervised Shared Response Model (SS-SRM)

Given multi-subject data, factorize it as a shared response S among all subjects and an orthogonal transform W per subject, using also labeled data to train a Multinomial Logistic Regression (MLR) classifier (with l2 regularization) in a semi-supervised manner:

(1)\[(1-\alpha) Loss_{SRM}(W_i,S;X_i) + \alpha/\gamma Loss_{MLR}(\theta, bias; {(W_i^T \times Z_i, y_i}) + R(\theta)\]

(see Equations (1) and (4) in [Turek2016]).

Parameters:
  • n_iter (int, default: 10) – Number of iterations to run the algorithm.
  • features (int, default: 50) – Number of features to compute.
  • gamma (float, default: 1.0) – Regularization parameter for the classifier.
  • alpha (float, default: 0.5) – Balance parameter between the SRM term and the MLR term.
  • rand_seed (int, default: 0) – Seed for initializing the random number generator.
w_

list of array, element i has shape=[voxels_i, features] – The orthogonal transforms (mappings) for each subject.

s_

array, shape=[features, samples] – The shared response.

theta_

array, shape=[classes, features] – The MLR class plane parameters.

bias_

array, shape=[classes] – The MLR class biases.

classes_

array of int, shape=[classes] – Mapping table for each classes to original class label.

random_state_

RandomState – Random number generator initialized using rand_seed

Note

The number of voxels may be different between subjects. However, the number of samples for the alignment data must be the same across subjects. The number of labeled samples per subject can be different.

The Semi-Supervised Shared Response Model is approximated using the Block-Coordinate Descent (BCD) algorithm proposed in [Turek2016].

This is a single node version.

fit(X, y, Z)

Compute the Semi-Supervised Shared Response Model

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, n_align]) – Each element in the list contains the fMRI data for alignment of one subject. There are n_align samples for each subject.
  • y (list of arrays of int, element i has shape=[samples_i]) – Each element in the list contains the labels for the data samples in Z.
  • Z (list of 2D arrays, element i has shape=[voxels_i, samples_i]) – Each element in the list contains the fMRI data of one subject for training the MLR classifier.
predict(X)

Classify the output for given data

Parameters:X (list of 2D arrays, element i has shape=[voxels_i, samples_i]) – Each element in the list contains the fMRI data of one subject The number of voxels should be according to each subject at the moment of training the model.
Returns:p – Predictions for each data sample.
Return type:list of arrays, element i has shape=[samples_i]
transform(X, y=None)

Use the model to transform matrix to Shared Response space

Parameters:
  • X (list of 2D arrays, element i has shape=[voxels_i, samples_i]) – Each element in the list contains the fMRI data of one subject note that number of voxels and samples can vary across subjects.
  • y (not used as it only applies the mappings) –
Returns:

s – Shared responses from input data (X)

Return type:

list of 2D arrays, element i has shape=[features_i, samples_i]