-
Notifications
You must be signed in to change notification settings - Fork 59
/
streamlines.py
362 lines (320 loc) · 14.9 KB
/
streamlines.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
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
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
# -*- coding: utf-8 -*-
import itertools
from dipy.io.stateful_tractogram import StatefulTractogram
from dipy.tracking.metrics import length
from dipy.tracking.streamline import set_number_of_points
from dipy.tracking.vox2track import _streamlines_in_mask
from nibabel.affines import apply_affine
import numpy as np
from scilpy.tractanalysis.streamlines_metrics import compute_tract_counts_map
def streamlines_in_mask(sft, target_mask, all_in=False):
"""
Parameters
----------
sft : StatefulTractogram
StatefulTractogram containing the streamlines to segment.
target_mask : numpy.ndarray
Binary mask in which the streamlines should pass.
Returns
-------
ids : list
Ids of the streamlines passing through the mask.
"""
sft.to_vox()
sft.to_corner()
# Copy-Paste from Dipy to get indices
if all_in:
target_mask = np.array(target_mask, dtype=np.bool, copy=True)
target_mask = np.invert(target_mask)
tractogram_mask = compute_tract_counts_map(sft.streamlines,
target_mask.shape)
tractogram_mask[tractogram_mask > 0] = 1
tmp_mask = tractogram_mask.astype(
np.uint8)*target_mask.astype(np.uint8)
streamlines_case = _streamlines_in_mask(list(sft.streamlines),
tmp_mask,
np.eye(3), [0, 0, 0])
return np.where(streamlines_case == [0, 1][False])[0].tolist()
else:
target_mask = np.array(target_mask, dtype=np.uint8, copy=True)
streamlines_case = _streamlines_in_mask(list(sft.streamlines),
target_mask,
np.eye(3), [0, 0, 0])
return np.where(streamlines_case == [0, 1][True])[0].tolist()
def filter_grid_roi(sft, mask, filter_type, is_exclude):
"""
Parameters
----------
sft : StatefulTractogram
StatefulTractogram containing the streamlines to segment.
target_mask : numpy.ndarray
Binary mask in which the streamlines should pass.
filter_type: str
One of the 3 following choices, 'any', 'all', 'either_end', 'both_ends'.
is_exclude: bool
Value to indicate if the ROI is an AND (false) or a NOT (true).
Returns
-------
ids : tuple
Filtered sft.
Ids of the streamlines passing through the mask.
"""
line_based_indices = []
if filter_type in ['any', 'all']:
line_based_indices = streamlines_in_mask(sft, mask,
all_in=filter_type == 'all')
else:
sft.to_vox()
sft.to_corner()
streamline_vox = sft.streamlines
# For endpoint filtering, we need to keep 2 separately
# Could be faster for either end, but the code look cleaner like this
line_based_indices_1 = []
line_based_indices_2 = []
for i, line_vox in enumerate(streamline_vox):
voxel_1 = tuple(line_vox[0].astype(np.int16))
voxel_2 = tuple(line_vox[-1].astype(np.int16))
if mask[voxel_1]:
line_based_indices_1.append(i)
if mask[voxel_2]:
line_based_indices_2.append(i)
# Both endpoints need to be in the mask (AND)
if filter_type == 'both_ends':
line_based_indices = np.intersect1d(line_based_indices_1,
line_based_indices_2)
# Only one endpoint need to be in the mask (OR)
elif filter_type == 'either_end':
line_based_indices = np.union1d(line_based_indices_1,
line_based_indices_2)
# If the 'exclude' option is used, the selection is inverted
if is_exclude:
line_based_indices = np.setdiff1d(range(len(sft)),
np.unique(line_based_indices))
line_based_indices = np.asarray(line_based_indices).astype(np.int32)
# From indices to sft
streamlines = sft.streamlines[line_based_indices]
data_per_streamline = sft.data_per_streamline[line_based_indices]
data_per_point = sft.data_per_point[line_based_indices]
new_sft = StatefulTractogram.from_sft(streamlines, sft,
data_per_streamline=data_per_streamline,
data_per_point=data_per_point)
return new_sft, line_based_indices
def pre_filtering_for_geometrical_shape(sft, size,
center, filter_type,
is_in_vox):
"""
Parameters
----------
sft : StatefulTractogram
StatefulTractogram containing the streamlines to segment.
size : numpy.ndarray (3)
Size in mm, x/y/z of the ROI.
center: numpy.ndarray (3)
Center x/y/z of the ROI.
filter_type: str
One of the 3 following choices, 'any', 'all', 'either_end', 'both_ends'.
is_in_vox: bool
Value to indicate if the ROI is in voxel space.
Returns
-------
ids : tuple
Filtered sft.
Ids of the streamlines passing through the mask.
"""
transfo, dim, _, _ = sft.space_attributes
inv_transfo = np.linalg.inv(transfo)
# Create relevant info about the ellipsoid in vox/world space
if is_in_vox:
center = np.asarray(apply_affine(transfo, center))
bottom_corner = center - size
top_corner = center + size
x_val = [bottom_corner[0], top_corner[0]]
y_val = [bottom_corner[1], top_corner[1]]
z_val = [bottom_corner[2], top_corner[2]]
corner_world = list(itertools.product(x_val, y_val, z_val))
corner_vox = apply_affine(inv_transfo, corner_world)
# Since the filtering using a grid is so fast, we pre-filter
# using a BB around the ellipsoid
min_corner = np.min(corner_vox, axis=0) - 1.0
max_corner = np.max(corner_vox, axis=0) + 1.5
pre_mask = np.zeros(dim)
min_x, max_x = int(max(min_corner[0], 0)), int(min(max_corner[0], dim[0]))
min_y, max_y = int(max(min_corner[1], 0)), int(min(max_corner[1], dim[1]))
min_z, max_z = int(max(min_corner[2], 0)), int(min(max_corner[2], dim[2]))
pre_mask[min_x:max_x, min_y:max_y, min_z:max_z] = 1
return filter_grid_roi(sft, pre_mask, filter_type, False)
def filter_ellipsoid(sft, ellipsoid_radius, ellipsoid_center,
filter_type, is_exclude, is_in_vox=False):
"""
Parameters
----------
sft : StatefulTractogram
StatefulTractogram containing the streamlines to segment.
ellipsoid_radius : numpy.ndarray (3)
Size in mm, x/y/z of the ellipsoid.
ellipsoid_center: numpy.ndarray (3)
Center x/y/z of the ellipsoid.
filter_type: str
One of the 3 following choices, 'any', 'all', 'either_end', 'both_ends'.
is_exclude: bool
Value to indicate if the ROI is an AND (false) or a NOT (true).
is_in_vox: bool
Value to indicate if the ROI is in voxel space.
Returns
-------
ids : tuple
Filtered sft.
Ids of the streamlines passing through the mask.
"""
pre_filtered_sft, pre_filtered_indices = \
pre_filtering_for_geometrical_shape(sft, ellipsoid_radius,
ellipsoid_center, filter_type,
is_in_vox)
pre_filtered_sft.to_rasmm()
pre_filtered_sft.to_center()
pre_filtered_streamlines = pre_filtered_sft.streamlines
transfo, _, res, _ = sft.space_attributes
if is_in_vox:
ellipsoid_center = np.asarray(apply_affine(transfo,
ellipsoid_center))
selected_by_ellipsoid = []
line_based_indices_1 = []
line_based_indices_2 = []
# This is still point based (but resampled), I had a ton of problems trying
# to use something with intersection, but even if I could do it :
# The result won't be identical to MI-Brain since I am not using the
# vtkPolydata. Also it won't be identical to TrackVis either,
# because TrackVis is point-based for Spherical ROI...
ellipsoid_radius = np.asarray(ellipsoid_radius)
ellipsoid_center = np.asarray(ellipsoid_center)
for i, line in enumerate(pre_filtered_streamlines):
if filter_type in ['any', 'all']:
# Resample to 1/10 of the voxel size
nb_points = max(int(length(line) / np.average(res) * 10), 2)
line = set_number_of_points(line, nb_points)
points_in_ellipsoid = np.sum(
((line - ellipsoid_center) / ellipsoid_radius) ** 2,
axis=1)
if filter_type == 'any' \
and np.argwhere(points_in_ellipsoid <= 1).any():
# If at least one point was in the ellipsoid
selected_by_ellipsoid.append(pre_filtered_indices[i])
elif filter_type == 'all' \
and len(np.argwhere(points_in_ellipsoid <= 1)) == len(line):
# If all points were in the ellipsoid
selected_by_ellipsoid.append(pre_filtered_indices[i])
else:
points_in_ellipsoid = np.sum(
((line[0] - ellipsoid_center) / ellipsoid_radius) ** 2)
if points_in_ellipsoid <= 1.0:
line_based_indices_1.append(pre_filtered_indices[i])
points_in_ellipsoid = np.sum(
((line[-1] - ellipsoid_center) / ellipsoid_radius) ** 2)
if points_in_ellipsoid <= 1.0:
line_based_indices_2.append(pre_filtered_indices[i])
# Both endpoints need to be in the mask (AND)
if filter_type == 'both_ends':
selected_by_ellipsoid = np.intersect1d(line_based_indices_1,
line_based_indices_2)
# Only one endpoint needs to be in the mask (OR)
elif filter_type == 'either_end':
selected_by_ellipsoid = np.union1d(line_based_indices_1,
line_based_indices_2)
# If the 'exclude' option is used, the selection is inverted
if is_exclude:
selected_by_ellipsoid = np.setdiff1d(range(len(sft)),
np.unique(selected_by_ellipsoid))
line_based_indices = np.asarray(selected_by_ellipsoid).astype(np.int32)
# From indices to sft
streamlines = sft.streamlines[line_based_indices]
data_per_streamline = sft.data_per_streamline[line_based_indices]
data_per_point = sft.data_per_point[line_based_indices]
new_sft = StatefulTractogram.from_sft(streamlines, sft,
data_per_streamline=data_per_streamline,
data_per_point=data_per_point)
return new_sft, line_based_indices
def filter_cuboid(sft, cuboid_radius, cuboid_center,
filter_type, is_exclude):
"""
Parameters
----------
sft : StatefulTractogram
StatefulTractogram containing the streamlines to segment.
cuboid_radius : numpy.ndarray (3)
Size in mm, x/y/z of the cuboid.
cuboid_center: numpy.ndarray (3)
Center x/y/z of the cuboid.
filter_type: str
One of the 3 following choices, 'any', 'all', 'either_end', 'both_ends'.
is_exclude: bool
Value to indicate if the ROI is an AND (false) or a NOT (true).
is_in_vox: bool
Value to indicate if the ROI is in voxel space.
Returns
-------
ids : tuple
Filtered sft.
Ids of the streamlines passing through the mask.
"""
pre_filtered_sft, pre_filtered_indices = \
pre_filtering_for_geometrical_shape(sft, cuboid_radius,
cuboid_center, filter_type,
False)
pre_filtered_sft.to_rasmm()
pre_filtered_sft.to_center()
pre_filtered_streamlines = pre_filtered_sft.streamlines
_, _, res, _ = sft.space_attributes
selected_by_cuboid = []
line_based_indices_1 = []
line_based_indices_2 = []
# Also here I am not using a mathematical intersection and
# I am not using vtkPolyData like in MI-Brain, so not exactly the same
cuboid_radius = np.asarray(cuboid_radius)
cuboid_center = np.asarray(cuboid_center)
for i, line in enumerate(pre_filtered_streamlines):
if filter_type in ['any', 'all']:
# Resample to 1/10 of the voxel size
nb_points = max(int(length(line) / np.average(res) * 10), 2)
line = set_number_of_points(line, nb_points)
points_in_cuboid = np.abs(line - cuboid_center) / cuboid_radius
points_in_cuboid = np.sum(np.where(points_in_cuboid <= 1, 1, 0),
axis=1)
if filter_type == 'any' \
and np.argwhere(points_in_cuboid == 3).any():
# If at least one point was in the cuboid in x/y/z
selected_by_cuboid.append(pre_filtered_indices[i])
elif filter_type == 'all' \
and len(np.argwhere(points_in_cuboid == 3)) == len(line):
# If all points were in the cuboid in x/y/z
selected_by_cuboid.append(pre_filtered_indices[i])
else:
# Faster to do it twice than trying to do in using an array of 2
points_in_cuboid = np.abs(line[0] - cuboid_center) / cuboid_radius
points_in_cuboid = np.sum(np.where(points_in_cuboid <= 1, 1, 0))
if points_in_cuboid == 3:
line_based_indices_1.append(pre_filtered_indices[i])
points_in_cuboid = np.abs(line[-1] - cuboid_center) / cuboid_radius
points_in_cuboid = np.sum(np.where(points_in_cuboid <= 1, 1, 0))
if points_in_cuboid == 3:
line_based_indices_2.append(pre_filtered_indices[i])
# Both endpoints need to be in the mask (AND)
if filter_type == 'both_ends':
selected_by_cuboid = np.intersect1d(line_based_indices_1,
line_based_indices_2)
# Only one endpoint need to be in the mask (OR)
elif filter_type == 'either_end':
selected_by_cuboid = np.union1d(line_based_indices_1,
line_based_indices_2)
# If the 'exclude' option is used, the selection is inverted
if is_exclude:
selected_by_cuboid = np.setdiff1d(range(len(sft)),
np.unique(selected_by_cuboid))
line_based_indices = np.asarray(selected_by_cuboid).astype(np.int32)
# From indices to sft
streamlines = sft.streamlines[line_based_indices]
data_per_streamline = sft.data_per_streamline[line_based_indices]
data_per_point = sft.data_per_point[line_based_indices]
new_sft = StatefulTractogram.from_sft(streamlines, sft,
data_per_streamline=data_per_streamline,
data_per_point=data_per_point)
return new_sft, line_based_indices