NiBabel

Access a cacophony of neuro-imaging file formats

Table Of Contents

Next topic

nifti2

Reggie -- the one

nifti1

Read / write access to NIfTI1 image format

NIfTI1 format defined at http://nifti.nimh.nih.gov/nifti-1/

Nifti1Extension(code, content) Baseclass for NIfTI1 header extensions.
Nifti1Extensions Simple extension collection, implemented as a list-subclass.
Nifti1Header([binaryblock, endianness, ...]) Class for NIfTI1 header
Nifti1Image(dataobj, affine[, header, ...]) Class for single file NIfTI1 format image
Nifti1Pair(dataobj, affine[, header, extra, ...]) Class for NIfTI1 format image, header pair
Nifti1PairHeader([binaryblock, endianness, ...]) Class for NIfTI1 pair header
load(filename) Load NIfTI1 single or pair from filename
save(img, filename) Save NIfTI1 single or pair to filename

Nifti1Extension

class nibabel.nifti1.Nifti1Extension(code, content)

Bases: object

Baseclass for NIfTI1 header extensions.

This class is sufficient to handle very simple text-based extensions, such as comment. More sophisticated extensions should/will be supported by dedicated subclasses.

Parameters:

code : int|str

Canonical extension code as defined in the NIfTI standard, given either as integer or corresponding label (see extension_codes)

content : str

Extension content as read from the NIfTI file header. This content is converted into a runtime representation.

__init__(code, content)
Parameters:

code : int|str

Canonical extension code as defined in the NIfTI standard, given either as integer or corresponding label (see extension_codes)

content : str

Extension content as read from the NIfTI file header. This content is converted into a runtime representation.

get_code()

Return the canonical extension type code.

get_content()

Return the extension content in its runtime representation.

get_sizeondisk()

Return the size of the extension in the NIfTI file.

write_to(fileobj, byteswap)

Write header extensions to fileobj

Write starts at fileobj current file position.

Parameters:

fileobj : file-like object

Should implement write method

byteswap : boolean

Flag if byteswapping the data is required.

Returns:

None :

Nifti1Extensions

class nibabel.nifti1.Nifti1Extensions

Bases: list

Simple extension collection, implemented as a list-subclass.

__init__()

x.__init__(...) initializes x; see help(type(x)) for signature

count(ecode)

Returns the number of extensions matching a given ecode.

Parameters:

code : int | str

The ecode can be specified either literal or as numerical value.

classmethod from_fileobj(klass, fileobj, size, byteswap)

Read header extensions from a fileobj

Parameters:

fileobj : file-like object

We begin reading the extensions at the current file position

size : int

Number of bytes to read. If negative, fileobj will be read till its end.

byteswap : boolean

Flag if byteswapping the read data is required.

Returns:

An extension list. This list might be empty in case not extensions :

were present in fileobj. :

get_codes()

Return a list of the extension code of all available extensions

get_sizeondisk()

Return the size of the complete header extensions in the NIfTI file.

write_to(fileobj, byteswap)

Write header extensions to fileobj

Write starts at fileobj current file position.

Parameters:

fileobj : file-like object

Should implement write method

byteswap : boolean

Flag if byteswapping the data is required.

Returns:

None :

Nifti1Header

class nibabel.nifti1.Nifti1Header(binaryblock=None, endianness=None, check=True, extensions=())

Bases: nibabel.spm99analyze.SpmAnalyzeHeader

Class for NIfTI1 header

The NIfTI1 header has many more coded fields than the simpler Analyze variants. NIfTI1 headers also have extensions.

Nifti allows the header to be a separate file, as part of a nifti image / header pair, or to precede the data in a single file. The object needs to know which type it is, in order to manage the voxel offset pointing to the data, extension reading, and writing the correct magic string.

This class handles the header-preceding-data case.

Initialize header from binary data block and extensions

__init__(binaryblock=None, endianness=None, check=True, extensions=())

Initialize header from binary data block and extensions

copy()

Return copy of header

Take reference to extensions as well as copy of header contents

classmethod default_structarr(klass, endianness=None)

Create empty header binary block with given endianness

