Skip to content

Commit a62c841

Browse files
committed
NF: spatial processing module
Add: * resampling images to images / mapped voxels * resampling images to output space * smoothing over voxel axes * FWHM / sigma conversion
1 parent 82f0dfa commit a62c841

File tree

2 files changed

+629
-0
lines changed

2 files changed

+629
-0
lines changed

nibabel/processing.py

+296
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,296 @@
1+
# emacs: -*- mode: python-mode; py-indent-offset: 4; indent-tabs-mode: nil -*-
2+
# vi: set ft=python sts=4 ts=4 sw=4 et:
3+
### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ##
4+
#
5+
# See COPYING file distributed along with the NiBabel package for the
6+
# copyright and license terms.
7+
#
8+
### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ### ##
9+
""" Image processing functions for:
10+
11+
* smoothing
12+
* resampling
13+
* converting sd to and from FWHM
14+
15+
Smoothing and resampling routines need scipy
16+
"""
17+
from __future__ import print_function, division, absolute_import
18+
19+
import numpy as np
20+
import numpy.linalg as npl
21+
22+
from .optpkg import optional_package
23+
spnd, _, _ = optional_package('scipy.ndimage')
24+
25+
from .affines import AffineError, to_matvec, from_matvec, append_diag
26+
from .spaces import vox2out_vox
27+
from .nifti1 import Nifti1Image
28+
29+
SIGMA2FWHM = np.sqrt(8 * np.log(2))
30+
31+
32+
def fwhm2sigma(fwhm):
33+
""" Convert a FWHM value to sigma in a Gaussian kernel.
34+
35+
Parameters
36+
----------
37+
fwhm : array-like
38+
FWHM value or values
39+
40+
Returns
41+
-------
42+
sigma : array or float
43+
sigma values corresponding to `fwhm` values
44+
45+
Examples
46+
--------
47+
>>> sigma = fwhm2sigma(6)
48+
>>> sigmae = fwhm2sigma([6, 7, 8])
49+
>>> sigma == sigmae[0]
50+
True
51+
"""
52+
return np.asarray(fwhm) / SIGMA2FWHM
53+
54+
55+
def sigma2fwhm(sigma):
56+
""" Convert a sigma in a Gaussian kernel to a FWHM value
57+
58+
Parameters
59+
----------
60+
sigma : array-like
61+
sigma value or values
62+
63+
Returns
64+
-------
65+
fwhm : array or float
66+
fwhm values corresponding to `sigma` values
67+
68+
Examples
69+
--------
70+
>>> fwhm = sigma2fwhm(3)
71+
>>> fwhms = sigma2fwhm([3, 4, 5])
72+
>>> fwhm == fwhms[0]
73+
True
74+
"""
75+
return np.asarray(sigma) * SIGMA2FWHM
76+
77+
78+
def adapt_affine(affine, n_dim):
79+
""" Adapt input / output dimensions of spatial `affine` for `n_dims`
80+
81+
Adapts a spatial (4, 4) affine that is being applied to an image with fewer
82+
than 3 spatial dimensions, or more than 3 dimensions. If there are more
83+
than three dimensions, assume an identity transformation for these
84+
dimensions.
85+
86+
Parameters
87+
----------
88+
affine : array-like
89+
affine transform. Usually shape (4, 4). For what follows ``N, M =
90+
affine.shape``
91+
n_dims : int
92+
Number of dimensions of underlying array, and therefore number of input
93+
dimensions for affine.
94+
95+
Returns
96+
-------
97+
adapted : shape (M, n_dims+1) array
98+
Affine array adapted to number of input dimensions. Columns of the
99+
affine corresponding to missing input dimensions have been dropped,
100+
columns corresponding to extra input dimensions have an extra identity
101+
column added
102+
"""
103+
affine = np.asarray(affine)
104+
rzs, trans = to_matvec(affine)
105+
# For missing input dimensions, drop columns in rzs
106+
rzs = rzs[:, :n_dim]
107+
adapted = from_matvec(rzs, trans)
108+
n_extra_columns = n_dim - adapted.shape[1] + 1
109+
if n_extra_columns > 0:
110+
adapted = append_diag(adapted, np.ones((n_extra_columns,)))
111+
return adapted
112+
113+
114+
def resample_from_to(from_img,
115+
to_vox_map,
116+
order=3,
117+
mode='constant',
118+
cval = 0.,
119+
out_class=Nifti1Image):
120+
""" Resample image `from_img` to mapped voxel space `to_vox_map`
121+
122+
Parameters
123+
----------
124+
from_img : object
125+
Object having attributes ``dataobj``, ``affine``, ``header``. If
126+
`out_class` is not None, ``img.__class__`` should be able to construct
127+
an image from data, affine and header.
128+
to_vox_map : image object or length 2 sequence
129+
If object, has attributes ``shape`` giving input voxel shape, and
130+
``affine`` giving mapping of input voxels to output space. If length 2
131+
sequence, elements are (shape, affine) with same meaning as above. The
132+
affine is a (4, 4) array-like.
133+
order : int, optional
134+
The order of the spline interpolation, default is 3. The order has to
135+
be in the range 0-5 (see ``scipy.ndimage.affine_transform``)
136+
mode : str, optional
137+
Points outside the boundaries of the input are filled according
138+
to the given mode ('constant', 'nearest', 'reflect' or 'wrap').
139+
Default is 'constant' (see ``scipy.ndimage.affine_transform``)
140+
cval : scalar, optional
141+
Value used for points outside the boundaries of the input if
142+
``mode='constant'``. Default is 0.0 (see
143+
``scipy.ndimage.affine_transform``)
144+
out_class : None or SpatialImage class, optional
145+
Class of output image. If None, use ``from_img.__class__``.
146+
147+
Returns
148+
-------
149+
out_img : object
150+
Image of instance specified by `out_class`, containing data output from
151+
resampling `from_img` into axes aligned to the output space of
152+
``from_img.affine``
153+
"""
154+
try:
155+
to_shape, to_affine = to_vox_map.shape, to_vox_map.affine
156+
except AttributeError:
157+
to_shape, to_affine = to_vox_map
158+
a_to_affine = adapt_affine(to_affine, len(to_shape))
159+
if out_class is None:
160+
out_class = from_img.__class__
161+
from_n_dim = len(from_img.shape)
162+
if from_n_dim < 3:
163+
raise AffineError('from_img must be at least 3D')
164+
a_from_affine = adapt_affine(from_img.affine, from_n_dim)
165+
to_vox2from_vox = npl.inv(a_from_affine).dot(a_to_affine)
166+
rzs, trans = to_matvec(to_vox2from_vox)
167+
data = spnd.affine_transform(from_img.dataobj,
168+
rzs,
169+
trans,
170+
to_shape,
171+
order = order,
172+
mode = mode,
173+
cval = cval)
174+
return out_class(data, to_affine, from_img.header)
175+
176+
177+
def resample_to_output(in_img,
178+
voxel_sizes = None,
179+
order=3,
180+
mode='constant',
181+
cval = 0.,
182+
out_class=Nifti1Image):
183+
""" Resample image `in_img` to output voxel axes (world space)
184+
185+
Parameters
186+
----------
187+
in_img : object
188+
Object having attributes ``dataobj``, ``affine``, ``header``. If
189+
`out_class` is not None, ``img.__class__`` should be able to construct
190+
an image from data, affine and header.
191+
voxel_sizes : None or sequence
192+
Gives the diagonal entries of ``out_img.affine` (except the trailing 1
193+
for the homogenous coordinates) (``out_img.affine == np.diag(voxel_sizes
194+
+ [1])``). If None, return identity `out_img.affine`.
195+
order : int, optional
196+
The order of the spline interpolation, default is 3. The order has to
197+
be in the range 0-5 (see ``scipy.ndimage.affine_transform``).
198+
mode : str, optional
199+
Points outside the boundaries of the input are filled according
200+
to the given mode ('constant', 'nearest', 'reflect' or 'wrap').
201+
Default is 'constant' (see ``scipy.ndimage.affine_transform``).
202+
cval : scalar, optional
203+
Value used for points outside the boundaries of the input if
204+
``mode='constant'``. Default is 0.0 (see
205+
``scipy.ndimage.affine_transform``).
206+
out_class : None or SpatialImage class, optional
207+
Class of output image. If None, use ``in_img.__class__``.
208+
209+
Returns
210+
-------
211+
out_img : object
212+
Image of instance specified by `out_class`, containing data output from
213+
resampling `in_img` into axes aligned to the output space of
214+
``in_img.affine``
215+
"""
216+
if out_class is None:
217+
out_class = in_img.__class__
218+
# Allow 2D images by promoting to 3D. We might want to see what a slice
219+
# looks like when resampled into world coordinates
220+
in_shape = in_img.shape
221+
n_dim = len(in_shape)
222+
if n_dim < 3: # Expand image to 3D, make voxel sizes match
223+
new_shape = in_shape + (1,) * (3 - n_dim)
224+
data = in_img.get_data().reshape(new_shape) # 2D data should be small
225+
in_img = out_class(data, in_img.affine, in_img.header)
226+
if not voxel_sizes is None and len(voxel_sizes) == n_dim:
227+
# Need to pad out voxel sizes to match new image dimensions
228+
voxel_sizes = tuple(voxel_sizes) + (1,) * (3 - n_dim)
229+
out_vox_map = vox2out_vox((in_img.shape, in_img.affine), voxel_sizes)
230+
return resample_from_to(in_img, out_vox_map, order, mode, cval, out_class)
231+
232+
233+
def smooth_image(img,
234+
fwhm,
235+
mode = 'nearest',
236+
cval = 0.,
237+
out_class=Nifti1Image):
238+
""" Smooth image `img` along voxel axes by FWHM `fwhm` millimeters
239+
240+
Parameters
241+
----------
242+
img : object
243+
Object having attributes ``dataobj``, ``affine``, ``header``. If
244+
`out_class` is not None, ``img.__class__`` should be able to construct
245+
an image from data, affine and header.
246+
fwhm : scalar or length 3 sequence
247+
FWHM *in mm* over which to smooth. The smoothing applies to the voxel
248+
axes, not to the output axes, but is in millimeters. The function
249+
adjusts the FWHM to voxels using the voxel sizes calculated from the
250+
affine. A scalar implies the same smoothing across the spatial
251+
dimensions of the image, but 0 smoothing over any further dimensions
252+
such as time. A vector should be the same length as the number of
253+
image dimensions.
254+
mode : str, optional
255+
Points outside the boundaries of the input are filled according
256+
to the given mode ('constant', 'nearest', 'reflect' or 'wrap').
257+
Default is 'nearest'. This is different from the default for
258+
``scipy.ndimage.affine_transform``, which is 'constant'. 'nearest'
259+
might be a better choice when smoothing to the edge of an image where
260+
there is still strong brain signal, otherwise this signal will get
261+
blurred towards zero.
262+
cval : scalar, optional
263+
Value used for points outside the boundaries of the input if
264+
``mode='constant'``. Default is 0.0 (see
265+
``scipy.ndimage.affine_transform``).
266+
out_class : None or SpatialImage class, optional
267+
Class of output image. If None, use ``img.__class__``.
268+
269+
Returns
270+
-------
271+
smoothed_img : object
272+
Image of instance specified by `out_class`, containing data output from
273+
smoothing `img` data by given FWHM kernel.
274+
"""
275+
if out_class is None:
276+
out_class = img.__class__
277+
n_dim = len(img.shape)
278+
# TODO: make sure time axis is last
279+
# Pad out fwhm from scalar, adding 0 for fourth etc (time etc) dimensions
280+
fwhm = np.asarray(fwhm)
281+
if fwhm.size == 1:
282+
fwhm_scalar = fwhm
283+
fwhm = np.zeros((n_dim,))
284+
fwhm[:3] = fwhm_scalar
285+
# Voxel sizes
286+
RZS = img.affine[:-1, :n_dim]
287+
vox = np.sqrt(np.sum(RZS ** 2, 0))
288+
# Smoothing in terms of voxels
289+
vox_fwhm = fwhm / vox
290+
vox_sd = fwhm2sigma(vox_fwhm)
291+
# Do the smoothing
292+
sm_data = spnd.gaussian_filter(img.dataobj,
293+
vox_sd,
294+
mode = mode,
295+
cval = cval)
296+
return out_class(sm_data, img.affine, img.header)

0 commit comments

Comments
 (0)