-
Notifications
You must be signed in to change notification settings - Fork 437
/
mapmri.py
2193 lines (1831 loc) · 72.2 KB
/
mapmri.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
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
import numpy as np
from scipy.special import gamma, genlaguerre, hermite
from dipy.reconst.base import ReconstFit, ReconstModel
from dipy.reconst.cache import Cache
from dipy.reconst.multi_voxel import multi_voxel_fit
try: # preferred scipy >= 0.14, required scipy >= 1.0
from scipy.special import factorial as sfactorial, factorial2
except ImportError:
from scipy.misc import factorial as sfactorial, factorial2
from math import factorial as mfactorial
from warnings import warn
from dipy.core.geometry import cart2sphere
from dipy.core.gradients import gradient_table
from dipy.core.optimize import Optimizer, PositiveDefiniteLeastSquares
from dipy.data import load_sdp_constraints
import dipy.reconst.dti as dti
from dipy.reconst.shm import real_sh_descoteaux_from_index, sph_harm_ind_list
from dipy.testing.decorators import warning_for_keywords
from dipy.utils.optpkg import optional_package
cvxpy, have_cvxpy, _ = optional_package("cvxpy", min_version="1.4.1")
class MapmriModel(ReconstModel, Cache):
r"""Mean Apparent Propagator MRI (MAPMRI) of the diffusion signal.
The main idea in MAPMRI footcite:p:`Ozarslan2013` is to model the diffusion
signal as a linear combination of the continuous functions presented in
footcite:p:`Ozarslan2008` but extending it in three dimensions.
The main difference with the SHORE proposed in footcite:p:`Merlet2013` is
that MAPMRI 3D extension is provided using a set of three basis functions
for the radial part, one for the signal along x, one for y and one for z,
while footcite:p:`Merlet2013` uses one basis function to model the radial
part and real Spherical Harmonics to model the angular part.
From the MAPMRI coefficients is possible to use the analytical formulae
to estimate the ODF.
See :footcite:p:`Avram2015` for additional tissue microstructure insights
provided by MAPMRI.
See also footcite:p:`Fick2016b`, footcite:p:`Cheng2012`,
footcite:p:`Hosseinbor2013`, footcite:p:`Craven1979`, and
footcite:p:`DelaHaije2020` for additional insight into to the model.
References
----------
.. footbibliography::
"""
@warning_for_keywords()
def __init__(
self,
gtab,
*,
radial_order=6,
laplacian_regularization=True,
laplacian_weighting=0.2,
positivity_constraint=False,
global_constraints=False,
pos_grid=15,
pos_radius="adaptive",
anisotropic_scaling=True,
eigenvalue_threshold=1e-04,
bval_threshold=np.inf,
dti_scale_estimation=True,
static_diffusivity=0.7e-3,
cvxpy_solver=None,
):
r"""Analytical and continuous modeling of the diffusion signal with
respect to the MAPMRI basis.
The main idea of the MAPMRI :footcite:p:`Ozarslan2013` is to model the
diffusion signal as a linear combination of the continuous functions
presented in :footcite:p:`Ozarslan2008` but extending it in three
dimensions.
The main difference with the SHORE proposed in
:footcite:p:`Ozarslan2009` is that MAPMRI 3D extension is provided using
a set of three basis functions for the radial part, one for the signal
along x, one for y and one for z, while :footcite:p:`Ozarslan2009` uses
one basis function to model the radial part and real Spherical Harmonics
to model the angular part.
From the MAPMRI coefficients it is possible to estimate various
q-space indices, the PDF and the ODF.
The fitting procedure can be constrained using the positivity
constraint proposed in :footcite:p:`Ozarslan2013` or
:footcite:p:`DelaHaije2020` and/or the laplacian regularization proposed
in :footcite:p:`Fick2016b`.
For the estimation of q-space indices we recommend using the 'regular'
anisotropic implementation of MAPMRI. However, it has been shown that
the ODF estimation in this implementation has a bias which
'squeezes together' the ODF peaks when there is a crossing at an angle
smaller than 90 degrees :footcite:p:`Fick2016b`. When you want to
estimate ODFs for tractography we therefore recommend using the
isotropic implementation (which is equivalent to
:footcite:p:`Ozarslan2009`).
The switch between isotropic and anisotropic can be easily made through
the anisotropic_scaling option.
Parameters
----------
gtab : GradientTable,
gradient directions and bvalues container class.
the gradient table has to include b0-images.
radial_order : unsigned int,
an even integer that represent the order of the basis
laplacian_regularization: bool,
Regularize using the Laplacian of the MAP-MRI basis.
laplacian_weighting: string or scalar,
The string 'GCV' makes it use generalized cross-validation
:footcite:p:`Craven1979` to find the regularization weight
:footcite:p:`DelaHaije2020`. A scalar sets the regularization
weight to that value and an array will make it selected the optimal
weight from the values in the array.
positivity_constraint : bool,
Constrain the propagator to be positive.
global_constraints : bool, optional
If set to False, positivity is enforced on a grid determined by
pos_grid and pos_radius. If set to True, positivity is enforced
everywhere using the constraints of :footcite:p:`Merlet2013`. Global
constraints are currently supported for anisotropic_scaling=True and
for radial_order <= 10.
pos_grid : int, optional
The number of points in the grid that is used in the local
positivity constraint.
pos_radius : float or string, optional
If set to a float, the maximum distance the local positivity
constraint constrains to posivity is that value. If set to
'adaptive', the maximum distance is dependent on the estimated
tissue diffusivity. If 'infinity', semidefinite programming
constraints are used :footcite:p:`DelaHaije2020`.
anisotropic_scaling : bool, optional
If True, uses the standard anisotropic MAP-MRI basis. If False,
uses the isotropic MAP-MRI basis (equal to 3D-SHORE).
eigenvalue_threshold : float, optional
Sets the minimum of the tensor eigenvalues in order to avoid
stability problem.
bval_threshold : float, optional
Sets the b-value threshold to be used in the scale factor
estimation. In order for the estimated non-Gaussianity to have
meaning this value should set to a lower value (b<2000 s/mm^2)
such that the scale factors are estimated on signal points that
reasonably represent the spins at Gaussian diffusion.
dti_scale_estimation : bool, optional
Whether or not DTI fitting is used to estimate the isotropic scale
factor for isotropic MAP-MRI.
When set to False the algorithm presets the isotropic tissue
diffusivity to static_diffusivity. This vastly increases fitting
speed but at the cost of slightly reduced fitting quality. Can
still be used in combination with regularization and constraints.
static_diffusivity : float, optional
the tissue diffusivity that is used when dti_scale_estimation is
set to False. The default is that of typical white matter
D=0.7e-3 :footcite:p:`Fick2016b`.
cvxpy_solver : str, optional
cvxpy solver name. Optionally optimize the positivity constraint
with a particular cvxpy solver. See https://www.cvxpy.org/ for
details.
Default: None (cvxpy chooses its own solver)
References
----------
.. footbibliography::
Examples
--------
In this example, where the data, gradient table and sphere tessellation
used for reconstruction are provided, we model the diffusion signal
with respect to the SHORE basis and compute the real and analytical
ODF.
>>> from dipy.data import dsi_voxels, default_sphere
>>> from dipy.core.gradients import gradient_table
>>> _, gtab_ = dsi_voxels()
>>> gtab = gradient_table(gtab_.bvals, bvecs=gtab_.bvecs,
... b0_threshold=gtab_.bvals.min())
>>> from dipy.sims.voxel import sticks_and_ball
>>> data, golden_directions = sticks_and_ball(gtab, d=0.0015, S0=1,
... angles=[(0, 0),
... (90, 0)],
... fractions=[50, 50],
... snr=None)
>>> from dipy.reconst.mapmri import MapmriModel
>>> radial_order = 4
>>> map_model = MapmriModel(gtab, radial_order=radial_order)
>>> mapfit = map_model.fit(data)
>>> odf = mapfit.odf(default_sphere)
"""
if np.sum(gtab.b0s_mask) == 0:
raise ValueError(
"gtab does not have any b0s, check in the "
"gradient_table if b0_threshold needs to be "
"increased."
)
self.gtab = gtab
if radial_order < 0 or radial_order % 2:
raise ValueError("radial_order must be a positive, even number.")
self.radial_order = radial_order
self.bval_threshold = bval_threshold
self.dti_scale_estimation = dti_scale_estimation
if laplacian_regularization:
msg = (
"Laplacian Regularization weighting must be 'GCV',"
" a positive float or an array of positive floats."
)
if isinstance(laplacian_weighting, str):
if not laplacian_weighting == "GCV":
raise ValueError(msg)
elif isinstance(laplacian_weighting, (float, np.ndarray)):
if np.sum(laplacian_weighting < 0) > 0:
raise ValueError(msg)
self.laplacian_weighting = laplacian_weighting
self.laplacian_regularization = laplacian_regularization
if positivity_constraint:
if not have_cvxpy:
raise ImportError("CVXPY package needed to enforce constraints.")
if cvxpy_solver is not None:
if cvxpy_solver not in cvxpy.installed_solvers():
installed_solvers = ", ".join(cvxpy.installed_solvers())
raise ValueError(
f"Input `cvxpy_solver` was set to"
f" {cvxpy_solver}. One of"
f" {installed_solvers} was expected."
)
self.cvxpy_solver = cvxpy_solver
if global_constraints:
if not anisotropic_scaling:
raise ValueError(
"Global constraints only available for"
" anistropic_scaling=True."
)
if radial_order > 10:
self.sdp_constraints = load_sdp_constraints("hermite", order=10)
warn(
"Global constraints are currently supported for"
" radial_order <= 10.",
stacklevel=2,
)
else:
self.sdp_constraints = load_sdp_constraints(
"hermite", order=radial_order
)
m = (2 + radial_order) * (4 + radial_order) * (3 + 2 * radial_order)
m = m // 24
self.sdp = PositiveDefiniteLeastSquares(m, A=self.sdp_constraints)
else:
msg = "pos_radius must be 'adaptive' or a positive float."
if isinstance(pos_radius, str):
if pos_radius != "adaptive":
raise ValueError(msg)
elif isinstance(pos_radius, (float, int)):
if pos_radius <= 0:
raise ValueError(msg)
self.constraint_grid = create_rspace(pos_grid, pos_radius)
if not anisotropic_scaling:
self.pos_K_independent = mapmri_isotropic_K_mu_independent(
radial_order, self.constraint_grid
)
else:
raise ValueError(msg)
self.pos_grid = pos_grid
self.pos_radius = pos_radius
self.global_constraints = global_constraints
self.positivity_constraint = positivity_constraint
self.anisotropic_scaling = anisotropic_scaling
if (gtab.big_delta is None) or (gtab.small_delta is None):
self.tau = 1 / (4 * np.pi**2)
else:
self.tau = gtab.big_delta - gtab.small_delta / 3.0
self.eigenvalue_threshold = eigenvalue_threshold
self.cutoff = gtab.bvals < self.bval_threshold
gtab_cutoff = gradient_table(
bvals=self.gtab.bvals[self.cutoff], bvecs=self.gtab.bvecs[self.cutoff]
)
self.tenmodel = dti.TensorModel(gtab_cutoff)
if self.anisotropic_scaling:
self.ind_mat = mapmri_index_matrix(self.radial_order)
self.Bm = b_mat(self.ind_mat)
self.S_mat, self.T_mat, self.U_mat = mapmri_STU_reg_matrices(radial_order)
else:
self.ind_mat = mapmri_isotropic_index_matrix(self.radial_order)
self.Bm = b_mat_isotropic(self.ind_mat)
self.laplacian_matrix = mapmri_isotropic_laplacian_reg_matrix(
radial_order, 1.0
)
qvals = np.sqrt(self.gtab.bvals / self.tau) / (2 * np.pi)
q = gtab.bvecs * qvals[:, None]
if self.dti_scale_estimation:
self.M_mu_independent = mapmri_isotropic_M_mu_independent(
self.radial_order, q
)
else:
D = static_diffusivity
mumean = np.sqrt(2 * D * self.tau)
self.mu = np.array([mumean, mumean, mumean])
self.M = mapmri_isotropic_phi_matrix(radial_order, mumean, q)
if (
self.laplacian_regularization
and isinstance(laplacian_weighting, float)
and not positivity_constraint
):
MMt = (
np.dot(self.M.T, self.M)
+ laplacian_weighting * mumean * self.laplacian_matrix
)
self.MMt_inv_Mt = np.dot(np.linalg.pinv(MMt), self.M.T)
@multi_voxel_fit
def fit(self, data, **kwargs):
errorcode = 0
tenfit = self.tenmodel.fit(data[self.cutoff])
evals = tenfit.evals
R = tenfit.evecs
evals = np.clip(evals, self.eigenvalue_threshold, evals.max())
qvals = np.sqrt(self.gtab.bvals / self.tau) / (2 * np.pi)
mu_max = max(np.sqrt(evals * 2 * self.tau)) # used for constraint
if self.anisotropic_scaling:
mu = np.sqrt(evals * 2 * self.tau)
qvecs = np.dot(self.gtab.bvecs, R)
q = qvecs * qvals[:, None]
M = mapmri_phi_matrix(self.radial_order, mu, q)
else:
try:
# self.MMt_inv_Mt
lopt = self.laplacian_weighting
coef = np.dot(self.MMt_inv_Mt, data)
coef = coef / sum(coef * self.Bm)
return MapmriFit(self, coef, self.mu, R, lopt, errorcode=errorcode)
except AttributeError:
try:
M = self.M
mu = self.mu
except AttributeError:
u0 = isotropic_scale_factor(evals * 2 * self.tau)
mu = np.array([u0, u0, u0])
M_mu_dependent = mapmri_isotropic_M_mu_dependent(
self.radial_order, mu[0], qvals
)
M = M_mu_dependent * self.M_mu_independent
if self.laplacian_regularization:
if self.anisotropic_scaling:
laplacian_matrix = mapmri_laplacian_reg_matrix(
self.ind_mat, mu, self.S_mat, self.T_mat, self.U_mat
)
else:
laplacian_matrix = self.laplacian_matrix * mu[0]
if (
isinstance(self.laplacian_weighting, str)
and self.laplacian_weighting.upper() == "GCV"
):
try:
lopt = generalized_crossvalidation(data, M, laplacian_matrix)
except np.linalg.linalg.LinAlgError:
# 1/0.
lopt = 0.05
errorcode = 1
elif np.isscalar(self.laplacian_weighting):
lopt = self.laplacian_weighting
else:
lopt = generalized_crossvalidation_array(
data, M, laplacian_matrix, weights_array=self.laplacian_weighting
)
else:
lopt = 0.0
laplacian_matrix = np.ones((self.ind_mat.shape[0], self.ind_mat.shape[0]))
if self.positivity_constraint:
data_norm = np.asarray(data / data[self.gtab.b0s_mask].mean())
if self.global_constraints:
coef = self.sdp.solve(M, data_norm, solver=self.cvxpy_solver)
else:
c = cvxpy.Variable(M.shape[1])
design_matrix = cvxpy.Constant(M) @ c
# workaround for the bug on cvxpy 1.0.15 when lopt = 0
# See https://github.com/cvxgrp/cvxpy/issues/672
if not lopt:
objective = cvxpy.Minimize(
cvxpy.sum_squares(design_matrix - data_norm)
)
else:
objective = cvxpy.Minimize(
cvxpy.sum_squares(design_matrix - data_norm)
+ lopt * cvxpy.quad_form(c, laplacian_matrix)
)
if self.pos_radius == "adaptive":
# custom constraint grid based on scale factor [Avram2015]
constraint_grid = create_rspace(self.pos_grid, np.sqrt(5) * mu_max)
else:
constraint_grid = self.constraint_grid
if self.anisotropic_scaling:
K = mapmri_psi_matrix(self.radial_order, mu, constraint_grid)
else:
if self.pos_radius == "adaptive":
# grid changes per voxel. Recompute entire K matrix.
K = mapmri_isotropic_psi_matrix(
self.radial_order, mu[0], constraint_grid
)
else:
# grid is static. Only compute mu-dependent part of K.
K_dependent = mapmri_isotropic_K_mu_dependent(
self.radial_order, mu[0], constraint_grid
)
K = K_dependent * self.pos_K_independent
M0 = M[self.gtab.b0s_mask, :]
constraints = [(M0[0] @ c) == 1, (K @ c) >= -0.1]
prob = cvxpy.Problem(objective, constraints)
try:
prob.solve(solver=self.cvxpy_solver)
coef = np.asarray(c.value).squeeze()
except Exception:
errorcode = 2
warn("Optimization did not find a solution", stacklevel=2)
try:
coef = np.dot(np.linalg.pinv(M), data) # least squares
except np.linalg.linalg.LinAlgError:
errorcode = 3
coef = np.zeros(M.shape[1])
return MapmriFit(self, coef, mu, R, lopt, errorcode=errorcode)
else:
try:
pseudoInv = np.dot(
np.linalg.inv(np.dot(M.T, M) + lopt * laplacian_matrix), M.T
)
coef = np.dot(pseudoInv, data)
except np.linalg.linalg.LinAlgError:
errorcode = 1
coef = np.zeros(M.shape[1])
return MapmriFit(self, coef, mu, R, lopt, errorcode=errorcode)
coef = coef / sum(coef * self.Bm)
return MapmriFit(self, coef, mu, R, lopt, errorcode=errorcode)
class MapmriFit(ReconstFit):
@warning_for_keywords()
def __init__(self, model, mapmri_coef, mu, R, lopt, *, errorcode=0):
"""Calculates diffusion properties for a single voxel
Parameters
----------
model : object,
AnalyticalModel
mapmri_coef : 1d ndarray,
mapmri coefficients
mu : array, shape (3,)
scale parameters vector for x, y and z
R : array, shape (3,3)
rotation matrix
lopt : float,
regularization weight used for laplacian regularization
errorcode : int
provides information on whether errors occurred in the fitting
of each voxel. 0 means no problem, 1 means a LinAlgError
occurred when trying to invert the design matrix. 2 means the
positivity constraint was unable to solve the problem. 3 means
that after positivity constraint failed, also matrix inversion
failed.
"""
self.model = model
self._mapmri_coef = mapmri_coef
self.gtab = model.gtab
self.radial_order = model.radial_order
self.mu = mu
self.R = R
self.lopt = lopt
self.errorcode = errorcode
@property
def mapmri_mu(self):
"""The MAPMRI scale factors"""
return self.mu
@property
def mapmri_R(self):
"""The MAPMRI rotation matrix"""
return self.R
@property
def mapmri_coeff(self):
"""The MAPMRI coefficients"""
return self._mapmri_coef
@warning_for_keywords()
def odf(self, sphere, *, s=2):
"""Calculates the analytical Orientation Distribution Function (ODF)
from the signal.
See :footcite:p:`Ozarslan2013` Eq. (32).
Parameters
----------
sphere : Sphere
A Sphere instance with vertices, edges and faces attributes.
s : unsigned int
radial moment of the ODF
References
----------
.. footbibliography::
"""
if self.model.anisotropic_scaling:
v_ = sphere.vertices
v = np.dot(v_, self.R)
I_s = mapmri_odf_matrix(self.radial_order, self.mu, s, v)
odf = np.dot(I_s, self._mapmri_coef)
else:
_I = self.model.cache_get("ODF_matrix", key=(sphere, s))
if _I is None:
_I = mapmri_isotropic_odf_matrix(
self.radial_order, 1, s, sphere.vertices
)
self.model.cache_set("ODF_matrix", (sphere, s), _I)
odf = self.mu[0] ** s * np.dot(_I, self._mapmri_coef)
return odf
@warning_for_keywords()
def odf_sh(self, *, s=2):
"""Calculates the real analytical odf for a given discrete sphere.
Computes the design matrix of the ODF for the given sphere vertices
and radial moment :footcite:p:`Ozarslan2013` eq. (32). The radial moment
s acts as a sharpening method. The analytical equation for the spherical
ODF basis is given in :footcite:p:`Fick2016b` eq. (C8).
References
----------
.. footbibliography::
"""
if self.model.anisotropic_scaling:
raise ValueError(
"odf in spherical harmonics not yet implemented "
"for anisotropic implementation"
)
_I = self.model.cache_get("ODF_sh_matrix", key=(self.radial_order, s))
if _I is None:
_I = mapmri_isotropic_odf_sh_matrix(self.radial_order, 1, s)
self.model.cache_set("ODF_sh_matrix", (self.radial_order, s), _I)
odf = self.mu[0] ** s * np.dot(_I, self._mapmri_coef)
return odf
def rtpp(self):
"""Calculates the analytical return to the plane probability (RTPP).
RTPP is defined in :footcite:p:`Ozarslan2013` eq. (42). The analytical
formula for the isotropic MAP-MRI basis was derived in
:footcite:p:`Fick2016b` eq. (C11).
References
----------
.. footbibliography::
"""
Bm = self.model.Bm
ind_mat = self.model.ind_mat
if self.model.anisotropic_scaling:
sel = Bm > 0.0 # select only relevant coefficients
const = 1 / (np.sqrt(2 * np.pi) * self.mu[0])
ind_sum = (-1.0) ** (ind_mat[sel, 0] / 2.0)
rtpp_vec = const * Bm[sel] * ind_sum * self._mapmri_coef[sel]
rtpp = rtpp_vec.sum()
return rtpp
else:
rtpp_vec = np.zeros((ind_mat.shape[0]))
count = 0
for n in range(0, self.model.radial_order + 1, 2):
for j in range(1, 2 + n // 2):
ell = n + 2 - 2 * j
const = (-1 / 2.0) ** (ell / 2) / np.sqrt(np.pi)
matsum = 0
for k in range(0, j):
matsum += (
(-1) ** k
* binomialfloat(j + ell - 0.5, j - k - 1)
* gamma(ell / 2 + k + 1 / 2.0)
/ (sfactorial(k) * 0.5 ** (ell / 2 + 1 / 2.0 + k))
)
for _ in range(-ell, ell + 1):
rtpp_vec[count] = const * matsum
count += 1
direction = np.array(self.R[:, 0], ndmin=2)
r, theta, phi = cart2sphere(
direction[:, 0], direction[:, 1], direction[:, 2]
)
rtpp = (
self._mapmri_coef
* (1 / self.mu[0])
* rtpp_vec
* real_sh_descoteaux_from_index(
ind_mat[:, 2], ind_mat[:, 1], theta, phi
)
)
return rtpp.sum()
def rtap(self):
"""Calculates the analytical return to the axis probability (RTAP).
RTAP is defined in :footcite:p:`Ozarslan2013` eq. (40, 44a). The
analytical formula for the isotropic MAP-MRI basis was derived in
:footcite:p:`Fick2016b` eq. (C11).
References
----------
.. footbibliography::
"""
Bm = self.model.Bm
ind_mat = self.model.ind_mat
if self.model.anisotropic_scaling:
sel = Bm > 0.0 # select only relevant coefficients
const = 1 / (2 * np.pi * np.prod(self.mu[1:]))
ind_sum = (-1.0) ** (np.sum(ind_mat[sel, 1:], axis=1) / 2.0)
rtap_vec = const * Bm[sel] * ind_sum * self._mapmri_coef[sel]
rtap = np.sum(rtap_vec)
else:
rtap_vec = np.zeros((ind_mat.shape[0]))
count = 0
for n in range(0, self.model.radial_order + 1, 2):
for j in range(1, 2 + n // 2):
ell = n + 2 - 2 * j
kappa = ((-1) ** (j - 1) * 2 ** (-(ell + 3) / 2.0)) / np.pi
matsum = 0
for k in range(0, j):
matsum += (
(-1) ** k
* binomialfloat(j + ell - 0.5, j - k - 1)
* gamma((ell + 1) / 2.0 + k)
) / (sfactorial(k) * 0.5 ** ((ell + 1) / 2.0 + k))
for _ in range(-ell, ell + 1):
rtap_vec[count] = kappa * matsum
count += 1
rtap_vec *= 2
direction = np.array(self.R[:, 0], ndmin=2)
r, theta, phi = cart2sphere(
direction[:, 0], direction[:, 1], direction[:, 2]
)
rtap_vec = (
self._mapmri_coef
* (1 / self.mu[0] ** 2)
* rtap_vec
* real_sh_descoteaux_from_index(
ind_mat[:, 2], ind_mat[:, 1], theta, phi
)
)
rtap = rtap_vec.sum()
return rtap
def rtop(self):
"""Calculates the analytical return to the origin probability (RTOP).
RTOP is defined in :footcite:p:`Ozarslan2013` eq. (36, 43). The
analytical formula for the isotropic MAP-MRI basis was derived in
:footcite:p:`Fick2016b` eq. (C11).
References
----------
.. footbibliography::
"""
Bm = self.model.Bm
if self.model.anisotropic_scaling:
const = 1 / (np.sqrt(8 * np.pi**3) * np.prod(self.mu))
ind_sum = (-1.0) ** (np.sum(self.model.ind_mat, axis=1) / 2)
rtop_vec = const * ind_sum * Bm * self._mapmri_coef
rtop = rtop_vec.sum()
else:
const = 1 / (2 * np.sqrt(2.0) * np.pi ** (3 / 2.0))
rtop_vec = const * (-1.0) ** (self.model.ind_mat[:, 0] - 1) * Bm
rtop = (1 / self.mu[0] ** 3) * rtop_vec * self._mapmri_coef
rtop = rtop.sum()
return rtop
def msd(self):
"""Calculates the analytical Mean Squared Displacement (MSD).
It is defined as the Laplacian of the origin of the estimated signal
:footcite:p:`Cheng2012`. The analytical formula for the MAP-MRI basis
was derived in :footcite:p:`Fick2016b` eq. (C13, D1).
References
----------
.. footbibliography::
"""
mu = self.mu
ind_mat = self.model.ind_mat
Bm = self.model.Bm
sel = self.model.Bm > 0.0 # select only relevant coefficients
mapmri_coef = self._mapmri_coef[sel]
if self.model.anisotropic_scaling:
ind_sum = np.sum(ind_mat[sel], axis=1)
nx, ny, nz = ind_mat[sel].T
numerator = (
(-1) ** (0.5 * (-ind_sum))
* np.pi ** (3 / 2.0)
* (
(1 + 2 * nx) * mu[0] ** 2
+ (1 + 2 * ny) * mu[1] ** 2
+ (1 + 2 * nz) * mu[2] ** 2
)
)
denominator = (
np.sqrt(
2.0 ** (-ind_sum) * sfactorial(nx) * sfactorial(ny) * sfactorial(nz)
)
* gamma(0.5 - 0.5 * nx)
* gamma(0.5 - 0.5 * ny)
* gamma(0.5 - 0.5 * nz)
)
msd_vec = self._mapmri_coef[sel] * (numerator / denominator)
msd = msd_vec.sum()
else:
msd_vec = (4 * ind_mat[sel, 0] - 1) * Bm[sel]
msd = self.mu[0] ** 2 * msd_vec * mapmri_coef
msd = msd.sum()
return msd
def qiv(self):
"""Calculates the analytical Q-space Inverse Variance (QIV).
It is defined as the inverse of the Laplacian of the origin of the
estimated propagator :footcite:p:`Hosseinbor2013` eq. (22). The
analytical formula for the MAP-MRI basis was derived in
:footcite:p:`Fick2016b` eq. (C14, D2).
References
----------
.. footbibliography::
"""
ux, uy, uz = self.mu
ind_mat = self.model.ind_mat
if self.model.anisotropic_scaling:
sel = self.model.Bm > 0 # select only relevant coefficients
nx, ny, nz = ind_mat[sel].T
numerator = (
8
* np.pi**2
* (ux * uy * uz) ** 3
* np.sqrt(sfactorial(nx) * sfactorial(ny) * sfactorial(nz))
* gamma(0.5 - 0.5 * nx)
* gamma(0.5 - 0.5 * ny)
* gamma(0.5 - 0.5 * nz)
)
denominator = np.sqrt(2.0 ** (-1 + nx + ny + nz)) * (
(1 + 2 * nx) * uy**2 * uz**2
+ ux**2 * ((1 + 2 * nz) * uy**2 + (1 + 2 * ny) * uz**2)
)
qiv_vec = self._mapmri_coef[sel] * (numerator / denominator)
qiv = qiv_vec.sum()
else:
sel = self.model.Bm > 0.0 # select only relevant coefficients
j = ind_mat[sel, 0]
qiv_vec = (8 * (-1.0) ** (1 - j) * np.sqrt(2) * np.pi ** (7 / 2.0)) / (
(4.0 * j - 1) * self.model.Bm[sel]
)
qiv = ux**5 * qiv_vec * self._mapmri_coef[sel]
qiv = qiv.sum()
return qiv
def ng(self):
r"""Calculates the analytical non-Gaussiannity (NG).
For the NG to be meaningful the mapmri scale factors must be estimated
only on data representing Gaussian diffusion of spins, i.e., bvals
smaller than about 2000 s/mm^2 :footcite:p:`Avram2015`.
See :footcite:p:`Ozarslan2013` for a definition of the metric.
References
----------
.. footbibliography::
"""
if self.model.bval_threshold > 2000.0:
warn(
"model bval_threshold must be lower than 2000 for the "
"non_Gaussianity to be physically meaningful [2].",
stacklevel=2,
)
if not self.model.anisotropic_scaling:
raise ValueError(
"Parallel non-Gaussianity is not defined using isotropic scaling."
)
coef = self._mapmri_coef
return np.sqrt(1 - coef[0] ** 2 / np.sum(coef**2))
def ng_parallel(self):
r"""Calculates the analytical parallel non-Gaussiannity (NG).
For the NG to be meaningful the mapmri scale factors must be estimated
only on data representing Gaussian diffusion of spins, i.e., bvals
smaller than about 2000 s/mm^2 :footcite:p:`Avram2015`.
See :footcite:p:`Ozarslan2013` for a definition of the metric.
References
----------
.. footbibliography::
"""
if self.model.bval_threshold > 2000.0:
warn(
"Model bval_threshold must be lower than 2000 for the "
"non_Gaussianity to be physically meaningful [2].",
stacklevel=2,
)
if not self.model.anisotropic_scaling:
raise ValueError(
"Parallel non-Gaussianity is not defined using isotropic scaling."
)
ind_mat = self.model.ind_mat
coef = self._mapmri_coef
a_par = np.zeros_like(coef)
a0 = np.zeros_like(coef)
for i in range(coef.shape[0]):
n1, n2, n3 = ind_mat[i]
if (n2 % 2 + n3 % 2) == 0:
a_par[i] = (
coef[i]
* (-1) ** ((n2 + n3) / 2)
* np.sqrt(sfactorial(n2) * sfactorial(n3))
/ (factorial2(n2) * factorial2(n3))
)
if n1 == 0:
a0[i] = a_par[i]
return np.sqrt(1 - np.sum(a0**2) / np.sum(a_par**2))
def ng_perpendicular(self):
r"""Calculates the analytical perpendicular non-Gaussiannity (NG)
For the NG to be meaningful the mapmri scale factors must be estimated
only on data representing Gaussian diffusion of spins, i.e., bvals
smaller than about 2000 s/mm^2 :footcite:p:`Avram2015`.
See :footcite:p:`Ozarslan2013` for a definition of the metric.
References
----------
.. footbibliography::
"""
if self.model.bval_threshold > 2000.0:
warn(
"model bval_threshold must be lower than 2000 for the "
"non_Gaussianity to be physically meaningful [2].",
stacklevel=2,
)
if not self.model.anisotropic_scaling:
raise ValueError(
"Parallel non-Gaussianity is not defined using isotropic scaling."
)
ind_mat = self.model.ind_mat
coef = self._mapmri_coef
a_perp = np.zeros_like(coef)
a00 = np.zeros_like(coef)
for i in range(coef.shape[0]):
n1, n2, n3 = ind_mat[i]
if n1 % 2 == 0:
if n2 % 2 == 0 and n3 % 2 == 0:
a_perp[i] = (
coef[i]
* (-1) ** (n1 / 2)
* np.sqrt(sfactorial(n1))
/ factorial2(n1)
)
if n2 == 0 and n3 == 0:
a00[i] = a_perp[i]
return np.sqrt(1 - np.sum(a00**2) / np.sum(a_perp**2))
def norm_of_laplacian_signal(self):
"""Calculates the norm of the laplacian of the fitted signal.
This information could be useful to assess if the extrapolation of the
fitted signal contains spurious oscillations. A high laplacian may
indicate that these are present, and any q-space indices that use
integrals of the signal may be corrupted (e.g. RTOP, RTAP, RTPP, QIV).
See :footcite:p:`Fick2016b` for a definition of the metric.
References
----------
.. footbibliography::
"""
if self.model.anisotropic_scaling:
laplacian_matrix = mapmri_laplacian_reg_matrix(
self.model.ind_mat,
self.mu,
self.model.S_mat,
self.model.T_mat,
self.model.U_mat,
)
else:
laplacian_matrix = self.mu[0] * self.model.laplacian_matrix
norm_of_laplacian = np.linalg.multi_dot(
[self._mapmri_coef, laplacian_matrix, self._mapmri_coef]
)
return norm_of_laplacian
@warning_for_keywords()
def fitted_signal(self, *, gtab=None):
"""Recovers the fitted signal for the given gradient table. If no
gradient table is given it recovers the signal for the gtab of the model
object.
"""
if gtab is None:
E = self.predict(self.model.gtab, S0=1.0)
else:
E = self.predict(gtab, S0=1.0)
return E
@warning_for_keywords()
def predict(self, qvals_or_gtab, *, S0=100.0):
"""Recovers the reconstructed signal for any qvalue array or gradient
table.
"""
if isinstance(qvals_or_gtab, np.ndarray):
q = qvals_or_gtab
# qvals = np.linalg.norm(q, axis=1)
else:
gtab = qvals_or_gtab
qvals = np.sqrt(gtab.bvals / self.model.tau) / (2 * np.pi)
q = qvals[:, None] * gtab.bvecs
if self.model.anisotropic_scaling:
q_rot = np.dot(q, self.R)
M = mapmri_phi_matrix(self.radial_order, self.mu, q_rot)
else:
M = mapmri_isotropic_phi_matrix(self.radial_order, self.mu[0], q)
E = S0 * np.dot(M, self._mapmri_coef)
return E
def pdf(self, r_points):
"""Diffusion propagator on a given set of real points.
if the array r_points is non writeable, then intermediate
results are cached for faster recalculation
"""
if self.model.anisotropic_scaling:
r_point_rotated = np.dot(r_points, self.R)
K = mapmri_psi_matrix(self.radial_order, self.mu, r_point_rotated)
EAP = np.dot(K, self._mapmri_coef)
else:
if not r_points.flags.writeable:
K_independent = self.model.cache_get(
"mapmri_matrix_pdf_independent", key=hash(r_points.data)
)
if K_independent is None:
K_independent = mapmri_isotropic_K_mu_independent(
self.radial_order, r_points
)
self.model.cache_set(
"mapmri_matrix_pdf_independent",
hash(r_points.data),
K_independent,
)
K_dependent = mapmri_isotropic_K_mu_dependent(
self.radial_order, self.mu[0], r_points