exts_klass

alias of Nifti1Extensions

classmethod from_fileobj(klass, fileobj, endianness=None, check=True)
classmethod from_header(klass, header=None, check=True)

Class method to create header from another header

Extend Analyze header copy by copying extensions from other Nifti types.

Parameters:

header : Header instance or mapping

a header of this class, or another class of header for conversion to this type

check : {True, False}

whether to check header for integrity

Returns:

hdr : header instance

fresh header instance of our own class

get_best_affine()

Select best of available transforms

get_data_shape()

Get shape of data

Notes

Allows for freesurfer hack for large vectors described in https://github.com/nipy/nibabel/issues/100 and https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/save_nifti.m?spec=svn5022&r=5022#77

Allows for freesurfer hack for 7th order icosahedron surface described in https://github.com/nipy/nibabel/issues/309 https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/load_nifti.m?r=8776#86 https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/save_nifti.m?r=8776#50

Examples

>>> hdr = Nifti1Header()
>>> hdr.get_data_shape()
(0,)
>>> hdr.set_data_shape((1,2,3))
>>> hdr.get_data_shape()
(1, 2, 3)

Expanding number of dimensions gets default zooms

>>> hdr.get_zooms()
(1.0, 1.0, 1.0)
get_dim_info()

Gets NIfTI MRI slice etc dimension information

Returns:

freq : {None,0,1,2}

Which data array axis is frequency encode direction

phase : {None,0,1,2}

Which data array axis is phase encode direction

slice : {None,0,1,2}

Which data array axis is slice encode direction

where ``data array`` is the array returned by ``get_data`` :

Because NIfTI1 files are natively Fortran indexed: :

0 is fastest changing in file 1 is medium changing in file 2 is slowest changing in file

``None`` means the axis appears not to be specified. :

Examples

See set_dim_info function

get_intent(code_repr='label')

Get intent code, parameters and name

Parameters:

code_repr : string

string giving output form of intent code representation. Default is ‘label’; use ‘code’ for integer representation.

Returns:

code : string or integer

intent code, or string describing code

parameters : tuple

parameters for the intent

name : string

intent name

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_intent('t test', (10,), name='some score')
>>> hdr.get_intent()
('t test', (10.0,), 'some score')
>>> hdr.get_intent('code')
(3, (10.0,), 'some score')
get_n_slices()

Return the number of slices

get_qform(coded=False)

Return 4x4 affine matrix from qform parameters in header

Parameters:

coded : bool, optional

If True, return {affine or None}, and qform code. If False, just return affine. {affine or None} means, return None if qform code == 0, and affine otherwise.

Returns:

affine : None or (4,4) ndarray

If coded is False, always return affine reconstructed from qform quaternion. If coded is True, return None if qform code is 0, else return the affine.

code : int

Qform code. Only returned if coded is True.

get_qform_quaternion()

Compute quaternion from b, c, d of quaternion

Fills a value by assuming this is a unit quaternion

get_sform(coded=False)

Return 4x4 affine matrix from sform parameters in header

Parameters:

coded : bool, optional

If True, return {affine or None}, and sform code. If False, just return affine. {affine or None} means, return None if sform code == 0, and affine otherwise.

Returns:

affine : None or (4,4) ndarray

If coded is False, always return affine from sform fields. If coded is True, return None if sform code is 0, else return the affine.

code : int

Sform code. Only returned if coded is True.

get_slice_duration()

Get slice duration

Returns:

slice_duration : float

time to acquire one slice

Notes

The NIfTI1 spec appears to require the slice dimension to be defined for slice_duration to have meaning.

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_dim_info(slice=2)
>>> hdr.set_slice_duration(0.3)
>>> print("%0.1f" % hdr.get_slice_duration())
0.3
get_slice_times()

Get slice times from slice timing information

Returns:

slice_times : tuple

