Skip to content

Commit bb1a257

Browse files
committed
enh: add a pure-Python interface to resample to specific resolutions
1 parent 736306b commit bb1a257

File tree

3 files changed

+133
-1
lines changed

3 files changed

+133
-1
lines changed

niworkflows/interfaces/images.py

Lines changed: 33 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -18,6 +18,39 @@
1818
LOGGER = logging.getLogger('nipype.interface')
1919

2020

21+
class _RegridToZoomsInputSpec(BaseInterfaceInputSpec):
22+
in_file = File(exists=True, mandatory=True,
23+
desc="a file whose resolution is to change")
24+
zooms = traits.Tuple(traits.Float, traits.Float, traits.Float,
25+
mandatory=True, desc="the new resolution")
26+
order = traits.Int(3, usedefault=True, desc="order of interpolator")
27+
clip = traits.Bool(True, usedefault=True,
28+
desc="clip the data array within the original image's range")
29+
30+
31+
class _RegridToZoomsOutputSpec(TraitedSpec):
32+
out_file = File(exists=True, dec="the regridded file")
33+
34+
35+
class RegridToZooms(SimpleInterface):
36+
"""Change the resolution of an image (regrid)."""
37+
38+
input_spec = _RegridToZoomsInputSpec
39+
output_spec = _RegridToZoomsOutputSpec
40+
41+
def _run_interface(self, runtime):
42+
from ..utils.images import resample_by_spacing
43+
self._results["out_file"] = fname_presuffix(
44+
self.inputs.in_file, suffix="_regrid", newpath=runtime.cwd)
45+
resample_by_spacing(
46+
self.inputs.in_file,
47+
self.inputs.zooms,
48+
order=self.inputs.order,
49+
clip=self.inputs.clip,
50+
).to_filename(self._results["out_file"])
51+
return runtime
52+
53+
2154
class _IntraModalMergeInputSpec(BaseInterfaceInputSpec):
2255
in_files = InputMultiPath(File(exists=True), mandatory=True,
2356
desc='input files')

niworkflows/utils/images.py

Lines changed: 74 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -127,3 +127,77 @@ def dseg_label(in_seg, label, newpath=None):
127127
new.set_data_dtype(np.uint8)
128128
new.to_filename(out_file)
129129
return out_file
130+
131+
132+
def resample_by_spacing(in_file, zooms, order=3, clip=True):
133+
"""Regrid the input image to match the new zooms."""
134+
from pathlib import Path
135+
import numpy as np
136+
import nibabel as nb
137+
from scipy.ndimage import map_coordinates
138+
139+
if isinstance(in_file, (str, Path)):
140+
in_file = nb.load(in_file)
141+
142+
hdr = in_file.header.copy()
143+
dtype = hdr.get_data_dtype()
144+
data = np.asanyarray(in_file.dataobj)
145+
zooms = np.array(zooms)
146+
147+
# Calculate the factors to normalize voxel size to the specific zooms
148+
pre_zooms = np.array(in_file.header.get_zooms()[:3])
149+
150+
# Calculate an affine aligned with cardinal axes, for simplicity
151+
card = nb.affines.from_matvec(np.diag(pre_zooms))
152+
extent = card[:3, :3].dot(np.array(in_file.shape[:3]))
153+
card[:3, 3] = -0.5 * extent
154+
155+
# Cover the FoV with the new grid
156+
new_size = np.ceil(extent / zooms).astype(int)
157+
offset = (extent - np.diag(zooms).dot(new_size)) * 0.5
158+
new_card = nb.affines.from_matvec(np.diag(zooms), card[:3, 3] + offset)
159+
160+
# Calculate the new indexes
161+
new_grid = np.array(np.meshgrid(
162+
np.arange(new_size[0]),
163+
np.arange(new_size[1]),
164+
np.arange(new_size[2]),
165+
indexing="ij")
166+
).reshape((3, -1))
167+
168+
# Calculate the locations of the new samples, w.r.t. the original grid
169+
ijk = np.linalg.inv(card).dot(new_card.dot(
170+
np.vstack((new_grid, np.ones((1, new_grid.shape[1]))))
171+
))
172+
173+
# Resample data in the new grid
174+
resampled = map_coordinates(
175+
data,
176+
ijk[:3, :],
177+
output=dtype,
178+
order=order,
179+
mode="constant",
180+
cval=0,
181+
prefilter=True,
182+
).reshape(new_size)
183+
if clip:
184+
resampled = np.clip(resampled, min=data.min(), max=data.max())
185+
186+
# Prepare output x-forms
187+
sform, scode = hdr.get_sform(coded=True)
188+
qform, qcode = hdr.get_qform(coded=True)
189+
190+
# Get the original image's affine
191+
affine = in_file.affine
192+
# Determine rotations w.r.t. cardinal axis and eccentricity
193+
rot = affine.dot(np.linalg.inv(card))
194+
# Apply to the new cardinal, so that the resampling is consistent
195+
new_affine = rot.dot(new_card)
196+
197+
if scode != 0:
198+
hdr.set_sform(new_affine.dot(np.linalg.inv(affine).dot(sform)), code=int(scode))
199+
if qcode != 0:
200+
hdr.set_qform(new_affine.dot(np.linalg.inv(affine).dot(qform)), code=int(qcode))
201+
202+
# Create a new x-form affine, aligned with cardinal axes, 1mm3 and centered.
203+
return nb.Nifti1Image(resampled, new_affine, hdr)

niworkflows/utils/tests/test_images.py

Lines changed: 26 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2,7 +2,7 @@
22
import numpy as np
33

44
import pytest
5-
from ..images import update_header_fields, overwrite_header, dseg_label
5+
from ..images import update_header_fields, overwrite_header, dseg_label, resample_by_spacing
66

77

88
def random_image():
@@ -86,3 +86,28 @@ def test_dseg_label(tmp_path):
8686

8787
new_im = nb.load(dseg_label(fname, label=2, newpath=tmp_path))
8888
assert np.all((data == 2).astype('int16') == np.int16(new_im.dataobj))
89+
90+
91+
def test_resample_by_spacing():
92+
"""Check the output zooms and data."""
93+
img_shape = (193, 229, 193) # Size from MNI152NLin2009cAsym, res=1
94+
img_affine = nb.affines.from_matvec(np.eye(3), (-96, -132, -78))
95+
96+
data = np.random.normal(size=img_shape).astype(float)
97+
new_affine = nb.affines.from_matvec(2.0 * np.eye(3), (-96.5, -132.5, -78.5))
98+
99+
nii = nb.Nifti1Image(data, img_affine, None)
100+
nii.set_qform(img_affine, code=4)
101+
resampled = resample_by_spacing(nii, (2.0, 2.0, 2.0), order=1, clip=False)
102+
assert resampled.header.get_zooms()[:3] == (2.0, 2.0, 2.0)
103+
assert np.all(resampled.affine == new_affine)
104+
105+
# Create a rotation matrix
106+
rot = nb.affines.from_matvec(nb.eulerangles.euler2mat(0.01, 0.02, 0.004))
107+
# Check this works with oblique images
108+
nii = nb.Nifti1Image(data, rot.dot(img_affine), None)
109+
nii.set_qform(img_affine, code=4)
110+
nii.set_sform(rot.dot(img_affine), code=4)
111+
resampled = resample_by_spacing(nii, (2.0, 2.0, 2.0), order=1, clip=False)
112+
assert resampled.header.get_zooms()[:3] == (2.0, 2.0, 2.0)
113+
assert np.allclose(resampled.affine, rot.dot(new_affine))

0 commit comments

Comments
 (0)