Newer
Older
Paul McCarthy
committed
#!/usr/bin/env python
#
# image.py - Provides the :class:`Image` class, for representing 3D/4D NIFTI
# images.
Paul McCarthy
committed
#
# Author: Paul McCarthy <pauldmccarthy@gmail.com>
#
Paul McCarthy
committed
"""This module provides the :class:`Nifti` and :class:`Image` classes, for
representing NIFTI1 and NIFTI2 images. The ``nibabel`` package is used
It is very easy to load a NIFTI image::
from fsl.data.image import Image
myimg = Image('MNI152_T1_2mm.nii.gz')
A handful of other functions are also provided for working with image files
and file names:
.. autosummary::
:nosignatures:
addExt
getExt
removeExt
defaultExt
Paul McCarthy
committed
import os
import os.path as op
import itertools as it

Paul McCarthy
committed
import json
Paul McCarthy
committed
import string
import logging
import warnings
Paul McCarthy
committed
import six
import numpy as np
import nibabel as nib
import nibabel.fileslice as fileslice
import fsl.utils.meta as meta
import fsl.transform.affine as affine
import fsl.utils.notifier as notifier
import fsl.utils.memoize as memoize
import fsl.utils.path as fslpath
import fsl.utils.deprecated as deprecated
import fsl.utils.bids as fslbids
import fsl.data.constants as constants
import fsl.data.imagewrapper as imagewrapper
Paul McCarthy
committed
Paul McCarthy
committed
log = logging.getLogger(__name__)
ALLOWED_EXTENSIONS = ['.nii.gz', '.nii', '.img', '.hdr', '.img.gz', '.hdr.gz']
"""The file extensions which we understand. This list is used as the default
if the ``allowedExts`` parameter is not passed to any of the ``*Ext``
functions, or the :func:`looksLikeImage` function.
"""
EXTENSION_DESCRIPTIONS = ['Compressed NIFTI images',
'NIFTI images',
'ANALYZE75 images',
'NIFTI/ANALYZE75 headers',
'Compressed NIFTI/ANALYZE75 images',
'Compressed NIFTI/ANALYZE75 headers']
"""Descriptions for each of the extensions in :data:`ALLOWED_EXTENSIONS`. """
FILE_GROUPS = [('.hdr', '.img'),
('.hdr.gz', '.img.gz')]
"""File suffix groups used by :func:`addExt` to resolve file path
ambiguities - see :func:`fsl.utils.path.addExt`.
"""
PathError = fslpath.PathError
"""Error raised by :mod:`fsl.utils.path` functions when an error occurs.
Made available in this module for convenience.
"""
class Nifti(notifier.Notifier, meta.Meta):
Paul McCarthy
committed
"""The ``Nifti`` class is intended to be used as a base class for
things which either are, or are associated with, a NIFTI image.
The ``Nifti`` class is intended to represent information stored in
the header of a NIFTI file - if you want to load the data from
a file, use the :class:`Image` class instead.
Paul McCarthy
committed
When a ``Nifti`` instance is created, it adds the following attributes
to itself:
================= ====================================================
``header`` The :mod:`nibabel` NIFTI1/NIFTI2/Analyze header
object.
``shape`` A list/tuple containing the number of voxels along
each image dimension.
``pixdim`` A list/tuple containing the length of one voxel
along each image dimension.
``voxToWorldMat`` A 4*4 array specifying the affine transformation
for transforming voxel coordinates into real world
coordinates.
``worldToVoxMat`` A 4*4 array specifying the affine transformation
for transforming real world coordinates into voxel
coordinates.
``intent`` The NIFTI intent code specified in the header (or
:attr:`.constants.NIFTI_INTENT_NONE` for Analyze
images).
================= ====================================================
The ``header`` field may either be a ``nifti1``, ``nifti2``, or
``analyze`` header object. Make sure to take this into account if you are
writing code that should work with all three. Use the :meth:`niftiVersion`
property if you need to know what type of image you are dealing with.
The ``shape`` attribute may not precisely match the image shape as
reported in the NIFTI header, because trailing dimensions of size 1 are
squeezed out. See the :meth:`__determineShape` and :meth:`mapIndices`
methods.
**Affine transformations**
The ``Nifti`` class is aware of three coordinate systems:
- The ``voxel`` coordinate system, used to access image data
- The ``world`` coordinate system, where voxel coordinates are
transformed into a millimetre coordinate system, defined by the
``sform`` and/or ``qform`` elements of the NIFTI header.
- The ``fsl`` coordinate system, where voxel coordinates are scaled by
the ``pixdim`` values in the NIFTI header, and the X axis is inverted
if the voxel-to-world affine has a positive determinant.
The :meth:`getAffine` method is a simple means of acquiring an affine
which will transform between any of these coordinate systems.
See `here <http://fsl.fmrib.ox.ac.uk/fsl/fslwiki/FLIRT/FAQ#What_is_the_format_of_the_matrix_used_by_FLIRT.2C_and_how_does_it_relate_to_the_transformation_parameters.3F>`_
for more details on the ``fsl`` coordinate system..
The ``Nifti`` class follows the same process as ``nibabel`` in determining
the ``voxel`` to ``world`` affine (see
http://nipy.org/nibabel/nifti_images.html#the-nifti-affines):
1. If ``sform_code != 0`` ("unknown") use the sform affine; else
2. If ``qform_code != 0`` ("unknown") use the qform affine; else
3. Use the fall-back affine.
However, the *fall-back* affine used by the ``Nifti`` class differs to
that used by ``nibabel``. In ``nibabel``, the origin (world coordinates
(0, 0, 0)) is set to the centre of the image. Here in the ``Nifti``
class, we set the world coordinate orign to be the corner of the image,
i.e. the corner of voxel (0, 0, 0).
You may change the ``voxToWorldMat`` of a ``Nifti`` instance (unless it
is an Analyze image). When you do so:
- Only the ``sform`` of the underlying ``Nifti1Header`` object is changed
- The ``qform`` is not modified.
- If the ``sform_code`` was previously set to ``NIFTI_XFORM_UNKNOWN``,
it is changed to ``NIFTI_XFORM_ALIGNED_ANAT``. Otherwise, the
``sform_code`` is not modified.
**ANALYZE support**
A ``Nifti`` instance expects to be passed either a
``nibabel.nifti1.Nifti1Header`` or a ``nibabel.nifti2.Nifti2Header``, but
can also encapsulate a ``nibabel.analyze.AnalyzeHeader``. In this case:
- The image voxel orientation is assumed to be R->L, P->A, I->S.
- The affine will be set to a diagonal matrix with the header pixdims as
its elements (with the X pixdim negated), and an offset specified by
the ANALYZE ``origin`` fields. Construction of the affine is handled
by ``nibabel``.
- The :meth:`niftiVersion` method will return ``0``.
- The :meth:`getXFormCode` method will return
:attr:`.constants.NIFTI_XFORM_ANALYZE`.
**Metadata**
The ``Image`` class inherits from the :class:`.Meta` class - its methods
can be used to store and query any meta-data associated with the image.
**Notification**
The ``Nifti`` class implements the :class:`.Notifier` interface -
listeners may register to be notified on the following topics:
=============== ========================================================
``'transform'`` The affine transformation matrix has changed. This topic
will occur when the :meth:`voxToWorldMat` is changed.
``'header'`` A header field has changed. This will occur when the
:meth:`intent` is changed.
=============== ========================================================
def __init__(self, header):
Paul McCarthy
committed
"""Create a ``Nifti`` object.
:arg header: A :class:`nibabel.nifti1.Nifti1Header`,
:class:`nibabel.nifti2.Nifti2Header`, or
``nibabel.analyze.AnalyzeHeader`` to be used as the
image header.
Paul McCarthy
committed
"""
# and Nifti1Header a sub-class of AnalyzeHeader,
# so we only need to test for the latter.
if not isinstance(header, nib.analyze.AnalyzeHeader):
raise ValueError('Unrecognised header: {}'.format(header))

Paul McCarthy
committed
origShape, shape, pixdim = Nifti.determineShape(header)
voxToWorldMat = Nifti.determineAffine(header)
affines, isneuro = Nifti.generateAffines(voxToWorldMat,
shape,
pixdim)
Paul McCarthy
committed
self.header = header
self.__shape = shape
self.__origShape = origShape
self.__pixdim = pixdim
self.__affines = affines
self.__isNeurological = isneuro

Paul McCarthy
committed
@staticmethod
def determineShape(header):
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
"""This method is called by :meth:`__init__`. It figures out the actual
shape of the image data, and the zooms/pixdims for each data axis. Any
empty trailing dimensions are squeezed, but the returned shape is
guaranteed to be at least 3 dimensions. Returns:
- A sequence/tuple containing the image shape, as reported in the
header.
- A sequence/tuple containing the effective image shape.
- A sequence/tuple containing the zooms/pixdims.
"""
# The canonicalShape function figures out
# the data shape that we should use.
origShape = list(header.get_data_shape())
shape = canonicalShape(origShape)
pixdims = list(header.get_zooms())
# if get_zooms() doesn't return at
# least len(shape) values, we'll
# fall back to the pixdim field in
# the header.
if len(pixdims) < len(shape):
pixdims = header['pixdim'][1:]
pixdims = pixdims[:len(shape)]

Paul McCarthy
committed
# should never happen, but if we only
# have zoom values for the original
# (< 3D) shape, pad them with 1s.
if len(pixdims) < len(shape):
pixdims = pixdims + [1] * (len(shape) - len(pixdims))
return origShape, shape, pixdims

Paul McCarthy
committed
@staticmethod
def determineAffine(header):
"""Called by :meth:`__init__`. Figures out the voxel-to-world
coordinate transformation matrix that is associated with this
Paul McCarthy
committed
``Nifti`` instance.
# We have to treat FSL/FNIRT images
# specially, as FNIRT clobbers the
# sform section of the NIFTI header
intent = header.get('intent_code', -1)
qform = header.get('qform_code', -1)
sform = header.get('sform_code', -1)

Paul McCarthy
committed
# FNIRT non-linear coefficient files
# clobber the sform/qform/intent
# and pixdims of the nifti header,
# so we can't correctly place it in
# the world coordinate system. See
# $FSLDIR/src/fnirt/fnirt_file_writer.cpp
# and fsl.transform.nonlinear for more
# details.
if intent in (constants.FSL_DCT_COEFFICIENTS,
constants.FSL_CUBIC_SPLINE_COEFFICIENTS,
constants.FSL_QUADRATIC_SPLINE_COEFFICIENTS,
constants.FSL_TOPUP_CUBIC_SPLINE_COEFFICIENTS,
constants.FSL_TOPUP_QUADRATIC_SPLINE_COEFFICIENTS):

Paul McCarthy
committed
log.debug('FNIRT coefficient field detected - generating affine')
# Knot spacing is stored in the pixdims
# (specified in terms of reference image
# voxels), and reference image pixdims
# are stored as intent code parameters.
# If we combine the two, we can at least
# get the shape/size of the coefficient
# field about right

Paul McCarthy
committed
knotpix = header.get_zooms()[:3]

Paul McCarthy
committed
refpix = (header.get('intent_p1', 1) or 1,
header.get('intent_p2', 1) or 1,
header.get('intent_p3', 1) or 1)
voxToWorldMat = affine.concat(
affine.scaleOffsetXform(refpix, 0),
affine.scaleOffsetXform(knotpix, 0))
# If the qform or sform codes are unknown,
# then we can't assume that the transform
# matrices are valid. So we fall back to a
# pixdim scaling matrix.
#
# n.b. For images like this, nibabel returns
# a scaling matrix where the centre voxel
# corresponds to world location (0, 0, 0).
Paul McCarthy
committed
# This goes against the NIFTI spec - it
# should just be a straight scaling matrix.
elif qform == 0 and sform == 0:
pixdims = header.get_zooms()
voxToWorldMat = affine.scaleOffsetXform(pixdims, 0)
# Otherwise we let nibabel decide
# which transform to use.
voxToWorldMat = np.array(header.get_best_affine())
return voxToWorldMat
Paul McCarthy
committed

Paul McCarthy
committed
@staticmethod
def generateAffines(voxToWorldMat, shape, pixdim):
"""Called by :meth:`__init__`, and the :meth:`voxToWorldMat` setter.
Generates and returns a dictionary containing affine transformations
between the ``voxel``, ``fsl``, and ``world`` coordinate
systems. These affines are accessible via the :meth:`getAffine`
method.
:arg voxToWorldMat: The voxel-to-world affine transformation

Paul McCarthy
committed
:arg shape: Image shape (number of voxels along each dimension
:arg pixdim: Image pixdims (size of one voxel along each
dimension)
:returns: A tuple containing:
- a dictionary of affine transformations between
each pair of coordinate systems
- ``True`` if the image is to be considered
"neurological", ``False`` otherwise - see the
:meth:`isNeurological` method.
import numpy.linalg as npla
affines = {}

Paul McCarthy
committed
shape = list(shape[ :3])
pixdim = list(pixdim[:3])
voxToScaledVoxMat = np.diag(pixdim + [1.0])
isneuro = npla.det(voxToWorldMat) > 0
if isneuro:
x = (shape[0] - 1) * pixdim[0]
flip = affine.scaleOffsetXform([-1, 1, 1],
[ x, 0, 0])
voxToScaledVoxMat = affine.concat(flip, voxToScaledVoxMat)
affines['fsl', 'fsl'] = np.eye(4)
affines['voxel', 'voxel'] = np.eye(4)
affines['world', 'world'] = np.eye(4)
affines['voxel', 'world'] = voxToWorldMat
affines['world', 'voxel'] = affine.invert(voxToWorldMat)
affines['voxel', 'fsl'] = voxToScaledVoxMat
affines['fsl', 'voxel'] = affine.invert(voxToScaledVoxMat)
affines['fsl', 'world'] = affine.concat(affines['voxel', 'world'],
affines['fsl', 'voxel'])
affines['world', 'fsl'] = affine.concat(affines['voxel', 'fsl'],
affines['world', 'voxel'])
return affines, isneuro
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
@staticmethod
def identifyAffine(image, xform, from_=None, to=None):
"""Attempt to identify the source or destination space for the given
affine.
``xform`` is assumed to be an affine transformation which can be used
to transform coordinates between two coordinate systems associated with
``image``.
If one of ``from_`` or ``to`` is provided, the other will be derived.
If neither are provided, both will be derived. See the
:meth:`.Nifti.getAffine` method for details on the valild values that
``from_`` and ``to`` may take.
:arg image: :class:`.Nifti` instance associated with the affine.
:arg xform: ``(4, 4)`` ``numpy`` array encoding an affine
transformation
:arg from_: Label specifying the coordinate system which ``xform``
takes as input
:arg to: Label specifying the coordinate system which ``xform``
produces as output
:returns: A tuple containing:
- A label for the ``from_`` coordinate system
- A label for the ``to`` coordinate system
"""
if (from_ is not None) and (to is not None):
return from_, to
if from_ is not None: froms = [from_]
else: froms = ['voxel', 'fsl', 'world']
if to is not None: tos = [to]
else: tos = ['voxel', 'fsl', 'world']
for from_, to in it.product(froms, tos):
candidate = image.getAffine(from_, to)
if np.all(np.isclose(candidate, xform)):
return from_, to
raise ValueError('Could not identify affine')
Paul McCarthy
committed
def strval(self, key):
"""Returns the specified NIFTI header field, converted to a python
string, correctly null-terminated, and with non-printable characters
removed.
This method is used to sanitise some NIFTI header fields. The default
Python behaviour for converting a sequence of bytes to a string is to
strip all termination characters (bytes with value of ``0x00``) from
the end of the sequence.
This default behaviour does not handle the case where a sequence of
bytes which did contain a long string is subsequently overwritten with
a shorter string - the short string will be terminated, but that
termination character will be followed by the remainder of the
original string.
"""
val = self.header[key]
try: val = bytes(val).partition(b'\0')[0]
except Exception: val = bytes(val)
Paul McCarthy
committed
val = val.decode('ascii')
return ''.join([c for c in val if c in string.printable]).strip()
Paul McCarthy
committed
@property
def niftiVersion(self):
"""Returns the NIFTI file version:
- ``0`` for ANALYZE
- ``1`` for NIFTI1
- ``2`` for NIFTI2
"""
# nib.Nifti2 is a subclass of Nifti1,
# and Nifti1 a subclass of Analyze,
# so we have to check in order
if isinstance(self.header, nib.nifti2.Nifti2Header): return 2
elif isinstance(self.header, nib.nifti1.Nifti1Header): return 1
elif isinstance(self.header, nib.analyze.AnalyzeHeader): return 0
else: raise RuntimeError('Unrecognised header: {}'.format(self.header))
@property
def shape(self):
"""Returns a tuple containing the image data shape. """
return tuple(self.__shape)
@property
def ndim(self):
"""Returns the number of dimensions in this image. This number may not
match the number of dimensions specified in the NIFTI header, as
trailing dimensions of length 1 are ignored. But it is guaranteed to be
at least 3.
"""
return len(self.__shape)
@property
def pixdim(self):
"""Returns a tuple containing the image pixdims (voxel sizes)."""
return tuple(self.__pixdim)
@property
def intent(self):
"""Returns the NIFTI intent code of this image. """
return self.header.get('intent_code', constants.NIFTI_INTENT_NONE)
@intent.setter
def intent(self, val):
"""Sets the NIFTI intent code of this image. """
# analyze has no intent
if (self.niftiVersion > 0) and (val != self.intent):
self.header.set_intent(val, allow_unknown=True)
self.notify(topic='header')
@property
def xyzUnits(self):
"""Returns the NIFTI XYZ dimension unit code. """
# analyze images have no unit field
if self.niftiVersion == 0:
return constants.NIFTI_UNITS_MM
# The nibabel get_xyzt_units returns labels,
# but we want the NIFTI codes. So we use
# the (undocumented) nifti1.unit_codes field
# to convert back to the raw codes.
units = self.header.get_xyzt_units()[0]
units = nib.nifti1.unit_codes[units]
return units
def timeUnits(self):
"""Returns the NIFTI time dimension unit code. """
# analyze images have no unit field
if self.niftiVersion == 0:
return constants.NIFTI_UNITS_SEC
# See xyzUnits
units = self.header.get_xyzt_units()[1]
units = nib.nifti1.unit_codes[units]
return units
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
def getAffine(self, from_, to):
"""Return an affine transformation which can be used to transform
coordinates from ``from_`` to ``to``.
Valid values for the ``from_`` and ``to`` arguments are:
- ``'voxel'``: The voxel coordinate system
- ``'world'``: The world coordinate system, as defined by the image
sform/qform
- ``'fsl'``: The FSL coordinate system (scaled voxels, with a
left-right flip if the sform/qform has a positive determinant)
:arg from_: Source coordinate system
:arg to: Destination coordinate system
:returns: A ``numpy`` array of shape ``(4, 4)``
"""
from_ = from_.lower()
to = to .lower()
if from_ not in ('voxel', 'fsl', 'world') or \
to not in ('voxel', 'fsl', 'world'):
raise ValueError('Invalid source/reference spaces: '
'{} -> {}'.format(from_, to))
return np.copy(self.__affines[from_, to])
@property
def worldToVoxMat(self):
"""Returns a ``numpy`` array of shape ``(4, 4)`` containing an
affine transformation from world coordinates to voxel coordinates.
"""
return self.getAffine('world', 'voxel')
@property
def voxToWorldMat(self):
"""Returns a ``numpy`` array of shape ``(4, 4)`` containing an
affine transformation from voxel coordinates to world coordinates.
"""
return self.getAffine('voxel', 'world')
@voxToWorldMat.setter
def voxToWorldMat(self, xform):
"""Update the ``voxToWorldMat``. The ``worldToVoxMat`` value is also
updated. This will result in notification on the ``'transform'``
topic.
"""
# Can't do much with
# an analyze image
if self.niftiVersion == 0:
raise Exception('voxToWorldMat cannot be '
'changed for an ANALYZE image')
header = self.header
sformCode = int(header['sform_code'])
if sformCode == constants.NIFTI_XFORM_UNKNOWN:
sformCode = constants.NIFTI_XFORM_ALIGNED_ANAT
header.set_sform(xform, code=sformCode)

Paul McCarthy
committed
affines, isneuro = Nifti.generateAffines(xform,
self.shape,
self.pixdim)
self.__affines = affines
self.__isNeurological = isneuro
'%s\nsform: %s\nqform: %s',
header.get_zooms(),
header.get_sform(),
header.get_qform())
self.notify(topic='transform')
@property
def voxToScaledVoxMat(self):
"""Returns a transformation matrix which transforms from voxel
coordinates into scaled voxel coordinates, with a left-right flip
if the image appears to be stored in neurological order.
See http://fsl.fmrib.ox.ac.uk/fsl/fslwiki/FLIRT/FAQ#What_is_the\
_format_of_the_matrix_used_by_FLIRT.2C_and_how_does_it_relate_to\
_the_transformation_parameters.3F
"""
return self.getAffine('voxel', 'fsl')
@property
def scaledVoxToVoxMat(self):
"""Returns a transformation matrix which transforms from scaled voxels
into voxels, the inverse of the :meth:`voxToScaledVoxMat` transform.
"""
return self.getAffine('fsl', 'voxel')
def mapIndices(self, sliceobj):
"""Adjusts the given slice object so that it may be used to index the
See the :meth:`__determineShape` method.
:arg sliceobj: Something that can be used to slice a
multi-dimensional array, e.g. ``arr[sliceobj]``.
"""
# How convenient - nibabel has a function
# that does the dirty work for us.
return fileslice.canonical_slicers(sliceobj, self.__origShape)
def getXFormCode(self, code=None):
Paul McCarthy
committed
"""This method returns the code contained in the NIFTI header,
indicating the space to which the (transformed) image is oriented.
The ``code`` parameter may be either ``sform`` (the default) or
``qform`` in which case the corresponding matrix is used.
:returns: one of the following codes:
- :data:`~.constants.NIFTI_XFORM_UNKNOWN`
- :data:`~.constants.NIFTI_XFORM_SCANNER_ANAT`
- :data:`~.constants.NIFTI_XFORM_ALIGNED_ANAT`
- :data:`~.constants.NIFTI_XFORM_TALAIRACH`
- :data:`~.constants.NIFTI_XFORM_MNI_152`
- :data:`~.constants.NIFTI_XFORM_ANALYZE`
"""
if self.niftiVersion == 0:
return constants.NIFTI_XFORM_ANALYZE
if code == 'sform' : code = 'sform_code'
elif code == 'qform' : code = 'qform_code'
elif code is not None:
raise ValueError('code must be None, sform, or qform')
if code is not None:
code = self.header[code]
# If the caller did not specify
# a code, we check both. If the
# sform is present, we return it.
# Otherwise, if the qform is
# present, we return that.
sform_code = self.header['sform_code']
qform_code = self.header['qform_code']
if sform_code != constants.NIFTI_XFORM_UNKNOWN: code = sform_code
elif qform_code != constants.NIFTI_XFORM_UNKNOWN: code = qform_code
Paul McCarthy
committed
# Invalid values
if code not in range(5):
code = constants.NIFTI_XFORM_UNKNOWN
Paul McCarthy
committed
return int(code)
# TODO Check what has worse performance - hashing
# a 4x4 array (via memoizeMD5), or the call
# to aff2axcodes. I'm guessing the latter,
# but am not 100% sure.
@memoize.Instanceify(memoize.memoizeMD5)
def axisMapping(self, xform):
"""Returns the (approximate) correspondence of each axis in the source
coordinate system to the axes in the destination coordinate system,
where the source and destinations are defined by the given affine
transformation matrix.
"""
inaxes = [[-1, 1], [-2, 2], [-3, 3]]
return nib.orientations.aff2axcodes(xform, inaxes)
def isNeurological(self):
Paul McCarthy
committed
"""Returns ``True`` if it looks like this ``Nifti`` object has a
neurological voxel orientation, ``False`` otherwise. This test is
purely based on the determinant of the voxel-to-mm transformation
matrix - if it has a positive determinant, the image is assumed to
be in neurological orientation, otherwise it is assumed to be in
radiological orientation.
..warning:: This method will return ``True`` for images with an
unknown orientation (e.g. the ``sform_code`` and
``qform_code`` are both set to ``0``). Therefore, you
must check the orientation via the :meth:`getXFormCode`
before trusting the result of this method.
See http://fsl.fmrib.ox.ac.uk/fsl/fslwiki/FLIRT/FAQ#What_is_the\
_format_of_the_matrix_used_by_FLIRT.2C_and_how_does_it_relate_to\
_the_transformation_parameters.3F
return self.__isNeurological
def sameSpace(self, other):
"""Returns ``True`` if the ``other`` image (assumed to be a
:class:`Nifti` instance) has the same dimensions and is in the
same space as this image.
"""
return np.all(np.isclose(self .shape[:3],
other.shape[:3])) and \
np.all(np.isclose(self .pixdim[:3],
other.pixdim[:3])) and \
np.all(np.isclose(self .voxToWorldMat,
other.voxToWorldMat))
def getOrientation(self, axis, xform):
"""Returns a code representing the orientation of the specified
Paul McCarthy
committed
axis in the input coordinate system of the given transformation
matrix.
:arg xform: A transformation matrix which is assumed to transform
Paul McCarthy
committed
coordinates from some coordinate system (the one
which you want an orientation for) into the image
world coordinate system.
For example, if you pass in the voxel-to-world
transformation matrix, you will get an orientation
for axes in the voxel coordinate system.
This method returns one of the following values, indicating the
direction in which coordinates along the specified axis increase:
- :attr:`~.constants.ORIENT_L2R`: Left to right
- :attr:`~.constants.ORIENT_R2L`: Right to left
- :attr:`~.constants.ORIENT_A2P`: Anterior to posterior
- :attr:`~.constants.ORIENT_P2A`: Posterior to anterior
- :attr:`~.constants.ORIENT_I2S`: Inferior to superior
- :attr:`~.constants.ORIENT_S2I`: Superior to inferior
- :attr:`~.constants.ORIENT_UNKNOWN`: Orientation is unknown
The returned value is dictated by the XForm code contained in the
image file header (see the :meth:`getXFormCode` method). Basically, if
the XForm code is *unknown*, this method will return
``ORIENT_UNKNOWN`` for all axes. Otherwise, it is assumed that the
image is in RAS orientation (i.e. the X axis increases from left to
right, the Y axis increases from posterior to anterior, and the Z axis
increases from inferior to superior).
"""
if self.getXFormCode() == constants.NIFTI_XFORM_UNKNOWN:
return constants.ORIENT_UNKNOWN
code = nib.orientations.aff2axcodes(
((constants.ORIENT_R2L, constants.ORIENT_L2R),
(constants.ORIENT_A2P, constants.ORIENT_P2A),
(constants.ORIENT_S2I, constants.ORIENT_I2S)))[axis]
return code
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
def adjust(self, pixdim=None, shape=None, origin=None):
"""Return a new ``Nifti`` object with the specified ``pixdim``
or ``shape``. The affine is of the new ``Nifti`` is adjusted
accordingly.
Only one of ``pixdim`` or ``shape`` can be specified.
The default behaviour of the ``origin`` argument (``'centre'``) causes
the corner voxel of the output to have the same centre as the corner
voxel of the input. If the origin is ``'corner'``, an offset is
applied which effectively causes the origin corners of the voxel grids
of the input and output to be aligned.
:arg pixdim: New voxel dimensions
:arg shape: New image shape
:arg origin: Voxel grid alignment - either ``'centre'`` (the default)
or ``'corner'``
:returns: A new ``Nifti`` object based on this one, with adjusted
pixdims, shape and affine.
"""
if ((pixdim is not None) and (shape is not None)) or \
((pixdim is None) and (shape is None)):
raise ValueError('Exactly one of pixdim or '
'shape must be specified')
newShape = shape
newPixdim = pixdim
# if pixdims were specified,
# convert them into a shape
if newPixdim is not None:
npixdim = len(newPixdim)
newPixdim = np.array(newPixdim)
oldShape = np.array(self.shape[ :npixdim])
oldPixdim = np.array(self.pixdim[:npixdim])
newShape = oldShape * (oldPixdim / newPixdim)
# pad shape to full dimensions
if len(newShape) != 3:
raise ValueError('Three dimensions must be specified')
# Rescale the voxel-to-world affine
xform = affine.rescale(oldShape, newShape, origin)
xform = affine.concat(self.getAffine('voxel', 'world'), xform)
# Now that we've got our spatial
# scaling/offset matrix, pad the
# new shape/pixdims with those
# from any non-spatial dimensions
newShape = list(newShape) + list(self.shape[ 3:])
newPixdim = list(newPixdim) + list(self.pixdim[3:])
# And create the new header
# and we're away
header = self.header.copy()
header.set_data_shape(newShape)
header.set_zooms( newPixdim)
header.set_sform( xform)
header.set_qform( xform)
return Nifti(header)
class Image(Nifti):
"""Class which represents a NIFTI image. Internally, the image is
loaded/stored using a :mod:`nibabel.nifti1.Nifti1Image` or
:mod:`nibabel.nifti2.Nifti2Image`, and data access managed by a
:class:`.ImageWrapper`.
Paul McCarthy
committed
In addition to the attributes added by the :meth:`Nifti.__init__` method,
the following attributes/properties are present on an ``Image`` instance
as properties (https://docs.python.org/2/library/functions.html#property):
============== ===========================================================
``name`` The name of this ``Image`` - defaults to the image
file name, sans-suffix.
``dataSource`` The data source of this ``Image`` - the name of the
file from where it was loaded, or some other string
describing its origin.
``nibImage`` A reference to the ``nibabel`` NIFTI image object.
``saveState`` A boolean value which is ``True`` if this image is
saved to disk, ``False`` if it is in-memory, or has
been edited.
``dataRange`` The minimum/maximum values in the image. Depending upon
the value of the ``calcRange`` parameter to
:meth:`__init__`, this may be calculated when the ``Image``
is created, or may be incrementally updated as more image
data is loaded from disk.
============== ===========================================================
The ``Image`` class adds some :class:`.Notifier` topics to those which are
already provided by the :class:`Nifti` class - listeners may register to
be notified of changes to the above properties, by registering on the
following _topic_ names (see the :class:`.Notifier` class documentation):
=============== ======================================================
``'data'`` This topic is notified whenever the image data changes
Paul McCarthy
committed
(via the :meth:`__setitem__` method). The indices/
slices of the portion of data that was modified is
passed to registered listeners as the notification
value (see :meth:`.Notifier.notify`).
``'saveState'`` This topic is notified whenever the saved state of the
image changes (i.e. data or ``voxToWorldMat`` is
edited, or the image saved to disk).
``'dataRange'`` This topic is notified whenever the image data range
is changed/adjusted.
=============== ======================================================
"""
def __init__(self,
image,
name=None,
header=None,
xform=None,
loadData=True,
calcRange=True,
indexed=False,
threaded=False,

Paul McCarthy
committed
dataSource=None,

Paul McCarthy
committed
loadMeta=False,
**kwargs):
"""Create an ``Image`` object with the given image data or file name.

Paul McCarthy
committed
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
:arg image: A string containing the name of an image file to load,
or a :mod:`numpy` array, or a :mod:`nibabel` image
object.
:arg name: A name for the image.
:arg header: If not ``None``, assumed to be a
:class:`nibabel.nifti1.Nifti1Header` or
:class:`nibabel.nifti2.Nifti2Header` to be used as the
image header. Not applied to images loaded from file,
or existing :mod:`nibabel` images.
:arg xform: A :math:`4\\times 4` affine transformation matrix
which transforms voxel coordinates into real world
coordinates. If not provided, and a ``header`` is
provided, the transformation in the header is used.
If neither a ``xform`` nor a ``header`` are provided,
an identity matrix is used. If both a ``xform`` and a
``header`` are provided, the ``xform`` is used in
preference to the header transformation.
:arg loadData: If ``True`` (the default) the image data is loaded
in to memory. Otherwise, only the image header
information is read, and the image data is kept
from disk. In either case, the image data is