Times of acquisition of slices, where 0 is the beginning of the acquisition, ordered by position in file. nifti allows slices at the top and bottom of the volume to be excluded from the standard slice timing specification, and calls these “padding slices”. We give padding slices None as a time of acquisition

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_dim_info(slice=2)
>>> hdr.set_data_shape((1, 1, 7))
>>> hdr.set_slice_duration(0.1)
>>> hdr['slice_code'] = slice_order_codes['sequential increasing']
>>> slice_times = hdr.get_slice_times()
>>> np.allclose(slice_times, [0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6])
True
get_slope_inter()

Get data scaling (slope) and DC offset (intercept) from header data

Returns:

slope : None or float

scaling (slope). None if there is no valid scaling from these fields

inter : None or float

offset (intercept). None if there is no valid scaling or if offset is not finite.

Examples

>>> hdr = Nifti1Header()
>>> hdr.get_slope_inter()
(1.0, 0.0)
>>> hdr['scl_slope'] = 0
>>> hdr.get_slope_inter()
(None, None)
>>> hdr['scl_slope'] = np.nan
>>> hdr.get_slope_inter()
(None, None)
>>> hdr['scl_slope'] = 1
>>> hdr['scl_inter'] = 1
>>> hdr.get_slope_inter()
(1.0, 1.0)
>>> hdr['scl_inter'] = np.inf
>>> hdr.get_slope_inter() 
Traceback (most recent call last):
    ...
HeaderDataError: Valid slope but invalid intercept inf
Traceback (most recent call last):
    ...
HeaderDataError: Valid slope but invalid intercept inf
get_xyzt_units()
has_data_intercept = True
has_data_slope = True
is_single = True
pair_magic = 'ni1'
pair_vox_offset = 0
quaternion_threshold = -3.5762786865234375e-07
set_data_shape(shape)

Set shape of data

If ndims == len(shape) then we set zooms for dimensions higher than ndims to 1.0

Parameters:

shape : sequence

sequence of integers specifying data array shape

Notes

Applies freesurfer hack for large vectors described in https://github.com/nipy/nibabel/issues/100 and https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/save_nifti.m?spec=svn5022&r=5022#77

Allows for freesurfer hack for 7th order icosahedron surface described in https://github.com/nipy/nibabel/issues/309 https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/load_nifti.m?r=8776#86 https://code.google.com/p/fieldtrip/source/browse/trunk/external/freesurfer/save_nifti.m?r=8776#50

set_dim_info(freq=None, phase=None, slice=None)

Sets nifti MRI slice etc dimension information

Parameters:

freq : {None, 0, 1, 2}

axis of data array referring to frequency encoding

phase : {None, 0, 1, 2}

axis of data array referring to phase encoding

slice : {None, 0, 1, 2}

axis of data array referring to slice encoding

``None`` means the axis is not specified. :

Notes

This is stored in one byte in the header

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_dim_info(1, 2, 0)
>>> hdr.get_dim_info()
(1, 2, 0)
>>> hdr.set_dim_info(freq=1, phase=2, slice=0)
>>> hdr.get_dim_info()
(1, 2, 0)
>>> hdr.set_dim_info()
>>> hdr.get_dim_info()
(None, None, None)
>>> hdr.set_dim_info(freq=1, phase=None, slice=0)
>>> hdr.get_dim_info()
(1, None, 0)
set_intent(code, params=(), name='')

Set the intent code, parameters and name

If parameters are not specified, assumed to be all zero. Each intent code has a set number of parameters associated. If you specify any parameters, then it will need to be the correct number (e.g the “f test” intent requires 2). However, parameters can also be set in the file data, so we also allow not setting any parameters (empty parameter tuple).

Parameters:

code : integer or string

code specifying nifti intent

params : list, tuple of scalars

parameters relating to intent (see intent_codes) defaults to (). Unspecified parameters are set to 0.0

name : string

intent name (description). Defaults to ‘’

Returns:

None :

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_intent(0) # unknown code
>>> hdr.set_intent('z score')
>>> hdr.get_intent()
('z score', (), '')
>>> hdr.get_intent('code')
(5, (), '')
>>> hdr.set_intent('t test', (10,), name='some score')
>>> hdr.get_intent()
('t test', (10.0,), 'some score')
>>> hdr.set_intent('f test', (2, 10), name='another score')
>>> hdr.get_intent()
('f test', (2.0, 10.0), 'another score')
>>> hdr.set_intent('f test')
>>> hdr.get_intent()
('f test', (0.0, 0.0), '')
set_qform(affine, code=None, strip_shears=True)

