{ "cells": [ { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n\n# Reconstruction with Constrained Spherical Deconvolution\n\nThis example shows how to use Constrained Spherical Deconvolution (CSD)\nintroduced by Tournier et al. [Tournier2007]_.\n\nThis method is mainly useful with datasets with gradient directions acquired on\na spherical grid.\n\nThe basic idea with this method is that if we could estimate the response\nfunction of a single fiber then we could deconvolve the measured signal and\nobtain the underlying fiber distribution.\n\nIn this way, the reconstruction of the fiber orientation distribution function\n(fODF) in CSD involves two steps:\n 1. Estimation of the fiber response function\n 2. Use the response function to reconstruct the fODF\n\nLet's first load the data. We will use a dataset with 10 b0s and 150 non-b0s\nwith b-value 2000.\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import numpy as np\n\nfrom dipy.core.gradients import gradient_table\nfrom dipy.data import get_fnames, default_sphere\nfrom dipy.io.gradients import read_bvals_bvecs\nfrom dipy.io.image import load_nifti\n\nhardi_fname, hardi_bval_fname, hardi_bvec_fname = get_fnames('stanford_hardi')\n\ndata, affine = load_nifti(hardi_fname)\n\nbvals, bvecs = read_bvals_bvecs(hardi_bval_fname, hardi_bvec_fname)\ngtab = gradient_table(bvals, bvecs)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "You can verify the b-values of the dataset by looking at the attribute\n``gtab.bvals``. Now that a dataset with multiple gradient directions is\nloaded, we can proceed with the two steps of CSD.\n\n## Step 1. Estimation of the fiber response function\n\nThere are many strategies to estimate the fiber response function. Here two\ndifferent strategies are presented.\n\n**Strategy 1 - response function estimates from a local brain region**\nOne simple way to estimate the fiber response function is to look for regions\nof the brain where it is known that there are single coherent fiber\npopulations. For example, if we use a ROI at the center of the brain, we will\nfind single fibers from the corpus callosum. The ``auto_response_ssst``\nfunction will calculate FA for a cuboid ROI of radii equal to ``roi_radii`` in\nthe center of the volume and return the response function estimated in that\nregion for the voxels with FA higher than 0.7.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from dipy.reconst.csdeconv import (auto_response_ssst,\n mask_for_response_ssst,\n response_from_mask_ssst)\n\nresponse, ratio = auto_response_ssst(gtab, data, roi_radii=10, fa_thr=0.7)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that the ``auto_response_ssst`` function calls two functions that can be\nused separately. First, the function ``mask_for_response_ssst`` creates a mask\nof voxels within the cuboid ROI that meet the FA threshold constraint. This\nmask can be used to calculate the number of voxels that were kept, or to also\napply an external mask (a WM mask for example). Second, the function\n``response_from_mask_ssst`` takes the mask and returns the response function\ncalculated within the mask. If no changes are made to the mask between the two\ncalls, the resulting responses should be identical.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "mask = mask_for_response_ssst(gtab, data, roi_radii=10, fa_thr=0.7)\nnvoxels = np.sum(mask)\nprint(nvoxels)\n\nresponse, ratio = response_from_mask_ssst(gtab, data, mask)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The ``response`` tuple contains two elements. The first is an array with\nthe eigenvalues of the response function and the second is the average S0 for\nthis response.\n\nIt is good practice to always validate the result of auto_response_ssst. For\nthis purpose we can print the elements of ``response`` and have a look at their\nvalues.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print(response)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "(array([ 0.0014, 0.00029, 0.00029]), 416.206)\n\nThe tensor generated from the response must be prolate (two smaller eigenvalues\nshould be equal) and look anisotropic with a ratio of second to first\neigenvalue of about 0.2. Or in other words, the axial diffusivity of this\ntensor should be around 5 times larger than the radial diffusivity.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "print(ratio)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "0.21197\n\nWe can double-check that we have a good response function by visualizing the\nresponse function's ODF. Here is how you would do that:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from dipy.viz import window, actor\nfrom dipy.sims.voxel import single_tensor_odf\n\n# Enables/disables interactive visualization\ninteractive = False\n\nscene = window.Scene()\nevals = response[0]\nevecs = np.array([[0, 1, 0], [0, 0, 1], [1, 0, 0]]).T\n\n\nresponse_odf = single_tensor_odf(default_sphere.vertices, evals, evecs)\n# transform our data from 1D to 4D\nresponse_odf = response_odf[None, None, None, :]\nresponse_actor = actor.odf_slicer(response_odf, sphere=default_sphere,\n colormap='plasma')\nscene.add(response_actor)\nprint('Saving illustration as csd_response.png')\nwindow.record(scene, out_path='csd_response.png', size=(200, 200))\nif interactive:\n window.show(scene)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: csd_response.png\n :align: center\n\n Estimated response function.\n\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "scene.rm(response_actor)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Strategy 2 - data-driven calibration of response function** Depending\non the dataset, FA threshold may not be the best way to find the best possible\nresponse function. For one, it depends on the diffusion tensor\n(FA and first eigenvector), which has lower accuracy at high\nb-values. Alternatively, the response function can be calibrated in a\ndata-driven manner [Tax2014]_.\n\nFirst, the data is deconvolved with a 'fat' response function. All voxels that\nare considered to contain only one peak in this deconvolution (as determined by\nthe peak threshold which gives an upper limit of the ratio of the second peak\nto the first peak) are maintained, and from these voxels a new response\nfunction is determined. This process is repeated until convergence is\nreached. Here we calibrate the response function on a small part of the data.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from dipy.reconst.csdeconv import recursive_response" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "A WM mask can shorten computation time for the whole dataset. Here it is\ncreated based on the DTI fit.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import dipy.reconst.dti as dti\ntenmodel = dti.TensorModel(gtab)\ntenfit = tenmodel.fit(data, mask=data[..., 0] > 200)\n\nfrom dipy.reconst.dti import fractional_anisotropy\nFA = fractional_anisotropy(tenfit.evals)\nMD = dti.mean_diffusivity(tenfit.evals)\nwm_mask = (np.logical_or(FA >= 0.4, (np.logical_and(FA >= 0.15, MD >= 0.0011))))\n\nresponse = recursive_response(gtab, data, mask=wm_mask, sh_order=8,\n peak_thr=0.01, init_fa=0.08,\n init_trace=0.0021, iter=8, convergence=0.001,\n parallel=True, num_processes=2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can check the shape of the signal of the response function, which should be\nlike a pancake:\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "response_signal = response.on_sphere(default_sphere)\n# transform our data from 1D to 4D\nresponse_signal = response_signal[None, None, None, :]\nresponse_actor = actor.odf_slicer(response_signal, sphere=default_sphere,\n colormap='plasma')\n\nscene = window.Scene()\n\nscene.add(response_actor)\nprint('Saving illustration as csd_recursive_response.png')\nwindow.record(scene, out_path='csd_recursive_response.png', size=(200, 200))\nif interactive:\n window.show(scene)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: csd_recursive_response.png\n :align: center\n\n Estimated response function using recursive calibration.\n\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "scene.rm(response_actor)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Step 2. fODF reconstruction\n\nAfter estimating a response function for one of the strategies shown above,\nwe are ready to start the deconvolution process. Let's import the CSD model\nand fit the datasets.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from dipy.reconst.csdeconv import ConstrainedSphericalDeconvModel\ncsd_model = ConstrainedSphericalDeconvModel(gtab, response)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "For illustration purposes we will fit only a small portion of the data.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "data_small = data[20:50, 55:85, 38:39]\ncsd_fit = csd_model.fit(data_small)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Show the CSD-based ODFs also known as FODFs (fiber ODFs).\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "csd_odf = csd_fit.odf(default_sphere)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here we visualize only a 30x30 region.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "fodf_spheres = actor.odf_slicer(csd_odf, sphere=default_sphere, scale=0.9,\n norm=False, colormap='plasma')\n\nscene.add(fodf_spheres)\n\nprint('Saving illustration as csd_odfs.png')\nwindow.record(scene, out_path='csd_odfs.png', size=(600, 600))\nif interactive:\n window.show(scene)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: csd_odfs.png\n :align: center\n\n CSD ODFs.\n\nIn DIPY we also provide tools for finding the peak directions (maxima) of the\nODFs. For this purpose we recommend using ``peaks_from_model``.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "from dipy.direction import peaks_from_model\n\ncsd_peaks = peaks_from_model(model=csd_model,\n data=data_small,\n sphere=default_sphere,\n relative_peak_threshold=.5,\n min_separation_angle=25,\n parallel=True,\n num_processes=2)\n\nscene.clear()\nfodf_peaks = actor.peak_slicer(csd_peaks.peak_dirs, csd_peaks.peak_values)\nscene.add(fodf_peaks)\n\nprint('Saving illustration as csd_peaks.png')\nwindow.record(scene, out_path='csd_peaks.png', size=(600, 600))\nif interactive:\n window.show(scene)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: csd_peaks.png\n :align: center\n\n CSD Peaks.\n\nWe can finally visualize both the ODFs and peaks in the same space.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "fodf_spheres.GetProperty().SetOpacity(0.4)\n\nscene.add(fodf_spheres)\n\nprint('Saving illustration as csd_both.png')\nwindow.record(scene, out_path='csd_both.png', size=(600, 600))\nif interactive:\n window.show(scene)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: csd_both.png\n :align: center\n\n CSD Peaks and ODFs.\n\n### References\n\n.. [Tournier2007] J-D. Tournier, F. Calamante and A. Connelly, \"Robust\n determination of the fibre orientation distribution in diffusion MRI:\n Non-negativity constrained super-resolved spherical deconvolution\",\n Neuroimage, vol. 35, no. 4, pp. 1459-1472, 2007.\n\n.. [Tax2014] C.M.W. Tax, B. Jeurissen, S.B. Vos, M.A. Viergever, A. Leemans,\n \"Recursive calibration of the fiber response function for spherical\n deconvolution of diffusion MRI data\", Neuroimage, vol. 86, pp. 67-80, 2014.\n\n.. include:: ../links_names.inc\n\n\n" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.9.13" } }, "nbformat": 4, "nbformat_minor": 0 }