Intravoxel incoherent motion¶
The intravoxel incoherent motion (IVIM) model describes diffusion and perfusion in the signal acquired with a diffusion MRI sequence that contains multiple low b-values. The IVIM model can be understood as an adaptation of the work of Stejskal and Tanner [Stejskal65] in biological tissue, and was proposed by Le Bihan [LeBihan84]. The model assumes two compartments: a slow moving compartment, where particles diffuse in a Brownian fashion as a consequence of thermal energy, and a fast moving compartment (the vascular compartment), where blood moves as a consequence of a pressure gradient. In the first compartment, the diffusion coefficient is \(\mathbf{D}\) while in the second compartment, a pseudo diffusion term \(\mathbf{D^*}\) is introduced that describes the displacement of the blood elements in an assumed randomly laid out vascular network, at the macroscopic level. According to [LeBihan84], \(\mathbf{D^*}\) is greater than \(\mathbf{D}\).The IVIM model expresses the MRI signal as follows:
where \(\mathbf{b}\) is the diffusion gradient weighing value (which is dependent on the measurement parameters), \(\mathbf{S_{0}}\) is the signal in the absence of diffusion gradient sensitization, \(\mathbf{f}\) is the perfusion fraction, \(\mathbf{D}\) is the diffusion coefficient and \(\mathbf{D^*}\) is the pseudo-diffusion constant, due to vascular contributions.\[S(b)=S_0(fe^{-bD^*}+(1-f)e^{-bD})\]
In the following example we show how to fit the IVIM model on a diffusion-weighteddataset and visualize the diffusion and pseudo diffusion coefficients. First, we import all relevant modules:
import matplotlib.pyplot as plt
from dipy.reconst.ivim import IvimModel
from dipy.data.fetcher import read_ivim
read_ivim
.
This dataset was acquired with 21 b-values in 3 different directions.
Volumes corresponding to different directions were registered to each
other, and averaged across directions. Thus, this dataset has 4 dimensions,
with the length of the last dimension corresponding to the number
of b-values. In order to use this model the data should contain signals
measured at 0 bvalue.img, gtab = read_ivim()
img
contains a nibabel
NIfTI image object (with the data)
and gtab contains a GradientTable object (information about
the gradients e.g. b-values and b-vectors). We get the
data from img using read_data
.data = img.get_data()
print('data.shape (%d, %d, %d, %d)' % data.shape)
z = 27
b = 20
plt.imshow(data[:, :, z, b].T, origin='lower', cmap='gray',
interpolation='nearest')
plt.axhline(y=100)
plt.axvline(x=170)
plt.savefig("ivim_data_slice.png")
plt.close()
x1, x2 = 160, 180
y1, y2 = 90, 110
data_slice = data[x1:x2, y1:y2, z, :]
plt.imshow(data[x1:x2, y1:y2, z, b].T, origin='lower',
cmap="gray", interpolation='nearest')
plt.savefig("CSF_slice.png")
plt.close()
_estimate_S0_D
is used as the starting point for the non-linear fit
of IVIM parameters using Scipy’s leastsq
or least_square
function
depending on which Scipy version you are using. All initializations for
the model such as split_b
are passed while creating the IvimModel
.
If you are using Scipy 0.17, you can also set bounds by setting
bounds=([0., 0., 0., 0.], [np.inf, 1., 1., 1.]))
while initializing
the IvimModel. It is recommeded that you upgrade to Scipy 0.17 since
the fitting results might at times return values which do not make sense
physically. (For example a negative \(\mathbf{f}\))ivimmodel = IvimModel(gtab)
ivimfit = ivimmodel.fit(data_slice)
ivimparams = ivimfit.model_params
print("ivimparams.shape : {}".format(ivimparams.shape))
i, j = 10, 10
estimated_params = ivimfit.model_params[i, j, :]
print(estimated_params)
estimated_signal = ivimfit.predict(gtab)[i, j, :]
plt.scatter(gtab.bvals, data_slice[i, j, :],
color="green", label="Actual signal")
plt.plot(gtab.bvals, estimated_signal, color="red", label="Estimated Signal")
plt.xlabel("bvalues")
plt.ylabel("Signals")
S0_est, f_est, D_star_est, D_est = estimated_params
text_fit = """Estimated \n S0={:06.3f} f={:06.4f}\n
D*={:06.5f} D={:06.5f}""".format(S0_est, f_est, D_star_est, D_est)
plt.text(0.65, 0.50, text_fit, horizontalalignment='center',
verticalalignment='center', transform=plt.gca().transAxes)
plt.legend(loc='upper right')
plt.savefig("ivim_voxel_plot.png")
plt.close()
variable
to out plotting function which gives the label for
the plot.def plot_map(raw_data, variable, limits, filename):
lower, upper = limits
plt.title('Map for {}'.format(variable))
plt.imshow(raw_data.T, origin='lower', clim=(lower, upper),
cmap="gray", interpolation='nearest')
plt.colorbar()
plt.savefig(filename)
plt.close()
plot_map(ivimfit.S0_predicted, "Predicted S0", (0, 10000), "predicted_S0.png")
plot_map(data_slice[:, :, 0], "Measured S0", (0, 10000), "measured_S0.png")
plot_map(ivimfit.perfusion_fraction, "f", (0, 1), "perfusion_fraction.png")
plot_map(ivimfit.D_star, "D*", (0, 0.01), "perfusion_coeff.png")
plot_map(ivimfit.D, "D", (0, 0.001), "diffusion_coeff.png")
[Stejskal65] | Stejskal, E. O.; Tanner, J. E. (1 January 1965). “Spin Diffusion Measurements: Spin Echoes in the Presence of a Time-Dependent Field Gradient”. The Journal of Chemical Physics 42 (1): 288. Bibcode: 1965JChPh..42..288S. doi:10.1063/1.1695690. |
[LeBihan84] | (1, 2) Le Bihan, Denis, et al. “Separation of diffusion and perfusion in intravoxel incoherent motion MR imaging.” Radiology 168.2 (1988): 497-505. |
Example source code
You can download
the full source code of this example
.
This same script is also included in the dipy source distribution under the
doc/examples/
directory.
No comments:
Post a Comment