Set qform header values from 4x4 affine

Parameters:

affine : None or 4x4 array

affine transform to write into sform. If None, only set code.

code : None, string or integer, optional

String or integer giving meaning of transform in affine. The default is None. If code is None, then:

  • If affine is None, code-> 0
  • If affine not None and existing qform code in header == 0, code-> 2 (aligned)
  • If affine not None and existing qform code in header != 0, code-> existing qform code in header

strip_shears : bool, optional

Whether to strip shears in affine. If True, shears will be silently stripped. If False, the presence of shears will raise a HeaderDataError

Notes

The qform transform only encodes translations, rotations and zooms. If there are shear components to the affine transform, and strip_shears is True (the default), the written qform gives the closest approximation where the rotation matrix is orthogonal. This is to allow quaternion representation. The orthogonal representation enforces orthogonal axes.

Examples

>>> hdr = Nifti1Header()
>>> int(hdr['qform_code']) # gives 0 - unknown
0
>>> affine = np.diag([1,2,3,1])
>>> np.all(hdr.get_qform() == affine)
False
>>> hdr.set_qform(affine)
>>> np.all(hdr.get_qform() == affine)
True
>>> int(hdr['qform_code']) # gives 2 - aligned
2
>>> hdr.set_qform(affine, code='talairach')
>>> int(hdr['qform_code'])
3
>>> hdr.set_qform(affine, code=None)
>>> int(hdr['qform_code'])
3
>>> hdr.set_qform(affine, code='scanner')
>>> int(hdr['qform_code'])
1
>>> hdr.set_qform(None)
>>> int(hdr['qform_code'])
0
set_sform(affine, code=None)

Set sform transform from 4x4 affine

Parameters:

affine : None or 4x4 array

affine transform to write into sform. If None, only set code

code : None, string or integer, optional

String or integer giving meaning of transform in affine. The default is None. If code is None, then:

  • If affine is None, code-> 0
  • If affine not None and existing sform code in header == 0, code-> 2 (aligned)
  • If affine not None and existing sform code in header != 0, code-> existing sform code in header

Examples

>>> hdr = Nifti1Header()
>>> int(hdr['sform_code']) # gives 0 - unknown
0
>>> affine = np.diag([1,2,3,1])
>>> np.all(hdr.get_sform() == affine)
False
>>> hdr.set_sform(affine)
>>> np.all(hdr.get_sform() == affine)
True
>>> int(hdr['sform_code']) # gives 2 - aligned
2
>>> hdr.set_sform(affine, code='talairach')
>>> int(hdr['sform_code'])
3
>>> hdr.set_sform(affine, code=None)
>>> int(hdr['sform_code'])
3
>>> hdr.set_sform(affine, code='scanner')
>>> int(hdr['sform_code'])
1
>>> hdr.set_sform(None)
>>> int(hdr['sform_code'])
0
set_slice_duration(duration)

Set slice duration

Parameters:

duration : scalar

time to acquire one slice

Examples

See get_slice_duration

set_slice_times(slice_times)

Set slice times into hdr

Parameters:

slice_times : tuple

tuple of slice times, one value per slice tuple can include None to indicate no slice time for that slice

Examples

>>> hdr = Nifti1Header()
>>> hdr.set_dim_info(slice=2)
>>> hdr.set_data_shape([1, 1, 7])
>>> hdr.set_slice_duration(0.1)
>>> times = [None, 0.2, 0.4, 0.1, 0.3, 0.0, None]
>>> hdr.set_slice_times(times)
>>> hdr.get_value_label('slice_code')
'alternating decreasing'
>>> int(hdr['slice_start'])
1
>>> int(hdr['slice_end'])
5
set_slope_inter(slope, inter=None)

Set slope and / or intercept into header

Set slope and intercept for image data, such that, if the image data is arr, then the scaled image data will be (arr * slope) + inter

