cifti2
¶
CIFTI-2 format IO
Read / write access to CIFTI-2 image format |
|
Defines |
Module: cifti2.cifti2
¶
Read / write access to CIFTI-2 image format
Format of the NIFTI2 container format described here:
Definition of the CIFTI-2 header format and file extensions can be found at:
|
Element representing a mapping of the dimension to vertex or voxels. |
|
Class for CIFTI-2 header extension |
Error in CIFTI-2 header |
|
|
Class for single file CIFTI-2 format image |
|
CIFTI-2 label: association of integer key with a name and RGBA values |
CIFTI-2 label table: a sequence of |
|
CIFTI-2 Matrix object |
|
|
Class for Matrix Indices Map |
|
A list of name-value pairs |
|
CIFTI-2 named map: association of name and optional data with a map index |
|
CIFTI-2 parcel: association of a name with vertices and/or voxels |
|
Cifti surface: association of brain structure and number of vertices |
Matrix that translates voxel indices to spatial coordinates |
|
|
CIFTI-2 vertex indices: vertex indices for an associated brain model |
|
CIFTI-2 vertices - association of brain structure and a list of vertices |
|
CIFTI-2 volume: information about a volume for mappings that use voxels |
|
CIFTI-2 VoxelIndicesIJK: Set of voxel indices contained in a structure |
|
Initialize header from binary data block and extensions |
Module: cifti2.cifti2_axes
¶
Defines Axis
objects to create, read, and manipulate CIFTI-2 files
These axes provide an alternative interface to the information in the CIFTI-2 header. Each type of CIFTI-2 axes describing the rows/columns in a CIFTI-2 matrix is given a unique class:
BrainModelAxis
: each row/column is a voxel or vertexParcelsAxis
: each row/column is a group of voxels and/or verticesScalarAxis
: each row/column has a unique name (with optional meta-data)LabelAxis
: each row/column has a unique name and label table (with optional meta-data)SeriesAxis
: each row/column is a timepoint, which increases monotonically
All of these classes are derived from the Axis class.
After loading a CIFTI-2 file a tuple of axes describing the rows and columns can be obtained
from the cifti2.Cifti2Header.get_axis()
method on the header object
(e.g. nibabel.load(<filename>).header.get_axis()
). Inversely, a new
cifti2.Cifti2Header
object can be created from existing Axis objects
using the cifti2.Cifti2Header.from_axes()
factory method.
CIFTI-2 Axis objects of the same type can be concatenated using the ‘+’-operator. Numpy indexing also works on axes (except for SeriesAxis objects, which have to remain monotonically increasing or decreasing).
Creating new CIFTI-2 axes¶
New Axis objects can be constructed by providing a description for what is contained in each row/column of the described tensor. For each Axis sub-class this descriptor is:
BrainModelAxis
: a CIFTI-2 structure name and a voxel or vertex indexParcelsAxis
: a name and a sequence of voxel and vertex indicesScalarAxis
: a name and optionally a dict of meta-dataLabelAxis
: a name, dict of label index to name and colour, and optionally a dict of meta-dataSeriesAxis
: the time-point of each row/column is set by setting the start, stop, size, and unit of the time-series
Several helper functions exist to create new BrainModelAxis
axes:
BrainModelAxis.from_mask()
creates a new BrainModelAxis volume covering the non-zero values of a maskBrainModelAxis.from_surface()
creates a new BrainModelAxis surface covering the provided indices of a surface
A ParcelsAxis
axis can be created from a sequence of BrainModelAxis
axes using
ParcelsAxis.from_brain_models()
.
Examples¶
We can create brain models covering the left cortex and left thalamus using:
>>> from nibabel import cifti2
>>> import numpy as np
>>> bm_cortex = cifti2.BrainModelAxis.from_mask([True, False, True, True],
... name='cortex_left')
>>> bm_thal = cifti2.BrainModelAxis.from_mask(np.ones((2, 2, 2)), affine=np.eye(4),
... name='thalamus_left')
In this very simple case bm_cortex
describes a left cortical surface skipping the second
out of four vertices. bm_thal
contains all voxels in a 2x2x2 volume.
Brain structure names automatically get converted to valid CIFTI-2 identifiers using
BrainModelAxis.to_cifti_brain_structure_name()
.
A 1-dimensional mask will be automatically interpreted as a surface element and a 3-dimensional
mask as a volume element.
These can be concatenated in a single brain model covering the left cortex and thalamus by simply adding them together
>>> bm_full = bm_cortex + bm_thal
Brain models covering the full HCP grayordinate space can be constructed by adding all the volumetric and surface brain models together like this (or by reading one from an already existing HCP file).
Getting a specific brain region from the full brain model is as simple as:
>>> assert bm_full[bm_full.name == 'CIFTI_STRUCTURE_CORTEX_LEFT'] == bm_cortex
>>> assert bm_full[bm_full.name == 'CIFTI_STRUCTURE_THALAMUS_LEFT'] == bm_thal
You can also iterate over all brain structures in a brain model:
>>> for idx, (name, slc, bm) in enumerate(bm_full.iter_structures()):
... print((str(name), slc))
... assert bm == bm_full[slc]
... assert bm == bm_cortex if idx == 0 else bm_thal
('CIFTI_STRUCTURE_CORTEX_LEFT', slice(0, 3, None))
('CIFTI_STRUCTURE_THALAMUS_LEFT', slice(3, None, None))
In this case there will be two iterations, namely: (‘CIFTI_STRUCTURE_CORTEX_LEFT’, slice(0, <size of cortex mask>), bm_cortex) and (‘CIFTI_STRUCTURE_THALAMUS_LEFT’, slice(<size of cortex mask>, None), bm_thal)
ParcelsAxis can be constructed from selections of these brain models:
>>> parcel = cifti2.ParcelsAxis.from_brain_models([
... ('surface_parcel', bm_cortex[:2]), # contains first 2 cortical vertices
... ('volume_parcel', bm_thal), # contains thalamus
... ('combined_parcel', bm_full[[1, 8, 10]]), # contains selected voxels/vertices
... ])
Time series are represented by their starting time (typically 0), step size (i.e. sampling time or TR), and number of elements:
>>> series = cifti2.SeriesAxis(start=0, step=100, size=5000)
So a header for fMRI data with a TR of 100 ms covering the left cortex and thalamus with 5000 timepoints could be created with
>>> type(cifti2.Cifti2Header.from_axes((series, bm_cortex + bm_thal)))
<class 'nibabel.cifti2.cifti2.Cifti2Header'>
Similarly the curvature and cortical thickness on the left cortex could be stored using a header like:
>>> type(cifti2.Cifti2Header.from_axes((cifti2.ScalarAxis(['curvature', 'thickness']),
... bm_cortex)))
<class 'nibabel.cifti2.cifti2.Cifti2Header'>
|
Abstract class for any object describing the rows or columns of a CIFTI-2 vector/matrix |
|
Each row/column in the CIFTI-2 vector/matrix represents a single vertex or voxel |
|
Defines CIFTI-2 axis for label array. |
|
Each row/column in the CIFTI-2 vector/matrix represents a parcel of voxels/vertices |
|
Along this axis of the CIFTI-2 vector/matrix each row/column has been given a unique name and optionally metadata |
|
Along this axis of the CIFTI-2 vector/matrix the rows/columns increase monotonously in time |
|
Parses the MatrixIndicesMap to find the appropriate CIFTI-2 axis describing the rows or columns |
|
Converts the axes describing the rows/columns of a CIFTI-2 vector/matrix to a Cifti2Header |
Module: cifti2.parse_cifti2
¶
|
|
|
Class to parse an XML string into a CIFTI-2 header object |
Cifti2BrainModel
¶
- class nibabel.cifti2.cifti2.Cifti2BrainModel(index_offset=None, index_count=None, model_type=None, brain_structure=None, n_surface_vertices=None, voxel_indices_ijk=None, vertex_indices=None)¶
Bases:
XmlSerializable
Element representing a mapping of the dimension to vertex or voxels.
Mapping to vertices of voxels must be specified.
Description - Maps a range of indices to surface vertices or voxels when IndicesMapToDataType is “CIFTI_INDEX_TYPE_BRAIN_MODELS.”
Attributes
IndexOffset - The matrix index of the first brainordinate of this BrainModel. Note that matrix indices are zero-based.
IndexCount - Number of surface vertices or voxels in this brain model, must be positive.
ModelType - Type of model representing the brain structure (surface or voxels). Valid values are listed in the table below.
BrainStructure - Identifies the brain structure. Valid values for BrainStructure are listed in the table below. However, if the needed structure is not listed in the table, a message should be posted to the CIFTI Forum so that a standardized name can be created for the structure and added to the table.
SurfaceNumberOfVertices - When ModelType is CIFTI_MODEL_TYPE_SURFACE this attribute contains the actual (or true) number of vertices in the surface that is associated with this BrainModel. When this BrainModel represents all vertices in the surface, this value is the same as IndexCount. When this BrainModel represents only a subset of the surface’s vertices, IndexCount will be less than this value.
Child Elements
VertexIndices (0…1)
VoxelIndicesIJK (0…1)
Text Content: [NA]
Parent Element - MatrixIndicesMap
For ModelType values, see CIFTI_MODEL_TYPES module attribute.
For BrainStructure values, see CIFTI_BRAIN_STRUCTURES model attribute.
- Attributes:
- index_offsetint
Start of the mapping
- index_countint
Number of elements in the array to be mapped
- model_typestr
One of CIFTI_MODEL_TYPES
- brain_structurestr
One of CIFTI_BRAIN_STRUCTURES
- surface_number_of_verticesint
Number of vertices in the surface. Use only for surface-type structure
- voxel_indices_ijkCifti2VoxelIndicesIJK, optional
Indices on the image towards where the array indices are mapped
- vertex_indicesCifti2VertexIndices, optional
Indices of the vertices towards where the array indices are mapped
- __init__(index_offset=None, index_count=None, model_type=None, brain_structure=None, n_surface_vertices=None, voxel_indices_ijk=None, vertex_indices=None)¶
- property vertex_indices¶
- property voxel_indices_ijk¶
Cifti2Header
¶
- class nibabel.cifti2.cifti2.Cifti2Header(matrix=None, version='2.0')¶
Bases:
FileBasedHeader
,XmlSerializable
Class for CIFTI-2 header extension
- __init__(matrix=None, version='2.0')¶
- classmethod from_axes(axes)¶
Creates a new Cifti2 header based on the Cifti2 axes
- Parameters:
- axestuple of :class`.cifti2_axes.Axis`
sequence of Cifti2 axes describing each row/column of the matrix to be stored
- Returns:
- headerCifti2Header
new header describing the rows/columns in a format consistent with Cifti2
- get_axis(index)¶
Generates the Cifti2 axis for a given dimension
- Parameters:
- indexint
Dimension for which we want to obtain the mapping.
- Returns:
- axis
cifti2_axes.Axis
- axis
- get_index_map(index)¶
Cifti2 Mapping class for a given index
- Parameters:
- indexint
Index for which we want to obtain the mapping. Must be in the mapped_indices sequence.
- Returns:
- cifti2_mapCifti2MatrixIndicesMap
Returns the Cifti2MatrixIndicesMap corresponding to the given index.
- property mapped_indices¶
List of matrix indices that are mapped
- classmethod may_contain_header(binaryblock)¶
- property number_of_mapped_indices¶
Number of mapped indices
Cifti2HeaderError
¶
Cifti2Image
¶
- class nibabel.cifti2.cifti2.Cifti2Image(dataobj=None, header=None, nifti_header=None, extra=None, file_map=None, dtype=None)¶
Bases:
DataobjImage
,SerializableImage
Class for single file CIFTI-2 format image
Initialize image
The image is a combination of (dataobj, header), with optional metadata in nifti_header (a NIfTI2 header). There may be more metadata in the mapping extra. Filename / file-like objects can also go in the file_map mapping.
- Parameters:
- dataobjobject
Object containing image data. It should be some object that returns an array from
np.asanyarray
. It should have ashape
attribute or property.- headerCifti2Header instance or sequence of
cifti2_axes.Axis
Header with data for / from XML part of CIFTI-2 format. Alternatively a sequence of cifti2_axes.Axis objects can be provided describing each dimension of the array.
- nifti_headerNone or mapping or NIfTI2 header instance, optional
Metadata for NIfTI2 component of this format.
- extraNone or mapping
Extra metadata not captured by header or nifti_header.
- file_mapmapping, optional
Mapping giving file information for this image format.
- __init__(dataobj=None, header=None, nifti_header=None, extra=None, file_map=None, dtype=None)¶
Initialize image
The image is a combination of (dataobj, header), with optional metadata in nifti_header (a NIfTI2 header). There may be more metadata in the mapping extra. Filename / file-like objects can also go in the file_map mapping.
- Parameters:
- dataobjobject
Object containing image data. It should be some object that returns an array from
np.asanyarray
. It should have ashape
attribute or property.- headerCifti2Header instance or sequence of
cifti2_axes.Axis
Header with data for / from XML part of CIFTI-2 format. Alternatively a sequence of cifti2_axes.Axis objects can be provided describing each dimension of the array.
- nifti_headerNone or mapping or NIfTI2 header instance, optional
Metadata for NIfTI2 component of this format.
- extraNone or mapping
Extra metadata not captured by header or nifti_header.
- file_mapmapping, optional
Mapping giving file information for this image format.
- classmethod from_file_map(file_map, *, mmap=True, keep_file_open=None)¶
Load a CIFTI-2 image from a file_map
- Parameters:
- file_mapfile_map
- Returns:
- imgCifti2Image
Returns a Cifti2Image
- classmethod from_image(img)¶
Class method to create new instance of own class from img
- Parameters:
- imginstance
In fact, an object with the API of
DataobjImage
.
- Returns:
- cimginstance
Image, of our own class
- get_data_dtype()¶
- header_class¶
alias of
Cifti2Header
- property nifti_header¶
- set_data_dtype(dtype)¶
- to_file_map(file_map=None, dtype=None)¶
Write image to file_map or contained
self.file_map
- Parameters:
- file_mapNone or mapping, optional
files mapping. If None (default) use object’s
file_map
attribute instead.
- Returns:
- None
- update_headers()¶
Harmonize NIfTI headers with image data
Ensures that the NIfTI-2 header records the data shape in the last three
dim
fields. Per the spec:Because the first four dimensions in NIfTI are reserved for space and time, the CIFTI dimensions are stored in the NIfTI header in dim[5] and up, where dim[5] is the length of the first CIFTI dimension (number of values in a row), dim[6] is the length of the second CIFTI dimension, and dim[7] is the length of the third CIFTI dimension, if applicable. The fields dim[1] through dim[4] will be 1; dim[0] will be 6 or 7, depending on whether a third matrix dimension exists.
>>> import numpy as np >>> data = np.zeros((2,3,4)) >>> img = Cifti2Image(data) >>> img.shape == (2, 3, 4) True >>> img.update_headers() >>> img.nifti_header.get_data_shape() == (1, 1, 1, 1, 2, 3, 4) True >>> img.shape == (2, 3, 4) True
Cifti2Label
¶
- class nibabel.cifti2.cifti2.Cifti2Label(key=0, label='', red=0.0, green=0.0, blue=0.0, alpha=0.0)¶
Bases:
XmlSerializable
CIFTI-2 label: association of integer key with a name and RGBA values
For all color components, value is floating point with range 0.0 to 1.0.
Description - Associates a label key value with a name and a display color.
Attributes
Key - Integer, data value which is assigned this name and color.
Red - Red color component for label. Value is floating point with range 0.0 to 1.0.
Green - Green color component for label. Value is floating point with range 0.0 to 1.0.
Blue - Blue color component for label. Value is floating point with range 0.0 to 1.0.
Alpha - Alpha color component for label. Value is floating point with range 0.0 to 1.0.
Child Elements: [NA]
Text Content - Name of the label.
Parent Element - LabelTable
- Attributes:
- keyint, optional
Integer, data value which is assigned this name and color.
- labelstr, optional
Name of the label.
- redfloat, optional
Red color component for label (between 0 and 1).
- greenfloat, optional
Green color component for label (between 0 and 1).
- bluefloat, optional
Blue color component for label (between 0 and 1).
- alphafloat, optional
Alpha color component for label (between 0 and 1).
- __init__(key=0, label='', red=0.0, green=0.0, blue=0.0, alpha=0.0)¶
- property rgba¶
Returns RGBA as tuple
Cifti2LabelTable
¶
- class nibabel.cifti2.cifti2.Cifti2LabelTable¶
Bases:
XmlSerializable
,MutableMapping
CIFTI-2 label table: a sequence of
Cifti2Label
sDescription - Used by NamedMap when IndicesMapToDataType is “CIFTI_INDEX_TYPE_LABELS” in order to associate names and display colors with label keys. Note that LABELS is the only mapping type that uses a LabelTable. Display coloring of continuous-valued data is not specified by CIFTI-2.
Attributes: [NA]
Child Elements
Label (0…N)
Text Content: [NA]
Parent Element - NamedMap
- __init__()¶
- append(label)¶
Cifti2Matrix
¶
- class nibabel.cifti2.cifti2.Cifti2Matrix¶
Bases:
XmlSerializable
,MutableSequence
CIFTI-2 Matrix object
This is a list-like container where the elements are instances of
Cifti2MatrixIndicesMap
.Description: contains child elements that describe the meaning of the values in the matrix.
Attributes: [NA]
Child Elements
MetaData (0 .. 1)
MatrixIndicesMap (1 .. N)
Text Content: [NA]
Parent Element: CIFTI
For each matrix (data) dimension, exactly one MatrixIndicesMap element must list it in the AppliesToMatrixDimension attribute.
- __init__()¶
- get_axis(index)¶
Generates the Cifti2 axis for a given dimension
- Parameters:
- indexint
Dimension for which we want to obtain the mapping.
- Returns:
- axis
cifti2_axes.Axis
- axis
- get_data_shape()¶
Returns data shape expected based on the CIFTI-2 header
Any dimensions omitted in the CIFTI-2 header will be given a default size of None.
- get_index_map(index)¶
Cifti2 Mapping class for a given index
- Parameters:
- indexint
Index for which we want to obtain the mapping. Must be in the mapped_indices sequence.
- Returns:
- cifti2_mapCifti2MatrixIndicesMap
Returns the Cifti2MatrixIndicesMap corresponding to the given index.
- insert(index, value)¶
S.insert(index, value) – insert value before index
- property mapped_indices¶
List of matrix indices that are mapped
- property metadata¶
Cifti2MatrixIndicesMap
¶
- class nibabel.cifti2.cifti2.Cifti2MatrixIndicesMap(applies_to_matrix_dimension, indices_map_to_data_type, number_of_series_points=None, series_exponent=None, series_start=None, series_step=None, series_unit=None, maps=[])¶
Bases:
XmlSerializable
,MutableSequence
Class for Matrix Indices Map
Description - Provides a mapping between matrix indices and their interpretation.
Attributes
AppliesToMatrixDimension - Lists the dimension(s) of the matrix to which this MatrixIndicesMap applies. The dimensions of the matrix start at zero (dimension 0 describes the indices along the first dimension, dimension 1 describes the indices along the second dimension, etc.). If this MatrixIndicesMap applies to more than one matrix dimension, the values are separated by a comma.
IndicesMapToDataType - Type of data to which the MatrixIndicesMap applies.
NumberOfSeriesPoints - Indicates how many samples there are in a series mapping type. For example, this could be the number of timepoints in a timeseries.
SeriesExponent - Integer, SeriesStart and SeriesStep must be multiplied by 10 raised to the power of the value of this attribute to give the actual values assigned to indices (e.g., if SeriesStart is “5” and SeriesExponent is “-3”, the value of the first series point is 0.005).
SeriesStart - Indicates what quantity should be assigned to the first series point.
SeriesStep - Indicates amount of change between each series point.
SeriesUnit - Indicates the unit of the result of multiplying SeriesStart and SeriesStep by 10 to the power of SeriesExponent.
Child Elements
BrainModel (0…N)
NamedMap (0…N)
Parcel (0…N)
Surface (0…N)
Volume (0…1)
Text Content: [NA]
Parent Element - Matrix
- Attributes:
- applies_to_matrix_dimensionlist of ints
Dimensions of this matrix that follow this mapping
- indices_map_to_data_typestr one of CIFTI_MAP_TYPES
Type of mapping to the matrix indices
- number_of_series_pointsint, optional
If it is a series, number of points in the series
- series_exponentint, optional
If it is a series the exponent of the increment
- series_startfloat, optional
If it is a series, starting time
- series_stepfloat, optional
If it is a series, step per element
- series_unitstr, optional
If it is a series, units
- __init__(applies_to_matrix_dimension, indices_map_to_data_type, number_of_series_points=None, series_exponent=None, series_start=None, series_step=None, series_unit=None, maps=[])¶
- property brain_models¶
- insert(index, value)¶
S.insert(index, value) – insert value before index
- property named_maps¶
- property parcels¶
- property surfaces¶
- property volume¶
Cifti2MetaData
¶
- class nibabel.cifti2.cifti2.Cifti2MetaData(*args, **kwargs)¶
Bases:
CaretMetaData
A list of name-value pairs
Description - Provides a simple method for user-supplied metadata that associates names with values.
Attributes: [NA]
Child Elements
MD (0…N)
Text Content: [NA]
Parent Elements - Matrix, NamedMap
MD elements are a single metadata entry consisting of a name and a value.
- Attributes:
- datalist of (name, value) tuples
- __init__(*args, **kwargs)¶
- property data¶
- difference_update(metadata)¶
Remove metadata key-value pairs
- Parameters:
- metadatadict-like datatype
- Returns:
- None
Cifti2NamedMap
¶
- class nibabel.cifti2.cifti2.Cifti2NamedMap(map_name=None, metadata=None, label_table=None)¶
Bases:
XmlSerializable
CIFTI-2 named map: association of name and optional data with a map index
Associates a name, optional metadata, and possibly a LabelTable with an index in a map.
Description - Associates a name, optional metadata, and possibly a LabelTable with an index in a map.
Attributes: [NA]
Child Elements
MapName (1)
LabelTable (0…1)
MetaData (0…1)
Text Content: [NA]
Parent Element - MatrixIndicesMap
- Attributes:
- map_namestr
Name of map
- metadataNone or Cifti2MetaData
Metadata associated with named map
- label_tableNone or Cifti2LabelTable
Label table associated with named map
- __init__(map_name=None, metadata=None, label_table=None)¶
- property label_table¶
- property metadata¶
Cifti2Parcel
¶
- class nibabel.cifti2.cifti2.Cifti2Parcel(name=None, voxel_indices_ijk=None, vertices=None)¶
Bases:
XmlSerializable
CIFTI-2 parcel: association of a name with vertices and/or voxels
Description - Associates a name, plus vertices and/or voxels, with an index.
Attributes
Name - The name of the parcel
Child Elements
Vertices (0…N)
VoxelIndicesIJK (0…1)
Text Content: [NA]
Parent Element - MatrixIndicesMap
- Attributes:
- namestr
Name of parcel
- voxel_indices_ijkNone or Cifti2VoxelIndicesIJK
Voxel indices associated with parcel
- verticeslist of Cifti2Vertices
Vertices associated with parcel
- __init__(name=None, voxel_indices_ijk=None, vertices=None)¶
- append_cifti_vertices(vertices)¶
Appends a Cifti2Vertices element to the Cifti2Parcel
- Parameters:
- verticesCifti2Vertices
- pop_cifti2_vertices(ith)¶
Pops the ith vertices element from the Cifti2Parcel
- property voxel_indices_ijk¶
Cifti2Surface
¶
- class nibabel.cifti2.cifti2.Cifti2Surface(brain_structure=None, surface_number_of_vertices=None)¶
Bases:
XmlSerializable
Cifti surface: association of brain structure and number of vertices
Description - Specifies the number of vertices for a surface, when IndicesMapToDataType is “CIFTI_INDEX_TYPE_PARCELS.” This is separate from the Parcel element because there can be multiple parcels on one surface, and one parcel may involve multiple surfaces.
Attributes
BrainStructure - A string from the BrainStructure list to identify what surface structure this element refers to (usually left cortex, right cortex, or cerebellum).
SurfaceNumberOfVertices - The number of vertices that this structure’s surface contains.
Child Elements: [NA]
Text Content: [NA]
Parent Element - MatrixIndicesMap
- Attributes:
- brain_structurestr
Name of brain structure
- surface_number_of_verticesint
Number of vertices on surface
- __init__(brain_structure=None, surface_number_of_vertices=None)¶
Cifti2TransformationMatrixVoxelIndicesIJKtoXYZ
¶
- class nibabel.cifti2.cifti2.Cifti2TransformationMatrixVoxelIndicesIJKtoXYZ(meter_exponent=None, matrix=None)¶
Bases:
XmlSerializable
Matrix that translates voxel indices to spatial coordinates
Description - Contains a matrix that translates Voxel IJK Indices to spatial XYZ coordinates (+X=>right, +Y=>anterior, +Z=> superior). The resulting coordinate is the center of the voxel.
Attributes
MeterExponent - Integer, specifies that the coordinate result from the transformation matrix should be multiplied by 10 to this power to get the spatial coordinates in meters (e.g., if this is “-3”, then the transformation matrix is in millimeters).
Child Elements: [NA]
Text Content - Sixteen floating-point values, in row-major order, that form a 4x4 homogeneous transformation matrix.
Parent Element - Volume
- Attributes:
- meter_exponentint
See attribute description above.
- matrixarray-like shape (4, 4)
Affine transformation matrix from voxel indices to RAS space.
- __init__(meter_exponent=None, matrix=None)¶
Cifti2VertexIndices
¶
- class nibabel.cifti2.cifti2.Cifti2VertexIndices(indices=None)¶
Bases:
XmlSerializable
,MutableSequence
CIFTI-2 vertex indices: vertex indices for an associated brain model
The vertex indices (which are independent for each surface, and zero-based) that are used in this brain model[.] The parent BrainModel’s
index_count
indicates the number of indices.Description - Contains a list of vertex indices for a BrainModel with ModelType equal to CIFTI_MODEL_TYPE_SURFACE.
Attributes: [NA]
Child Elements: [NA]
Text Content - The vertex indices (which are independent for each surface, and zero-based) that are used in this brain model, with each index separated by a whitespace character. The parent BrainModel’s IndexCount attribute indicates the number of indices in this element’s content.
Parent Element - BrainModel
- __init__(indices=None)¶
- insert(index, value)¶
S.insert(index, value) – insert value before index
Cifti2Vertices
¶
- class nibabel.cifti2.cifti2.Cifti2Vertices(brain_structure=None, vertices=None)¶
Bases:
XmlSerializable
,MutableSequence
CIFTI-2 vertices - association of brain structure and a list of vertices
Description - Contains a BrainStructure type and a list of vertex indices within a Parcel.
Attributes
BrainStructure - A string from the BrainStructure list to identify what surface this vertex list is from (usually left cortex, right cortex, or cerebellum).
Child Elements: [NA]
Text Content - Vertex indices (which are independent for each surface, and zero-based) separated by whitespace characters.
Parent Element - Parcel
The class behaves like a list of Vertex indices (which are independent for each surface, and zero-based)
- Attributes:
- brain_structurestr
A string from the BrainStructure list to identify what surface this vertex list is from (usually left cortex, right cortex, or cerebellum).
- __init__(brain_structure=None, vertices=None)¶
- insert(index, value)¶
S.insert(index, value) – insert value before index
Cifti2Volume
¶
- class nibabel.cifti2.cifti2.Cifti2Volume(volume_dimensions=None, transform_matrix=None)¶
Bases:
XmlSerializable
CIFTI-2 volume: information about a volume for mappings that use voxels
Description - Provides information about the volume for any mappings that use voxels.
Attributes
VolumeDimensions - Three integer values separated by commas, the lengths of the three volume file dimensions that are related to spatial coordinates, in number of voxels. Voxel indices (which are zero-based) that are used in the mapping that this element applies to must be within these dimensions.
Child Elements
TransformationMatrixVoxelIndicesIJKtoXYZ (1)
Text Content: [NA]
Parent Element - MatrixIndicesMap
- Attributes:
- volume_dimensionsarray-like shape (3,)
See attribute description above.
- transformation_matrix_voxel_indices_ijk_to_xyzCifti2TransformationMatrixVoxelIndicesIJKtoXYZ
Matrix that translates voxel indices to spatial coordinates
- __init__(volume_dimensions=None, transform_matrix=None)¶
Cifti2VoxelIndicesIJK
¶
- class nibabel.cifti2.cifti2.Cifti2VoxelIndicesIJK(indices=None)¶
Bases:
XmlSerializable
,MutableSequence
CIFTI-2 VoxelIndicesIJK: Set of voxel indices contained in a structure
Description - Identifies the voxels that model a brain structure, or participate in a parcel. Note that when this is a child of BrainModel, the IndexCount attribute of the BrainModel indicates the number of voxels contained in this element.
Attributes: [NA]
Child Elements: [NA]
Text Content - IJK indices (which are zero-based) of each voxel in this brain model or parcel, with each index separated by a whitespace character. There are three indices per voxel. If the parent element is BrainModel, then the BrainModel element’s IndexCount attribute indicates the number of triplets (IJK indices) in this element’s content.
Parent Elements - BrainModel, Parcel
Each element of this sequence is a triple of integers.
- __init__(indices=None)¶
- insert(index, value)¶
S.insert(index, value) – insert value before index
LimitedNifti2Header
¶
- class nibabel.cifti2.cifti2.LimitedNifti2Header(binaryblock=None, endianness=None, check=True, extensions=())¶
Bases:
Nifti2Header
Initialize header from binary data block and extensions
- __init__(binaryblock=None, endianness=None, check=True, extensions=())¶
Initialize header from binary data block and extensions
Axis
¶
- class nibabel.cifti2.cifti2_axes.Axis¶
Bases:
ABC
Abstract class for any object describing the rows or columns of a CIFTI-2 vector/matrix
Mainly used for type checking.
Base class for the following concrete CIFTI-2 axes:
BrainModelAxis
: each row/column is a voxel or vertexParcelsAxis
: each row/column is a group of voxels and/or verticesScalarAxis
: each row/column has a unique name with optional meta-dataLabelAxis
: each row/column has a unique name and label table with optional meta-dataSeriesAxis
: each row/column is a timepoint, which increases monotonically
- __init__(*args, **kwargs)¶
- property size¶
BrainModelAxis
¶
- class nibabel.cifti2.cifti2_axes.BrainModelAxis(name, voxel=None, vertex=None, affine=None, volume_shape=None, nvertices=None)¶
Bases:
Axis
Each row/column in the CIFTI-2 vector/matrix represents a single vertex or voxel
This Axis describes which vertex/voxel is represented by each row/column.
New BrainModelAxis axes can be constructed by passing on the greyordinate brain-structure names and voxel/vertex indices to the constructor or by one of the factory methods:
from_mask()
: creates surface or volumetric BrainModelAxis axis from respectively 1D or 3D masksfrom_surface()
: creates a surface BrainModelAxis axis
The resulting BrainModelAxis axes can be concatenated by adding them together.
- Parameters:
- namearray_like
brain structure name or (N, ) string array with the brain structure names
- voxelarray_like, optional
(N, 3) array with the voxel indices (can be omitted for CIFTI-2 files only covering the surface)
- vertexarray_like, optional
(N, ) array with the vertex indices (can be omitted for volumetric CIFTI-2 files)
- affinearray_like, optional
(4, 4) array mapping voxel indices to mm space (not needed for CIFTI-2 files only covering the surface)
- volume_shapetuple of three integers, optional
shape of the volume in which the voxels were defined (not needed for CIFTI-2 files only covering the surface)
- nverticesdict from string to integer, optional
maps names of surface elements to integers (not needed for volumetric CIFTI-2 files)
- __init__(name, voxel=None, vertex=None, affine=None, volume_shape=None, nvertices=None)¶
New BrainModelAxis axes can be constructed by passing on the greyordinate brain-structure names and voxel/vertex indices to the constructor or by one of the factory methods:
from_mask()
: creates surface or volumetric BrainModelAxis axis from respectively 1D or 3D masksfrom_surface()
: creates a surface BrainModelAxis axis
The resulting BrainModelAxis axes can be concatenated by adding them together.
- Parameters:
- namearray_like
brain structure name or (N, ) string array with the brain structure names
- voxelarray_like, optional
(N, 3) array with the voxel indices (can be omitted for CIFTI-2 files only covering the surface)
- vertexarray_like, optional
(N, ) array with the vertex indices (can be omitted for volumetric CIFTI-2 files)
- affinearray_like, optional
(4, 4) array mapping voxel indices to mm space (not needed for CIFTI-2 files only covering the surface)
- volume_shapetuple of three integers, optional
shape of the volume in which the voxels were defined (not needed for CIFTI-2 files only covering the surface)
- nverticesdict from string to integer, optional
maps names of surface elements to integers (not needed for volumetric CIFTI-2 files)
- property affine¶
Affine of the volumetric image in which the greyordinate voxels were defined
- classmethod from_index_mapping(mim)¶
Creates a new BrainModel axis based on a CIFTI-2 dataset
- Parameters:
- Returns:
- BrainModelAxis
- classmethod from_mask(mask, name='other', affine=None)¶
Creates a new BrainModelAxis axis describing the provided mask
- Parameters:
- maskarray_like
all non-zero voxels will be included in the BrainModelAxis axis should be (Nx, Ny, Nz) array for volume mask or (Nvertex, ) array for surface mask
- namestr, optional
Name of the brain structure (e.g. ‘CortexRight’, ‘thalamus_left’ or ‘brain_stem’)
- affinearray_like, optional
(4, 4) array with the voxel to mm transformation (defaults to identity matrix) Argument will be ignored for surface masks
- Returns:
- BrainModelAxis which covers the provided mask
- classmethod from_surface(vertices, nvertex, name='Other')¶
Creates a new BrainModelAxis axis describing the vertices on a surface
- Parameters:
- verticesarray_like
indices of the vertices on the surface
- nvertexint
total number of vertices on the surface
- namestr
Name of the brain structure (e.g. ‘CortexLeft’ or ‘CortexRight’)
- Returns:
- BrainModelAxis which covers (part of) the surface
- get_element(index)¶
Describes a single element from the axis
- Parameters:
- indexint
Indexes the row/column of interest
- Returns:
- tuple with 3 elements
- str, ‘CIFTI_MODEL_TYPE_SURFACE’ for vertex or ‘CIFTI_MODEL_TYPE_VOXELS’ for voxel
- vertex index if it is a surface element, otherwise array with 3 voxel indices
- structure.BrainStructure object describing the brain structure the element was taken from
- iter_structures()¶
Iterates over all brain structures in the order that they appear along the axis
- Yields:
- tuple with 3 elements:
- CIFTI-2 brain structure name
- slice to select the data associated with the brain structure from the tensor
- brain model covering that specific brain structure
- property name¶
The brain structure to which the voxel/vertices of belong
- property surface_mask¶
(N, ) boolean array which is true for any element on the surface
- static to_cifti_brain_structure_name(name)¶
Attempts to convert the name of an anatomical region in a format recognized by CIFTI-2
This function returns:
the name if it is in the CIFTI-2 format already
if the name is a tuple the first element is assumed to be the structure name while the second is assumed to be the hemisphere (left, right or both). The latter will default to both.
names like left_cortex, cortex_left, LeftCortex, or CortexLeft will be converted to CIFTI_STRUCTURE_CORTEX_LEFT
see
nibabel.cifti2.tests.test_name()
for examples of which conversions are possible- Parameters:
- name: iterable of 2-element tuples of integer and string
input name of an anatomical region
- Returns:
- CIFTI-2 compatible name
- Raises:
- ValueError: raised if the input name does not match a known anatomical structure in CIFTI-2
- to_mapping(dim)¶
Converts the brain model axis to a MatrixIndicesMap for storage in CIFTI-2 format
- Parameters:
- dimint
which dimension of the CIFTI-2 vector/matrix is described by this dataset (zero-based)
- Returns:
- property volume_mask¶
(N, ) boolean array which is true for any element on the surface
- property volume_shape¶
Shape of the volumetric image in which the greyordinate voxels were defined
LabelAxis
¶
- class nibabel.cifti2.cifti2_axes.LabelAxis(name, label, meta=None)¶
Bases:
Axis
Defines CIFTI-2 axis for label array.
Along this axis of the CIFTI-2 vector/matrix each row/column has been given a unique name, label table, and optionally metadata
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- labelarray_like
single dictionary or (N, ) object array with dictionaries mapping from integers to (name, (R, G, B, A)), where name is a string and R, G, B, and A are floats between 0 and 1 giving the colour and alpha (i.e., transparency)
- metaarray_like, optional
(N, ) object array with a dictionary of metadata for each row/column
- __init__(name, label, meta=None)¶
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- labelarray_like
single dictionary or (N, ) object array with dictionaries mapping from integers to (name, (R, G, B, A)), where name is a string and R, G, B, and A are floats between 0 and 1 giving the colour and alpha (i.e., transparency)
- metaarray_like, optional
(N, ) object array with a dictionary of metadata for each row/column
- classmethod from_index_mapping(mim)¶
Creates a new Label axis based on a CIFTI-2 dataset
- Parameters:
- Returns:
- LabelAxis
- get_element(index)¶
Describes a single element from the axis
- Parameters:
- indexint
Indexes the row/column of interest
- Returns:
- tuple with 2 elements
- unicode name of the row/column
- dictionary with the label table
- dictionary with the element metadata
- to_mapping(dim)¶
Converts the hcp_labels to a MatrixIndicesMap for storage in CIFTI-2 format
- Parameters:
- dimint
which dimension of the CIFTI-2 vector/matrix is described by this dataset (zero-based)
- Returns:
ParcelsAxis
¶
- class nibabel.cifti2.cifti2_axes.ParcelsAxis(name, voxels, vertices, affine=None, volume_shape=None, nvertices=None)¶
Bases:
Axis
Each row/column in the CIFTI-2 vector/matrix represents a parcel of voxels/vertices
This Axis describes which parcel is represented by each row/column.
Individual parcels can be accessed based on their name, using
parcel = parcel_axis[name]
Use of this constructor is not recommended. New ParcelsAxis axes can be constructed more easily from a sequence of BrainModelAxis axes using
from_brain_models()
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- voxelsarray_like
(N, ) object array each containing a sequence of voxels. For each parcel the voxels are represented by a (M, 3) index array
- verticesarray_like
(N, ) object array each containing a sequence of vertices. For each parcel the vertices are represented by a mapping from brain structure name to (M, ) index array
- affinearray_like, optional
(4, 4) array mapping voxel indices to mm space (not needed for CIFTI-2 files only covering the surface)
- volume_shapetuple of three integers, optional
shape of the volume in which the voxels were defined (not needed for CIFTI-2 files only covering the surface)
- nverticesdict from string to integer, optional
maps names of surface elements to integers (not needed for volumetric CIFTI-2 files)
- __init__(name, voxels, vertices, affine=None, volume_shape=None, nvertices=None)¶
Use of this constructor is not recommended. New ParcelsAxis axes can be constructed more easily from a sequence of BrainModelAxis axes using
from_brain_models()
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- voxelsarray_like
(N, ) object array each containing a sequence of voxels. For each parcel the voxels are represented by a (M, 3) index array
- verticesarray_like
(N, ) object array each containing a sequence of vertices. For each parcel the vertices are represented by a mapping from brain structure name to (M, ) index array
- affinearray_like, optional
(4, 4) array mapping voxel indices to mm space (not needed for CIFTI-2 files only covering the surface)
- volume_shapetuple of three integers, optional
shape of the volume in which the voxels were defined (not needed for CIFTI-2 files only covering the surface)
- nverticesdict from string to integer, optional
maps names of surface elements to integers (not needed for volumetric CIFTI-2 files)
- property affine¶
Affine of the volumetric image in which the greyordinate voxels were defined
- classmethod from_brain_models(named_brain_models)¶
Creates a Parcel axis from a list of BrainModelAxis axes with names
- Parameters:
- named_brain_modelsiterable of 2-element tuples of string and BrainModelAxis
list of (parcel name, brain model representation) pairs defining each parcel
- Returns:
- ParcelsAxis
- classmethod from_index_mapping(mim)¶
Creates a new Parcels axis based on a CIFTI-2 dataset
- Parameters:
- mim
cifti2.Cifti2MatrixIndicesMap
- mim
- Returns:
- ParcelsAxis
- get_element(index)¶
Describes a single element from the axis
- Parameters:
- indexint
Indexes the row/column of interest
- Returns:
- tuple with 3 elements
- unicode name of the parcel
- (M, 3) int array with voxel indices
- dict from string to (K, ) int array with vertex indices
for a specific surface brain structure
- to_mapping(dim)¶
Converts the Parcel to a MatrixIndicesMap for storage in CIFTI-2 format
- Parameters:
- dimint
which dimension of the CIFTI-2 vector/matrix is described by this dataset (zero-based)
- Returns:
cifti2.Cifti2MatrixIndicesMap
- property volume_shape¶
Shape of the volumetric image in which the greyordinate voxels were defined
ScalarAxis
¶
- class nibabel.cifti2.cifti2_axes.ScalarAxis(name, meta=None)¶
Bases:
Axis
Along this axis of the CIFTI-2 vector/matrix each row/column has been given a unique name and optionally metadata
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- metaarray_like
(N, ) object array with a dictionary of metadata for each row/column. Defaults to empty dictionary
- __init__(name, meta=None)¶
- Parameters:
- namearray_like
(N, ) string array with the parcel names
- metaarray_like
(N, ) object array with a dictionary of metadata for each row/column. Defaults to empty dictionary
- classmethod from_index_mapping(mim)¶
Creates a new Scalar axis based on a CIFTI-2 dataset
- Parameters:
- Returns:
- ScalarAxis
- get_element(index)¶
Describes a single element from the axis
- Parameters:
- indexint
Indexes the row/column of interest
- Returns:
- tuple with 2 elements
- unicode name of the row/column
- dictionary with the element metadata
- to_mapping(dim)¶
Converts the hcp_labels to a MatrixIndicesMap for storage in CIFTI-2 format
- Parameters:
- dimint
which dimension of the CIFTI-2 vector/matrix is described by this dataset (zero-based)
- Returns:
SeriesAxis
¶
- class nibabel.cifti2.cifti2_axes.SeriesAxis(start, step, size, unit='SECOND')¶
Bases:
Axis
Along this axis of the CIFTI-2 vector/matrix the rows/columns increase monotonously in time
This Axis describes the time point of each row/column.
Creates a new SeriesAxis axis
- Parameters:
- startfloat
starting time point
- stepfloat
sampling time (TR)
- sizeint
number of time points
- unitstr
Unit of the step size (one of ‘second’, ‘hertz’, ‘meter’, or ‘radian’)
- __init__(start, step, size, unit='SECOND')¶
Creates a new SeriesAxis axis
- Parameters:
- startfloat
starting time point
- stepfloat
sampling time (TR)
- sizeint
number of time points
- unitstr
Unit of the step size (one of ‘second’, ‘hertz’, ‘meter’, or ‘radian’)
- classmethod from_index_mapping(mim)¶
Creates a new SeriesAxis axis based on a CIFTI-2 dataset
- Parameters:
- Returns:
- SeriesAxis
- get_element(index)¶
Gives the time point of a specific row/column
- Parameters:
- indexint
Indexes the row/column of interest
- Returns:
- float
- size = None¶
- property time¶
- to_mapping(dim)¶
Converts the SeriesAxis to a MatrixIndicesMap for storage in CIFTI-2 format
- Parameters:
- dimint
which dimension of the CIFTI-2 vector/matrix is described by this dataset (zero-based)
- Returns:
cifti2.Cifti2MatrixIndicesMap
- property unit¶
from_index_mapping¶
to_header¶
- nibabel.cifti2.cifti2_axes.to_header(axes)¶
Converts the axes describing the rows/columns of a CIFTI-2 vector/matrix to a Cifti2Header
- Parameters:
- axesiterable of
Axis
objects one or more axes describing each dimension in turn
- axesiterable of
- Returns:
- header
cifti2.Cifti2Header
- header
Cifti2Extension
¶
- class nibabel.cifti2.parse_cifti2.Cifti2Extension(code: int | str, content: bytes = b'', object: T | None = None)¶
Bases:
Nifti1Extension
[Cifti2Header
]- Parameters:
- codeint or str
Canonical extension code as defined in the NIfTI standard, given either as integer or corresponding label (see
extension_codes
)- contentbytes, optional
Extension content as read from the NIfTI file header.
- objectoptional
Extension content in runtime form.
- __init__(code: int | str, content: bytes = b'', object: T | None = None) None ¶
- Parameters:
- codeint or str
Canonical extension code as defined in the NIfTI standard, given either as integer or corresponding label (see
extension_codes
)- contentbytes, optional
Extension content as read from the NIfTI file header.
- objectoptional
Extension content in runtime form.
Cifti2Parser
¶
- class nibabel.cifti2.parse_cifti2.Cifti2Parser(encoding=None, buffer_size=3500000, verbose=0)¶
Bases:
XmlParser
Class to parse an XML string into a CIFTI-2 header object
- Parameters:
- encodingstr
string containing xml document
- buffer_size: None or int, optional
size of read buffer. None uses default buffer_size from xml.parsers.expat.
- verboseint, optional
amount of output during parsing (0=silent, by default).
- __init__(encoding=None, buffer_size=3500000, verbose=0)¶
- Parameters:
- encodingstr
string containing xml document
- buffer_size: None or int, optional
size of read buffer. None uses default buffer_size from xml.parsers.expat.
- verboseint, optional
amount of output during parsing (0=silent, by default).
- CharacterDataHandler(data)¶
Collect character data chunks pending collation
The parser breaks the data up into chunks of size depending on the buffer_size of the parser. A large bit of character data, with standard parser buffer_size (such as 8K) can easily span many calls to this function. We thus collect the chunks and process them when we hit start or end tags.
- EndElementHandler(name)¶
- StartElementHandler(name, attrs)¶
- flush_chardata()¶
Collate and process collected character data
- property pending_data¶
True if there is character data pending for processing