algorithms.group.parcel_analysis

Module: algorithms.group.parcel_analysis

Inheritance diagram for nipy.algorithms.group.parcel_analysis:

digraph inheritancef8e0cf5d6c { rankdir=LR; size="8.0, 12.0"; "group.parcel_analysis.ParcelAnalysis" [URL="#nipy.algorithms.group.parcel_analysis.ParcelAnalysis",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5)",target="_top"]; }

Parcel-based group analysis of multi-subject image data.

Routines implementing Bayesian inference on group-level effects assumed to be constant within given brain parcels. The model accounts for both estimation errors and localization uncertainty in reference space of first-level images.

See:

Keller, Merlin et al (2008). Dealing with Spatial Normalization Errors in fMRI Group Inference using Hierarchical Modeling. Statistica Sinica; 18(4).

Keller, Merlin et al (2009). Anatomically Informed Bayesian Model Selection for fMRI Group Data Analysis. In MICCAI‘09, Lecture Notes in Computer Science; 5762:450–457.

Roche, Alexis (2012). OHBM‘12 talk, slides at: https://sites.google.com/site/alexisroche/slides/Talk_Beijing12.pdf

ParcelAnalysis

class nipy.algorithms.group.parcel_analysis.ParcelAnalysis(con_imgs, parcel_img, parcel_info=None, msk_img=None, vcon_imgs=None, design_matrix=None, cvect=None, fwhm=8, smooth_method='default', res_path=None, write_smoothed_images=False)[source]

Bases: object

__init__(con_imgs, parcel_img, parcel_info=None, msk_img=None, vcon_imgs=None, design_matrix=None, cvect=None, fwhm=8, smooth_method='default', res_path=None, write_smoothed_images=False)[source]

Bayesian parcel-based analysis.

Given a sequence of independent images registered to a common space (for instance, a set of contrast images from a first-level fMRI analysis), perform a second-level analysis assuming constant effects throughout parcels defined from a given label image in reference space. Specifically, a model of the following form is assumed:

Y = X * beta + variability,

where Y denotes the input image sequence, X is a design matrix, and beta are parcel-wise parameter vectors. The algorithm computes the Bayesian posterior probability of beta in each parcel using an expectation propagation scheme.

Parameters

con_imgs: sequence of nipy-like images

Images input to the group analysis.

parcel_img: nipy-like image

Label image where each label codes for a parcel.

parcel_info: sequence of arrays, optional

A sequence of two arrays with same length equal to the number of distinct parcels consistently with the parcel_img argument. The first array gives parcel names and the second, parcel values, i.e., corresponding intensities in the associated parcel image. By default, parcel values are taken as np.unique(parcel_img.get_data()) and parcel names are these values converted to strings.

msk_img: nipy-like image, optional

Binary mask to restrict analysis. By default, analysis is carried out on all parcels with nonzero value.

vcon_imgs: sequece of nipy-like images, optional

First-level variance estimates corresponding to con_imgs. This is useful if the input images are “noisy”. By default, first-level variances are assumed to be zero.

design_matrix: array, optional

If None, a one-sample analysis model is used. Otherwise, an array with shape (n, p) where n matches the number of input scans, and p is the number of regressors.

cvect: array, optional

Contrast vector of interest. The method makes an inference on the contrast defined as the dot product cvect’*beta, where beta are the unknown parcel-wise effects. If None, cvect is assumed to be np.array((1,)). However, the cvect argument is mandatory if design_matrix is provided.

fwhm: float, optional

A parameter that represents the localization uncertainty in reference space in terms of the full width at half maximum of an isotropic Gaussian kernel.

smooth_method: str, optional

One of ‘default’ and ‘spm’. Setting smooth_method=spm results in simply smoothing the input images using a Gaussian kernel, while the default method involves more complex smoothing in order to propagate spatial uncertainty into the inference process.

res_path: str, optional

An existing path to write output images. If None, no output is written.

write_smoothed_images: bool, optional

Specify whether smoothed images computed throughout the inference process are to be written on disk in res_path.

dump_results(path=None)[source]

Save parcel analysis information in NPZ file.

t_map()[source]

Compute voxel-wise t-statistic map. This map is different from what you would get from an SPM-style mass univariate analysis because the method accounts for both spatial uncertainty in reference space and possibly errors on first-level inputs (if variance images are provided).

Returns

tmap_img: nipy image

t-statistic map.

parcel_maps(full_res=True)[source]

Compute parcel-based posterior contrast means and positive contrast probabilities.

Parameters

full_res: boolean

If True, the output images will be at the same resolution as the parcel image. Otherwise, resolution will match the first-level images.

Returns

pmap_mu_img: nipy image

Image of posterior contrast means for each parcel.

pmap_prob_img: nipy image

Corresponding image of posterior probabilities of positive contrast.

nipy.algorithms.group.parcel_analysis.parcel_analysis(con_imgs, parcel_img, msk_img=None, vcon_imgs=None, design_matrix=None, cvect=None, fwhm=8, smooth_method='default', res_path=None)[source]

Helper function for Bayesian parcel-based analysis.

Given a sequence of independent images registered to a common space (for instance, a set of contrast images from a first-level fMRI analysis), perform a second-level analysis assuming constant effects throughout parcels defined from a given label image in reference space. Specifically, a model of the following form is assumed:

Y = X * beta + variability,

where Y denotes the input image sequence, X is a design matrix, and beta are parcel-wise parameter vectors. The algorithm computes the Bayesian posterior probability of cvect’*beta, where cvect is a given contrast vector, in each parcel using an expectation propagation scheme.

Parameters

con_imgs: sequence of nipy-like images

Images input to the group analysis.

parcel_img: nipy-like image

Label image where each label codes for a parcel.

msk_img: nipy-like image, optional

Binary mask to restrict analysis. By default, analysis is carried out on all parcels with nonzero value.

vcon_imgs: sequece of nipy-like images, optional

First-level variance estimates corresponding to con_imgs. This is useful if the input images are “noisy”. By default, first-level variances are assumed to be zero.

design_matrix: array, optional

If None, a one-sample analysis model is used. Otherwise, an array with shape (n, p) where n matches the number of input scans, and p is the number of regressors.

cvect: array, optional

Contrast vector of interest. The method makes an inference on the contrast defined as the dot product cvect’*beta, where beta are the unknown parcel-wise effects. If None, cvect is assumed to be np.array((1,)). However, the cvect argument is mandatory if design_matrix is provided.

fwhm: float, optional

A parameter that represents the localization uncertainty in reference space in terms of the full width at half maximum of an isotropic Gaussian kernel.

smooth_method: str, optional

One of ‘default’ and ‘spm’. Setting smooth_method=spm results in simply smoothing the input images using a Gaussian kernel, while the default method involves more complex smoothing in order to propagate spatial uncertainty into the inference process.

res_path: str, optional

An existing path to write output images. If None, no output is written.

Returns

pmap_mu_img: nipy image

Image of posterior contrast means for each parcel.

pmap_prob_img: nipy image

Corresponding image of posterior probabilities of positive contrast.