(slope, inter) of (NaN, NaN) is a signal to a containing image to set slope, inter automatically on write.

Parameters:

slope : None or float

If None, implies slope of NaN. If slope is None or NaN then inter should be None or NaN. Values of 0, Inf or -Inf raise HeaderDataError

inter : None or float, optional

Intercept. If None, implies inter of NaN. If slope is None or NaN then inter should be None or NaN. Values of Inf or -Inf raise HeaderDataError

set_xyzt_units(xyz=None, t=None)
single_magic = 'n+1'
single_vox_offset = 352
template_dtype = dtype([('sizeof_hdr', '<i4'), ('data_type', 'S10'), ('db_name', 'S18'), ('extents', '<i4'), ('session_error', '<i2'), ('regular', 'S1'), ('dim_info', 'u1'), ('dim', '<i2', (8,)), ('intent_p1', '<f4'), ('intent_p2', '<f4'), ('intent_p3', '<f4'), ('intent_code', '<i2'), ('datatype', '<i2'), ('bitpix', '<i2'), ('slice_start', '<i2'), ('pixdim', '<f4', (8,)), ('vox_offset', '<f4'), ('scl_slope', '<f4'), ('scl_inter', '<f4'), ('slice_end', '<i2'), ('slice_code', 'u1'), ('xyzt_units', 'u1'), ('cal_max', '<f4'), ('cal_min', '<f4'), ('slice_duration', '<f4'), ('toffset', '<f4'), ('glmax', '<i4'), ('glmin', '<i4'), ('descrip', 'S80'), ('aux_file', 'S24'), ('qform_code', '<i2'), ('sform_code', '<i2'), ('quatern_b', '<f4'), ('quatern_c', '<f4'), ('quatern_d', '<f4'), ('qoffset_x', '<f4'), ('qoffset_y', '<f4'), ('qoffset_z', '<f4'), ('srow_x', '<f4', (4,)), ('srow_y', '<f4', (4,)), ('srow_z', '<f4', (4,)), ('intent_name', 'S16'), ('magic', 'S4')])
write_to(fileobj)

Nifti1Image

class nibabel.nifti1.Nifti1Image(dataobj, affine, header=None, extra=None, file_map=None)

Bases: nibabel.nifti1.Nifti1Pair

Class for single file NIfTI1 format image

__init__(dataobj, affine, header=None, extra=None, file_map=None)
files_types = (('image', '.nii'),)
header_class

alias of Nifti1Header

update_header()

Harmonize header with image data and affine

Nifti1Pair

class nibabel.nifti1.Nifti1Pair(dataobj, affine, header=None, extra=None, file_map=None)

Bases: nibabel.analyze.AnalyzeImage

Class for NIfTI1 format image, header pair

__init__(dataobj, affine, header=None, extra=None, file_map=None)
get_qform(coded=False)

Return 4x4 affine matrix from qform parameters in header

Parameters:

coded : bool, optional

If True, return {affine or None}, and qform code. If False, just return affine. {affine or None} means, return None if qform code == 0, and affine otherwise.

Returns:

affine : None or (4,4) ndarray

If coded is False, always return affine reconstructed from qform quaternion. If coded is True, return None if qform code is 0, else return the affine.

code : int

Qform code. Only returned if coded is True.

See also

set_qform, get_sform

get_sform(coded=False)

Return 4x4 affine matrix from sform parameters in header

Parameters:

coded : bool, optional

If True, return {affine or None}, and sform code. If False, just return affine. {affine or None} means, return None if sform code == 0, and affine otherwise.

Returns:

affine : None or (4,4) ndarray

If coded is False, always return affine from sform fields. If coded is True, return None if sform code is 0, else return the affine.

code : int

Sform code. Only returned if coded is True.

See also

set_sform, get_qform

header_class

alias of Nifti1PairHeader

set_qform(affine, code=None, strip_shears=True, **kwargs)

Set qform header values from 4x4 affine

Parameters:

affine : None or 4x4 array

affine transform to write into sform. If None, only set code.

code : None, string or integer

