Skip to content

Commit ea18e1e

Browse files
authored
Merge pull request #355 from nipreps/enh/calculate-bspline-dense-reference
ENH: Function to calculate reference grids aligned with the coefficients
2 parents 1edb800 + f01cee1 commit ea18e1e

File tree

2 files changed

+154
-1
lines changed

2 files changed

+154
-1
lines changed

sdcflows/utils/tests/test_tools.py

Lines changed: 82 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -22,12 +22,93 @@
2222
#
2323
"""Test EPI manipulation routines."""
2424
import numpy as np
25+
import pytest
2526
import nibabel as nb
26-
from ..tools import brain_masker
27+
from nitransforms.linear import Affine
28+
from sdcflows.utils.tools import brain_masker, deoblique_and_zooms
2729

2830

2931
def test_epi_mask(tmpdir, testdata_dir):
3032
"""Check mask algorithm."""
3133
tmpdir.chdir()
3234
mask = brain_masker(testdata_dir / "epi.nii.gz")[-1]
3335
assert abs(np.asanyarray(nb.load(mask).dataobj).sum() - 166476) < 10
36+
37+
38+
@pytest.mark.parametrize("padding", [0, 1, 4])
39+
@pytest.mark.parametrize("factor", [1, 4, 0.8])
40+
@pytest.mark.parametrize("centered", [True, False])
41+
@pytest.mark.parametrize(
42+
"rotate",
43+
[
44+
np.eye(4),
45+
# Rotate 90 degrees around x-axis
46+
np.array([[0, 1, 0, 0], [0, 0, 1, 0], [1, 0, 0, 0], [0, 0, 0, 1]]),
47+
# Rotate 30 degrees around x-axis
48+
nb.affines.from_matvec(
49+
nb.eulerangles.euler2mat(0, 0, 30 * np.pi / 180),
50+
(0, 0, 0),
51+
),
52+
# Rotate 48 degrees around y-axis and translation
53+
nb.affines.from_matvec(
54+
nb.eulerangles.euler2mat(0, 48 * np.pi / 180, 0),
55+
(2.0, 1.2, 0.7),
56+
),
57+
],
58+
)
59+
def test_deoblique_and_zooms(tmpdir, padding, factor, centered, rotate, debug=False):
60+
"""Check deoblique and denser."""
61+
tmpdir.chdir()
62+
63+
# Generate an example reference image
64+
ref_shape = (80, 128, 160) if factor < 1 else (20, 32, 40)
65+
ref_data = np.zeros(ref_shape, dtype=np.float32)
66+
ref_data[1:-2, 10:-11, 1:-2] = 1
67+
ref_affine = np.diag([2.0, 1.25, 1.0, 1.0])
68+
ref_zooms = nb.affines.voxel_sizes(ref_affine)
69+
if centered:
70+
ref_affine[:3, 3] -= nb.affines.apply_affine(
71+
ref_affine,
72+
0.5 * (np.array(ref_data.shape) - 1),
73+
)
74+
ref_img = nb.Nifti1Image(ref_data, ref_affine)
75+
ref_img.header.set_qform(ref_affine, 1)
76+
ref_img.header.set_sform(ref_affine, 0)
77+
78+
# Generate an example oblique image
79+
mov_affine = rotate @ ref_affine
80+
mov_img = nb.Nifti1Image(ref_data, mov_affine)
81+
mov_img.header.set_qform(mov_affine, 1)
82+
mov_img.header.set_sform(mov_affine, 0)
83+
84+
# Call function with default parameters
85+
out_img = deoblique_and_zooms(ref_img, mov_img, padding=padding, factor=factor)
86+
87+
# Check output shape and zooms
88+
assert np.allclose(out_img.header.get_zooms()[:3], ref_zooms / factor)
89+
90+
# Check idempotency with a lot of tolerance
91+
ref_resampled = Affine(reference=out_img).apply(ref_img, order=0)
92+
ref_back = Affine(reference=ref_img).apply(ref_resampled, order=0)
93+
resampled = Affine(reference=out_img).apply(mov_img, order=2)
94+
if debug:
95+
ref_img.to_filename("reference.nii.gz")
96+
ref_back.to_filename("reference_back.nii.gz")
97+
ref_resampled.to_filename("reference_resampled.nii.gz")
98+
mov_img.to_filename("moving.nii.gz")
99+
resampled.to_filename("resampled.nii.gz")
100+
101+
# Allow up to 3% pixels wrong after up(down)sampling and walk back.
102+
assert (
103+
np.abs(np.clip(ref_back.get_fdata(), 0, 1) - ref_data).sum()
104+
< ref_data.size * 0.03
105+
)
106+
107+
vox_vol_out = np.prod(out_img.header.get_zooms())
108+
vox_vol_mov = np.prod(mov_img.header.get_zooms())
109+
vol_factor = vox_vol_out / vox_vol_mov
110+
111+
ref_volume = ref_data.sum()
112+
res_volume = np.clip(resampled.get_fdata(), 0, 1).sum() * vol_factor
113+
# Tolerate up to 2% variation of the volume of the moving image
114+
assert abs(1 - res_volume / ref_volume) < 0.02

