-
Notifications
You must be signed in to change notification settings - Fork 35
/
track.py
485 lines (442 loc) · 17.7 KB
/
track.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
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
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
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
#!/usr/bin/env python
"""
ndmg.track
~~~~~~~~~~
Contains ndmg's fiber reconstruction and tractography functionality.
Theory described here: https://neurodata.io/talks/ndmg.pdf#page=21
"""
# system imports
import os
# external package imports
import numpy as np
import nibabel as nib
# dipy imports
from dipy.tracking.streamline import Streamlines
from dipy.tracking import utils
from dipy.tracking.local_tracking import LocalTracking
from dipy.tracking.local_tracking import ParticleFilteringTracking
from dipy.tracking.stopping_criterion import BinaryStoppingCriterion
from dipy.tracking.stopping_criterion import ActStoppingCriterion
from dipy.tracking.stopping_criterion import CmcStoppingCriterion
from dipy.reconst.dti import fractional_anisotropy, TensorModel, quantize_evecs
from dipy.reconst.shm import CsaOdfModel
from dipy.reconst.csdeconv import ConstrainedSphericalDeconvModel, recursive_response
from dipy.data import get_sphere
from dipy.direction import peaks_from_model, ProbabilisticDirectionGetter
from ndmg.utils.gen_utils import timer
from ndmg.stats import qa_tensor
def build_seed_list(mask_img_file, stream_affine, dens):
"""uses dipy tractography utilities in order to create a seed list for tractography
Parameters
----------
mask_img_file : str
path to mask of area to generate seeds for
stream_affine : ndarray
4x4 array with 1s diagonally and 0s everywhere else
dens : int
seed density
Returns
-------
ndarray
locations for the seeds
"""
mask_img = nib.load(mask_img_file)
mask_img_data = mask_img.get_data().astype("bool")
seeds = utils.random_seeds_from_mask(
mask_img_data,
affine=stream_affine,
seeds_count=int(dens),
seed_count_per_voxel=True,
)
return seeds
def tens_mod_fa_est(gtab, dwi_file, B0_mask):
"""Estimate a tensor FA image to use for registrations using dipy functions
Parameters
----------
gtab : GradientTable
gradient table created from bval and bvec file
dwi_file : str
Path to eddy-corrected and RAS reoriented dwi image
B0_mask : str
Path to nodif B0 mask (averaged b0 mask)
Returns
-------
str
Path to tensor_fa image file
"""
data = nib.load(dwi_file).get_fdata()
print("Generating simple tensor FA image to use for registrations...")
nodif_B0_img = nib.load(B0_mask)
B0_mask_data = nodif_B0_img.get_fdata().astype("bool")
nodif_B0_affine = nodif_B0_img.affine
model = TensorModel(gtab)
mod = model.fit(data, B0_mask_data)
FA = fractional_anisotropy(mod.evals)
FA[np.isnan(FA)] = 0
fa_img = nib.Nifti1Image(FA.astype(np.float32), nodif_B0_affine)
fa_path = f"{os.path.dirname(B0_mask)}/tensor_fa.nii.gz"
nib.save(fa_img, fa_path)
return fa_path
class RunTrack:
def __init__(
self,
dwi_in,
nodif_B0_mask,
gm_in_dwi,
vent_csf_in_dwi,
csf_in_dwi,
wm_in_dwi,
gtab,
mod_type,
track_type,
mod_func,
qa_tensor_out,
seeds,
stream_affine,
):
"""A class for deterministic tractography in native space
Parameters
----------
dwi_in : str
path to the input dwi image to perform tractography on.
Should be a nifti, gzipped nifti, or other image that nibabel
is capable of reading, with data as a 4D object.
nodif_B0_mask : str
path to the mask of the b0 mean volume. Should be a nifti,
gzipped nifti, or other image file that nibabel is capable of
reading, with data as a 3D object.
gm_in_dwi : str
Path to gray matter segmentation in EPI space. Should be a nifti,
gzipped nifti, or other image file that nibabel is capable of
reading, with data as a 3D object
vent_csf_in_dwi : str
Ventricular CSF Mask in EPI space. Should be a nifti,
gzipped nifti, or other image file that nibabel is capable of
reading, with data as a 3D object
csf_in_dwi : str
Path to CSF mask in EPI space. Should be a nifti, gzipped nifti, or other image file that nibabel compatable
wm_in_dwi : str
Path to white matter probabilities in EPI space. Should be a nifti,
gzipped nifti, or other image file that nibabel is capable of
reading, with data as a 3D object.
gtab : gradient table
gradient table created from bval and bvec files
mod_type : str
Determinstic (det) or probabilistic (prob) tracking
track_type : str
Tracking approach: local or particle
mod_func : str
Diffusion model: csd or csa
qa_tensor: str
path to store the qa for tensor/directions of model
seeds : ndarray
ndarray of seeds for tractography
stream_affine : ndarray
4x4 2D array with 1s diagonaly and 0s everywhere else
"""
self.dwi = dwi_in
self.nodif_B0_mask = nodif_B0_mask
self.gm_in_dwi = gm_in_dwi
self.vent_csf_in_dwi = vent_csf_in_dwi
self.csf_in_dwi = csf_in_dwi
self.wm_in_dwi = wm_in_dwi
self.gtab = gtab
self.mod_type = mod_type
self.track_type = track_type
self.qa_tensor_out = qa_tensor_out
self.seeds = seeds
self.mod_func = mod_func
self.stream_affine = stream_affine
@timer
def run(self):
"""Creates the tracktography tracks using dipy commands and the specified tracking type and approach
Returns
-------
ArraySequence
contains the tractography track raw data for further analysis
Raises
------
ValueError
Raised when no seeds are supplied or no valid seeds were found in white-matter interface
ValueError
Raised when no seeds are supplied or no valid seeds were found in white-matter interface
"""
self.tiss_classifier = self.prep_tracking()
if self.mod_type == "det":
if self.mod_func == "csa":
self.mod = self.odf_mod_est()
elif self.mod_func == "csd":
self.mod = self.csd_mod_est()
if self.track_type == "local":
tracks = self.local_tracking()
elif self.track_type == "particle":
tracks = self.particle_tracking()
else:
raise ValueError(
"Error: Either no seeds supplied, or no valid seeds found in white-matter interface"
)
elif self.mod_type == "prob":
if self.mod_func == "csa":
self.mod = self.odf_mod_est()
elif self.mod_func == "csd":
self.mod = self.csd_mod_est()
if self.track_type == "local":
tracks = self.local_tracking()
elif self.track_type == "particle":
tracks = self.particle_tracking()
else:
raise ValueError(
"Error: Either no seeds supplied, or no valid seeds found in white-matter interface"
)
tracks = Streamlines([track for track in tracks if len(track) > 60])
return tracks
@staticmethod
def make_hdr(streamlines, hdr):
trk_hdr = nib.streamlines.trk.TrkFile.create_empty_header()
trk_hdr["hdr_size"] = 1000
trk_hdr["dimensions"] = hdr["dim"][1:4].astype("float32")
trk_hdr["voxel_sizes"] = hdr["pixdim"][1:4]
trk_hdr["voxel_to_rasmm"] = np.eye(4)
trk_hdr["voxel_order"] = "RAS"
trk_hdr["pad2"] = "RAS"
trk_hdr["image_orientation_patient"] = np.array(
[0.0, 0.0, 0.0, 0.0, 0.0, 0.0]
).astype("float32")
trk_hdr["endianness"] = "<"
trk_hdr["_offset_data"] = 1000
trk_hdr["nb_streamlines"] = streamlines.total_nb_rows
return trk_hdr
def prep_tracking(self):
"""Uses nibabel and dipy functions in order to load the grey matter, white matter, and csf masks
and use a tissue classifier (act, cmc, or binary) on the include/exclude maps to make a tissueclassifier object
Returns
-------
ActStoppingCriterion, CmcStoppingCriterion, or BinaryStoppingCriterion
The resulting tissue classifier object, depending on which method you use (currently only does act)
"""
if self.track_type == "local":
tiss_class = "bin"
elif self.track_type == "particle":
tiss_class = "cmc"
self.dwi_img = nib.load(self.dwi)
self.data = self.dwi_img.get_data()
# Loads mask and ensures it's a true binary mask
self.mask_img = nib.load(self.nodif_B0_mask)
self.mask = self.mask_img.get_data() > 0
# Load tissue maps and prepare tissue classifier
self.gm_mask = nib.load(self.gm_in_dwi)
self.gm_mask_data = self.gm_mask.get_data()
self.wm_mask = nib.load(self.wm_in_dwi)
self.wm_mask_data = self.wm_mask.get_data()
self.wm_in_dwi_data = nib.load(self.wm_in_dwi).get_data().astype("bool")
if tiss_class == "act":
self.vent_csf_in_dwi = nib.load(self.vent_csf_in_dwi)
self.vent_csf_in_dwi_data = self.vent_csf_in_dwi.get_data()
self.background = np.ones(self.gm_mask.shape)
self.background[
(self.gm_mask_data + self.wm_mask_data + self.vent_csf_in_dwi_data) > 0
] = 0
self.include_map = self.wm_mask_data
self.include_map[self.background > 0] = 0
self.exclude_map = self.vent_csf_in_dwi_data
self.tiss_classifier = ActStoppingCriterion(
self.include_map, self.exclude_map
)
elif tiss_class == "bin":
self.tiss_classifier = BinaryStoppingCriterion(self.wm_in_dwi_data)
# self.tiss_classifier = BinaryStoppingCriterion(self.mask)
elif tiss_class == "cmc":
self.vent_csf_in_dwi = nib.load(self.vent_csf_in_dwi)
self.vent_csf_in_dwi_data = self.vent_csf_in_dwi.get_data()
voxel_size = np.average(self.wm_mask.get_header()["pixdim"][1:4])
step_size = 0.2
self.tiss_classifier = CmcStoppingCriterion.from_pve(
self.wm_mask_data,
self.gm_mask_data,
self.vent_csf_in_dwi_data,
step_size=step_size,
average_voxel_size=voxel_size,
)
else:
pass
return self.tiss_classifier
@timer
def tens_mod_est(self):
print("Fitting tensor model...")
self.model = TensorModel(self.gtab)
self.ten = self.model.fit(self.data, self.wm_in_dwi_data)
self.fa = self.ten.fa
self.fa[np.isnan(self.fa)] = 0
self.sphere = get_sphere("repulsion724")
self.ind = quantize_evecs(self.ten.evecs, self.sphere.vertices)
return self.ten
@timer
def odf_mod_est(self):
print("Fitting CSA ODF model...")
self.mod = CsaOdfModel(self.gtab, sh_order=6)
return self.mod
@timer
def csd_mod_est(self):
print("Fitting CSD model...")
try:
print("Attempting to use spherical harmonic basis first...")
self.mod = ConstrainedSphericalDeconvModel(self.gtab, None, sh_order=6)
except:
print("Falling back to estimating recursive response...")
self.response = recursive_response(
self.gtab,
self.data,
mask=self.wm_in_dwi_data,
sh_order=6,
peak_thr=0.01,
init_fa=0.08,
init_trace=0.0021,
iter=8,
convergence=0.001,
parallel=False,
)
print("CSD Reponse: " + str(self.response))
self.mod = ConstrainedSphericalDeconvModel(self.gtab, self.response,sh_order=6)
return self.mod
@timer
def local_tracking(self):
self.sphere = get_sphere("repulsion724")
if self.mod_type == "det":
print("Obtaining peaks from model...")
self.mod_peaks = peaks_from_model(
self.mod,
self.data,
self.sphere,
relative_peak_threshold=0.5,
min_separation_angle=25,
mask=self.wm_in_dwi_data,
npeaks=5,
normalize_peaks=True,
)
qa_tensor.create_qa_figure(self.mod_peaks.peak_dirs, self.mod_peaks.peak_values, self.qa_tensor_out, self.mod_func)
self.streamline_generator = LocalTracking(
self.mod_peaks,
self.tiss_classifier,
self.seeds,
self.stream_affine,
step_size=0.5,
return_all=True,
)
elif self.mod_type == "prob":
print("Preparing probabilistic tracking...")
print("Fitting model to data...")
self.mod_fit = self.mod.fit(self.data, self.wm_in_dwi_data)
print("Building direction-getter...")
self.mod_peaks = peaks_from_model(
self.mod,
self.data,
self.sphere,
relative_peak_threshold=0.5,
min_separation_angle=25,
mask=self.wm_in_dwi_data,
npeaks=5,
normalize_peaks=True,
)
qa_tensor.create_qa_figure(self.mod_peaks.peak_dirs, self.mod_peaks.peak_values, self.qa_tensor_out, self.mod_func)
try:
print(
"Proceeding using spherical harmonic coefficient from model estimation..."
)
self.pdg = ProbabilisticDirectionGetter.from_shcoeff(
self.mod_fit.shm_coeff, max_angle=60.0, sphere=self.sphere
)
except:
print("Proceeding using FOD PMF from model estimation...")
self.fod = self.mod_fit.odf(self.sphere)
self.pmf = self.fod.clip(min=0)
self.pdg = ProbabilisticDirectionGetter.from_pmf(
self.pmf, max_angle=60.0, sphere=self.sphere
)
self.streamline_generator = LocalTracking(
self.pdg,
self.tiss_classifier,
self.seeds,
self.stream_affine,
step_size=0.5,
return_all=True,
)
print("Reconstructing tractogram streamlines...")
self.streamlines = Streamlines(self.streamline_generator)
return self.streamlines
@timer
def particle_tracking(self):
self.sphere = get_sphere("repulsion724")
if self.mod_type == "det":
maxcrossing = 1
print("Obtaining peaks from model...")
self.mod_peaks = peaks_from_model(
self.mod,
self.data,
self.sphere,
relative_peak_threshold=0.5,
min_separation_angle=25,
mask=self.wm_in_dwi_data,
npeaks=5,
normalize_peaks=True,
)
qa_tensor.create_qa_figure(self.mod_peaks.peak_dirs, self.mod_peaks.peak_values, self.qa_tensor_out, self.mod_func)
self.streamline_generator = ParticleFilteringTracking(
self.mod_peaks,
self.tiss_classifier,
self.seeds,
self.stream_affine,
max_cross=maxcrossing,
step_size=0.5,
maxlen=1000,
pft_back_tracking_dist=2,
pft_front_tracking_dist=1,
particle_count=15,
return_all=True,
)
elif self.mod_type == "prob":
maxcrossing = 2
print("Preparing probabilistic tracking...")
print("Fitting model to data...")
self.mod_fit = self.mod.fit(self.data, self.wm_in_dwi_data)
print("Building direction-getter...")
self.mod_peaks = peaks_from_model(
self.mod,
self.data,
self.sphere,
relative_peak_threshold=0.5,
min_separation_angle=25,
mask=self.wm_in_dwi_data,
npeaks=5,
normalize_peaks=True,
)
qa_tensor.create_qa_figure(self.mod_peaks.peak_dirs, self.mod_peaks.peak_values, self.qa_tensor_out, self.mod_func)
try:
print(
"Proceeding using spherical harmonic coefficient from model estimation..."
)
self.pdg = ProbabilisticDirectionGetter.from_shcoeff(
self.mod_fit.shm_coeff, max_angle=60.0, sphere=self.sphere
)
except:
print("Proceeding using FOD PMF from model estimation...")
self.fod = self.mod_fit.odf(self.sphere)
self.pmf = self.fod.clip(min=0)
self.pdg = ProbabilisticDirectionGetter.from_pmf(
self.pmf, max_angle=60.0, sphere=self.sphere
)
self.streamline_generator = ParticleFilteringTracking(
self.pdg,
self.tiss_classifier,
self.seeds,
self.stream_affine,
max_cross=maxcrossing,
step_size=0.5,
maxlen=1000,
pft_back_tracking_dist=2,
pft_front_tracking_dist=1,
particle_count=15,
return_all=True,
)
print("Reconstructing tractogram streamlines...")
self.streamlines = Streamlines(self.streamline_generator)
return self.streamlines