String or integer giving meaning of transform in affine. The default is None. If code is None, then:

  • If affine is None, code-> 0
  • If affine not None and existing qform code in header == 0, code-> 2 (aligned)
  • If affine not None and existing qform code in header != 0, code-> existing qform code in header

strip_shears : bool, optional

Whether to strip shears in affine. If True, shears will be silently stripped. If False, the presence of shears will raise a HeaderDataError

update_affine : bool, optional

Whether to update the image affine from the header best affine after setting the qform. Must be keyword argument (because of different position in set_qform). Default is True

See also

get_qform, set_sform

Examples

>>> data = np.arange(24).reshape((2,3,4))
>>> aff = np.diag([2, 3, 4, 1])
>>> img = Nifti1Pair(data, aff)
>>> img.get_qform()
array([[ 2.,  0.,  0.,  0.],
       [ 0.,  3.,  0.,  0.],
       [ 0.,  0.,  4.,  0.],
       [ 0.,  0.,  0.,  1.]])
>>> img.get_qform(coded=True)
(None, 0)
>>> aff2 = np.diag([3, 4, 5, 1])
>>> img.set_qform(aff2, 'talairach')
>>> qaff, code = img.get_qform(coded=True)
>>> np.all(qaff == aff2)
True
>>> int(code)
3
set_sform(affine, code=None, **kwargs)

Set sform transform from 4x4 affine

Parameters:

affine : None or 4x4 array

affine transform to write into sform. If None, only set code

code : None, string or integer

String or integer giving meaning of transform in affine. The default is None. If code is None, then:

  • If affine is None, code-> 0
  • If affine not None and existing sform code in header == 0, code-> 2 (aligned)
  • If affine not None and existing sform code in header != 0, code-> existing sform code in header

update_affine : bool, optional

Whether to update the image affine from the header best affine after setting the qform. Must be keyword argument (because of different position in set_qform). Default is True

See also

get_sform, set_qform

Examples

>>> data = np.arange(24).reshape((2,3,4))
>>> aff = np.diag([2, 3, 4, 1])
>>> img = Nifti1Pair(data, aff)
>>> img.get_sform()
array([[ 2.,  0.,  0.,  0.],
       [ 0.,  3.,  0.,  0.],
       [ 0.,  0.,  4.,  0.],
       [ 0.,  0.,  0.,  1.]])
>>> saff, code = img.get_sform(coded=True)
>>> saff
array([[ 2.,  0.,  0.,  0.],
       [ 0.,  3.,  0.,  0.],
       [ 0.,  0.,  4.,  0.],
       [ 0.,  0.,  0.,  1.]])
>>> int(code)
2
>>> aff2 = np.diag([3, 4, 5, 1])
>>> img.set_sform(aff2, 'talairach')
>>> saff, code = img.get_sform(coded=True)
>>> np.all(saff == aff2)
True
>>> int(code)
3
update_header()

Harmonize header with image data and affine

See AnalyzeImage.update_header for more examples

Examples

>>> data = np.zeros((2,3,4))
>>> affine = np.diag([1.0,2.0,3.0,1.0])
>>> img = Nifti1Image(data, affine)
>>> hdr = img.header
>>> np.all(hdr.get_qform() == affine)
True
>>> np.all(hdr.get_sform() == affine)
True

Nifti1PairHeader

class nibabel.nifti1.Nifti1PairHeader(binaryblock=None, endianness=None, check=True, extensions=())

Bases: nibabel.nifti1.Nifti1Header

Class for NIfTI1 pair header

Initialize header from binary data block and extensions

__init__(binaryblock=None, endianness=None, check=True, extensions=())

Initialize header from binary data block and extensions

is_single = False

load

nibabel.nifti1.load(filename)

Load NIfTI1 single or pair from filename

Parameters:

filename : str

filename of image to be loaded

Returns:

img : Nifti1Image or Nifti1Pair

NIfTI1 single or pair image instance

Raises:

ImageFileError :

if filename doesn’t look like NIfTI1;

IOError :

if filename does not exist.

save

nibabel.nifti1.save(img, filename)

Save NIfTI1 single or pair to filename

Parameters:

filename : str

filename to which to save image