sdcflows/utils/tools.py

Lines changed: 72 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -21,6 +21,78 @@
2121
# https://www.nipreps.org/community/licensing/
2222
#
2323
"""Image processing tools."""
24+
import nibabel as nb
25+
26+
27+
def deoblique_and_zooms(
28+
in_reference: nb.spatialimages.SpatialImage,
29+
oblique: nb.spatialimages.SpatialImage,
30+
factor: int = 4,
31+
padding: int = 1,
32+
factor_tol: float = 1e-4,
33+
):
34+
"""
35+
Generate a sampling reference aligned with in_reference fully covering oblique.
36+
37+
Parameters
38+
----------
39+
in_reference : :obj:`~nibabel.spatialimages.SpatialImage`
40+
The sampling reference.
41+
oblique : :obj:`~nibabel.spatialimages.SpatialImage`
42+
The rotated coordinate system whose extent should be fully covered by the
43+
generated reference.
44+
factor : :obj:`int`
45+
A factor to increase the resolution of the generated reference.
46+
padding : :obj:`int`
47+
Number of additional voxels around the most extreme positions of the projection of
48+
oblique on to the reference.
49+
factor_tol : :obj:`float`
50+
Absolute tolerance to determine whether factor is one.
51+
52+
"""
53+
from itertools import product
54+
import numpy as np
55+
from nibabel.affines import apply_affine, rescale_affine
56+
57+
# Reference space metadata
58+
hdr = in_reference.header.copy()
59+
affine = in_reference.affine.copy()
60+
ref_shape = np.array(in_reference.shape[:3])
61+
ref_zooms = np.array(hdr.get_zooms()[:3])
62+
_, scode = in_reference.get_sform(coded=True)
63+
_, qcode = in_reference.get_qform(coded=True)
64+
65+
# Calculate the 8 most extreme coordinates of oblique in in_reference space
66+
corners = np.array(list(product((0, 1), repeat=3))) * (
67+
np.array(oblique.shape[:3]) - 1
68+
)
69+
extent_ijk = apply_affine(np.linalg.inv(affine) @ oblique.affine, corners)
70+
71+
underflow = np.clip(extent_ijk.min(0) - padding, None, 0).astype(int)
72+
overflow = np.ceil(
73+
np.clip(extent_ijk.max(0) + padding + 1 - ref_shape, 0, None)
74+
).astype(int)
75+
if np.any(underflow < 0) or np.any(overflow > 0):
76+
# Add under/overflow voxels
77+
ref_shape += overflow - underflow
78+
# Consistently update origin
79+
affine[:-1, -1] = apply_affine(affine, underflow)
80+
81+
# Make grid denser
82+
if abs(1.0 - factor) > factor_tol:
83+
new_shape = np.rint(ref_shape * factor)
84+
affine = rescale_affine(affine, ref_shape, ref_zooms / factor, new_shape)
85+
ref_shape = new_shape
86+
87+
# Generate new reference
88+
hdr.set_sform(affine, scode)
89+
hdr.set_qform(affine, qcode)
90+
91+
return in_reference.__class__(
92+
nb.fileslice.strided_scalar(ref_shape.astype(int)),
93+
affine,
94+
hdr,
95+
)
2496

2597

2698
def resample_to_zooms(in_file, zooms, order=3, prefilter=True):

0 commit comments

Comments
 (0)