labs.datasets.volumes.volume_img

Module: labs.datasets.volumes.volume_img

Inheritance diagram for nipy.labs.datasets.volumes.volume_img:

Inheritance diagram of nipy.labs.datasets.volumes.volume_img

An image that stores the data as an (x, y, z, …) array, with an affine mapping to the world space

VolumeImg

class nipy.labs.datasets.volumes.volume_img.VolumeImg(data, affine, world_space, metadata=None, interpolation='continuous')

Bases: VolumeGrid

A regularly-spaced image for embedding data in an x, y, z 3D world, for neuroimaging.

This object is an ndarray representing a volume, with the first 3 dimensions being spatial, and mapped to a named world space using an affine (4x4 matrix).

Notes

The data is stored in an undefined way: prescalings might need to be applied to it before using it, or the data might be loaded on demand. The best practice to access the data is not to access the _data attribute, but to use the get_fdata method.

Attributes:
affine4x4 ndarray

Affine mapping from indices to world coordinates.

world_spacestring

Name of the world space the data is embedded in. For instance mni152.

metadatadictionary

Optional, user-defined, dictionary used to carry around extra information about the data as it goes through transformations. The consistency of this information may not be maintained as the data is modified.

interpolation‘continuous’ or ‘nearest’

String giving the interpolation logic used when calculating values in different world spaces

_data

Private pointer to the data.

__init__(data, affine, world_space, metadata=None, interpolation='continuous')

Creates a new neuroimaging image with an affine mapping.

Parameters:
datandarray

ndarray representing the data.

affine4x4 ndarray

affine transformation to the reference world space

world_spacestring

name of the reference world space.

metadatadictionary

dictionary of user-specified information to store with the image.

affine = array([[1., 0., 0., 0.],        [0., 1., 0., 0.],        [0., 0., 1., 0.],        [0., 0., 0., 1.]])
as_volume_img(affine=None, shape=None, interpolation=None, copy=True)

Resample the image to be an image with the data points lying on a regular grid with an affine mapping to the word space (a nipy VolumeImg).

Parameters:
affine: 4x4 or 3x3 ndarray, optional

Affine of the new voxel grid or transform object pointing to the new voxel coordinate grid. If a 3x3 ndarray is given, it is considered to be the rotation part of the affine, and the best possible bounding box is calculated, in this case, the shape argument is not used. If None is given, a default affine is provided by the image.

shape: (n_x, n_y, n_z), tuple of integers, optional

The shape of the grid used for sampling, if None is given, a default affine is provided by the image.

interpolationNone, ‘continuous’ or ‘nearest’, optional

Interpolation type used when calculating values in different word spaces. If None, the image’s interpolation logic is used.

Returns:
resampled_imagenipy VolumeImg

New nipy VolumeImg with the data sampled on the grid defined by the affine and shape.

Notes

The coordinate system of the image is not changed: the returned image points to the same world space.

composed_with_transform(w2w_transform)

Return a new image embedding the same data in a different word space using the given world to world transform.

Parameters:
w2w_transformtransform object

The transform object giving the mapping between the current world space of the image, and the new word space.

Returns:
remapped_imagenipy image

An image containing the same data, expressed in the new world space.

get_affine()
get_fdata()

Return data as a numpy array.

get_transform()

Returns the transform object associated with the volumetric structure which is a general description of the mapping from the values to the world space.

Returns:
transformnipy.datasets.Transform object
get_world_coords()

Return the data points coordinates in the world space.

Returns:
x: ndarray

x coordinates of the data points in world space

y: ndarray

y coordinates of the data points in world space

z: ndarray

z coordinates of the data points in world space

interpolation = 'continuous'
like_from_data(data)

Returns an volumetric data structure with the same relationship between data and world space, and same metadata, but different data.

Parameters:
data: ndarray
metadata = {}
resampled_to_img(target_image, interpolation=None)

Resample the data to be on the same voxel grid than the target volume structure.

Parameters:
target_imagenipy image

Nipy image onto the voxel grid of which the data will be resampled. This can be any kind of img understood by Nipy (datasets, pynifti objects, nibabel object) or a string giving the path to a nifti of analyse image.

interpolationNone, ‘continuous’ or ‘nearest’, optional

Interpolation type used when calculating values in different word spaces. If None, the image’s interpolation logic is used.

Returns:
resampled_imagenipy_image

New nipy image with the data resampled.

Notes

Both the target image and the original image should be embedded in the same world space.

values_in_world(x, y, z, interpolation=None)

Return the values of the data at the world-space positions given by x, y, z

Parameters:
xnumber or ndarray

x positions in world space, in other words millimeters

ynumber or ndarray

y positions in world space, in other words millimeters. The shape of y should match the shape of x

znumber or ndarray

z positions in world space, in other words millimeters. The shape of z should match the shape of x

interpolationNone, ‘continuous’ or ‘nearest’, optional

Interpolation type used when calculating values in different word spaces. If None, the image’s interpolation logic is used.

Returns:
valuesnumber or ndarray

Data values interpolated at the given world position. This is a number or an ndarray, depending on the shape of the input coordinate.

world_space = ''
xyz_ordered(resample=False, copy=True)

Returns an image with the affine diagonal and positive in the world space it is embedded in.

Parameters:
resample: boolean, optional

If resample is False, no resampling is performed, the axis are only permuted. If it is impossible to get xyz ordering by permuting the axis, a ‘CompositionError’ is raised.

copy: boolean, optional

If copy is True, a deep copy of the image (including the data) is made.