{ "cells": [ { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Denoise images using Non-Local Means (NLMEANS)\n\nUsing the non-local means filter [Coupe08]_ and [Coupe11]_ and you can denoise\n3D or 4D images and boost the SNR of your datasets. You can also decide between\nmodeling the noise as Gaussian or Rician (default).\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import numpy as np\nimport matplotlib.pyplot as plt\nfrom time import time\nfrom dipy.denoise.nlmeans import nlmeans\nfrom dipy.denoise.noise_estimate import estimate_sigma\nfrom dipy.data import get_fnames\nfrom dipy.io.image import load_nifti, save_nifti\n\n\nt1_fname = get_fnames('stanford_t1')\ndata, affine = load_nifti(t1_fname)\n\nmask = data > 1500\n\nprint(\"vol size\", data.shape)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In order to call ``non_local_means`` first you need to estimate the standard\ndeviation of the noise. We have used N=32 since the Stanford dataset was\nacquired on a 3T GE scanner with a 32 array head coil.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "sigma = estimate_sigma(data, N=32)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Calling the main function ``non_local_means``\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "t = time()\n\nden = nlmeans(data, sigma=sigma, mask=mask, patch_radius=1,\n block_radius=2, rician=True)\n\nprint(\"total time\", time() - t)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let us plot the axial slice of the denoised output\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "axial_middle = data.shape[2] // 2\n\nbefore = data[:, :, axial_middle].T\nafter = den[:, :, axial_middle].T\n\ndifference = np.abs(after.astype(np.float64) - before.astype(np.float64))\n\ndifference[~mask[:, :, axial_middle].T] = 0\n\n\nfig, ax = plt.subplots(1, 3)\nax[0].imshow(before, cmap='gray', origin='lower')\nax[0].set_title('before')\nax[1].imshow(after, cmap='gray', origin='lower')\nax[1].set_title('after')\nax[2].imshow(difference, cmap='gray', origin='lower')\nax[2].set_title('difference')\n\nplt.savefig('denoised.png', bbox_inches='tight')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ ".. figure:: denoised.png\n :align: center\n\n **Showing axial slice before (left) and after (right) NLMEANS denoising**\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "save_nifti('denoised.nii.gz', den, affine)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "An improved version of non-local means denoising is adaptive soft coefficient\nmatching, please refer to `example_denoise_ascm` for more details.\n\n## References\n\n.. [Coupe08] P. Coupe, P. Yger, S. Prima, P. Hellier, C. Kervrann, C. Barillot,\n \"An Optimized Blockwise Non Local Means Denoising Filter for 3D Magnetic\n Resonance Images\", IEEE Transactions on Medical Imaging, 27(4):425-441, 2008\n\n.. [Coupe11] Pierrick Coupe, Jose Manjon, Montserrat Robles, Louis Collins.\n \"Adaptive Multiresolution Non-Local Means Filter for 3D MR Image Denoising\"\n IET Image Processing, Institution of Engineering and Technology, 2011\n\n\n.. include:: ../links_names.inc\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 }