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
|
"""
===============================================================================
Suppress Gibbs oscillations
===============================================================================
Magnetic Resonance (MR) images are reconstructed from the Fourier coefficients
of acquired k-space images. Since only a finite number of Fourier coefficients
can be acquired in practice, reconstructed MR images can be corrupted by Gibbs
artefacts, which is manifested by intensity oscillations adjacent to edges of
different tissues types :footcite:p:`Veraart2016a`. Although this artefact
affects MR images in general, in the context of diffusion-weighted imaging,
Gibbs oscillations can be magnified in derived diffusion-based estimates
:footcite:p:`Veraart2016a`, :footcite:p:`Perrone2015`.
In the following example, we show how to suppress Gibbs artefacts of MR images.
This algorithm is based on an adapted version of a sub-voxel Gibbs suppression
procedure :footcite:p:`Kellner2016`. Full details of the implemented algorithm
can be found in chapter 3 of :footcite:p:`NetoHenriques2018` (please cite
:footcite:p:`Kellner2016`, :footcite:p:`NetoHenriques2018` if you are using this
code).
The algorithm to suppress Gibbs oscillations can be imported from the denoise
module of dipy:
"""
import matplotlib.pyplot as plt
import numpy as np
from dipy.data import get_fnames, read_cenir_multib
from dipy.denoise.gibbs import gibbs_removal
from dipy.io.image import load_nifti_data
import dipy.reconst.msdki as msdki
from dipy.segment.mask import median_otsu
###############################################################################
# We first apply this algorithm to a T1-weighted dataset which can be fetched
# using the following code:
t1_fname, t1_denoised_fname, ap_fname = get_fnames(name="tissue_data")
t1 = load_nifti_data(t1_denoised_fname)
###############################################################################
# Let's plot a slice of this dataset.
axial_slice = 88
t1_slice = t1[..., axial_slice]
fig = plt.figure(figsize=(15, 4))
fig.subplots_adjust(wspace=0.2)
t1_slice = np.rot90(t1_slice)
plt.subplot(1, 2, 1)
plt.imshow(t1_slice, cmap="gray", vmin=100, vmax=400)
plt.colorbar()
fig.savefig("structural.png")
###############################################################################
# .. rst-class:: centered small fst-italic fw-semibold
#
# Representative slice of a T1-weighted structural image.
#
#
# Due to the high quality of the data, Gibbs artefacts are not visually
# evident in this dataset. Therefore, to analyse the benefits of the Gibbs
# suppression algorithm, Gibbs artefacts are artificially introduced by
# removing high frequencies of the image's Fourier transform.
c = np.fft.fft2(t1_slice)
c = np.fft.fftshift(c)
N = c.shape[0]
c_crop = c[64:192, 64:192]
N = c_crop.shape[0]
t1_gibbs = abs(np.fft.ifft2(c_crop) / 4)
###############################################################################
# Gibbs oscillation suppression of this single data slice can be performed by
# running the following command:
t1_unring = gibbs_removal(t1_gibbs, inplace=False)
###############################################################################
# Let’s plot the results:
fig1, ax = plt.subplots(1, 3, figsize=(12, 6), subplot_kw={"xticks": [], "yticks": []})
ax.flat[0].imshow(t1_gibbs, cmap="gray", vmin=100, vmax=400)
ax.flat[0].annotate(
"Rings",
fontsize=10,
xy=(81, 70),
color="red",
xycoords="data",
xytext=(30, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
ax.flat[0].annotate(
"Rings",
fontsize=10,
xy=(75, 50),
color="red",
xycoords="data",
xytext=(30, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
ax.flat[1].imshow(t1_unring, cmap="gray", vmin=100, vmax=400)
ax.flat[1].annotate(
"WM/GM",
fontsize=10,
xy=(75, 50),
color="green",
xycoords="data",
xytext=(30, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "green"},
)
ax.flat[2].imshow(t1_unring - t1_gibbs, cmap="gray", vmin=0, vmax=10)
ax.flat[2].annotate(
"Rings",
fontsize=10,
xy=(81, 70),
color="red",
xycoords="data",
xytext=(30, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
ax.flat[2].annotate(
"Rings",
fontsize=10,
xy=(75, 50),
color="red",
xycoords="data",
xytext=(30, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
plt.show()
fig1.savefig("Gibbs_suppression_structural.png")
###############################################################################
# .. rst-class:: centered small fst-italic fw-semibold
#
# Uncorrected and corrected structural images are shown in the left
# and middle panels, while the difference between these images is shown
# in the right panel.
#
#
# The image artificially corrupted with Gibb's artefacts is shown in the left
# panel. In this panel, the characteristic ringing profile of Gibbs artefacts
# can be visually appreciated (see intensity oscillations pointed by the red
# arrows). The corrected image is shown in the middle panel. One can appreciate
# that artefactual oscillations are visually suppressed without compromising
# the contrast between white and grey matter (e.g. details pointed by the green
# arrow). The difference between uncorrected and corrected data is plotted in
# the right panel which highlights the suppressed Gibbs ringing profile.
#
#
# Now let's show how to use the Gibbs suppression algorithm in
# diffusion-weighted images. We fetch the multi-shell diffusion-weighted
# dataset which was kindly supplied by Romain Valabrègue,
# CENIR, ICM, Paris :footcite:p:`Valabregue2015`.
bvals = [200, 400, 1000, 2000]
img, gtab = read_cenir_multib(bvals=bvals)
data = np.asarray(img.dataobj)
###############################################################################
# For illustration purposes, we select two slices of this dataset
data_slices = data[:, :, 40:42, :]
###############################################################################
# Gibbs oscillation suppression of all multi-shell data and all slices
# can be performed in the following way:
data_corrected = gibbs_removal(data_slices, slice_axis=2, num_processes=-1)
###############################################################################
# Due to the high dimensionality of diffusion-weighted data, we recommend
# that you specify which is the axis of data matrix that corresponds to
# different slices in the above step. This is done by using the optional
# parameter 'slice_axis'.
#
# Below we plot the results for an image acquired with b-value=0:
fig2, ax = plt.subplots(1, 3, figsize=(12, 6), subplot_kw={"xticks": [], "yticks": []})
ax.flat[0].imshow(
data_slices[:, :, 0, 0].T, cmap="gray", origin="lower", vmin=0, vmax=10000
)
ax.flat[0].set_title("Uncorrected b0")
ax.flat[1].imshow(
data_corrected[:, :, 0, 0].T, cmap="gray", origin="lower", vmin=0, vmax=10000
)
ax.flat[1].set_title("Corrected b0")
ax.flat[2].imshow(
data_corrected[:, :, 0, 0].T - data_slices[:, :, 0, 0].T,
cmap="gray",
origin="lower",
vmin=-500,
vmax=500,
)
ax.flat[2].set_title("Gibbs residuals")
plt.show()
fig2.savefig("Gibbs_suppression_b0.png")
###############################################################################
# .. rst-class:: centered small fst-italic fw-semibold
#
# Uncorrected (left panel) and corrected (middle panel) b-value=0 images. For
# reference, the difference between uncorrected and corrected images is shown
# in the right panel.
#
#
# The above figure shows that the benefits of suppressing Gibbs artefacts is
# hard to observe on b-value=0 data. Therefore, diffusion derived metrics for
# both uncorrected and corrected data are computed using the mean signal
# diffusion kurtosis image technique
# (:ref:`sphx_glr_examples_built_reconstruction_reconst_msdki.py`).
#
# To avoid unnecessary calculations on the background of the image, we also
# compute a brain mask.
# Create a brain mask
maskdata, mask = median_otsu(
data_slices,
vol_idx=range(10, 50),
median_radius=3,
numpass=1,
autocrop=False,
dilate=1,
)
# Define mean signal diffusion kurtosis model
dki_model = msdki.MeanDiffusionKurtosisModel(gtab)
# Fit the uncorrected data
dki_fit = dki_model.fit(data_slices, mask=mask)
MSKini = dki_fit.msk
# Fit the corrected data
dki_fit = dki_model.fit(data_corrected, mask=mask)
MSKgib = dki_fit.msk
###############################################################################
# Let's plot the results
fig3, ax = plt.subplots(1, 3, figsize=(12, 12), subplot_kw={"xticks": [], "yticks": []})
ax.flat[0].imshow(MSKini[:, :, 0].T, cmap="gray", origin="lower", vmin=0, vmax=1.5)
ax.flat[0].set_title("MSK (uncorrected)")
ax.flat[0].annotate(
"Rings",
fontsize=12,
xy=(59, 63),
color="red",
xycoords="data",
xytext=(45, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
ax.flat[1].imshow(MSKgib[:, :, 0].T, cmap="gray", origin="lower", vmin=0, vmax=1.5)
ax.flat[1].set_title("MSK (corrected)")
ax.flat[2].imshow(
MSKgib[:, :, 0].T - MSKini[:, :, 0].T,
cmap="gray",
origin="lower",
vmin=-0.2,
vmax=0.2,
)
ax.flat[2].set_title("MSK (uncorrected - corrected")
ax.flat[2].annotate(
"Rings",
fontsize=12,
xy=(59, 63),
color="red",
xycoords="data",
xytext=(45, 0),
textcoords="offset points",
arrowprops={"arrowstyle": "->", "color": "red"},
)
plt.show()
fig3.savefig("Gibbs_suppression_msdki.png")
###############################################################################
# .. rst-class:: centered small fst-italic fw-semibold
#
# Uncorrected and corrected mean signal kurtosis images are shown in the left
# and middle panels. The difference between uncorrected and corrected images
# are show in the right panel.
#
#
# In the left panel of the figure above, Gibbs artefacts can be appreciated by
# the negative values of mean signal kurtosis (black voxels) adjacent to the
# brain ventricle (red arrows). These negative values seem to be suppressed
# after the `gibbs_removal` function is applied. For a better visualization of
# Gibbs oscillations, the difference between corrected and uncorrected images
# are shown in the right panel.
#
#
# References
# ----------
#
# .. footbibliography::
#
|