denoise

Module: denoise.adaptive_soft_matching

adaptive_soft_matching(ima, fimau, fimao, sigma)

Adaptive Soft Coefficient Matching

Module: denoise.gibbs

gibbs_removal(vol[, slice_axis, n_points, ...])

Suppresses Gibbs ringing artefacts of images volumes.

Module: denoise.localpca

genpca(arr[, sigma, mask, patch_radius, ...])

General function to perform PCA-based denoising of diffusion datasets.

localpca(arr[, sigma, mask, patch_radius, ...])

Performs local PCA denoising according to Manjon et al. [1]_.

mppca(arr[, mask, patch_radius, pca_method, ...])

Performs PCA-based denoising using the Marcenko-Pastur distribution [1]_.

Module: denoise.nlmeans

nlmeans(arr, sigma[, mask, patch_radius, ...])

Non-local means for denoising 3D and 4D images

Module: denoise.noise_estimate

piesno(data, N[, alpha, l, itermax, eps, ...])

Probabilistic Identification and Estimation of Noise (PIESNO).

estimate_sigma(arr[, ...])

Standard deviation estimation from local patches

Module: denoise.non_local_means

non_local_means(arr, sigma[, mask, ...])

Non-local means for denoising 3D and 4D images, using

Module: denoise.patch2self

patch2self(data, bvals[, patch_radius, ...])

Patch2Self Denoiser.

adaptive_soft_matching

dipy.denoise.adaptive_soft_matching.adaptive_soft_matching(ima, fimau, fimao, sigma)

Adaptive Soft Coefficient Matching

Combines two filtered 3D-images at different resolutions and the original image. Returns the resulting combined image.

Parameters

ima : the original (not filtered) image fimau : 3D double array,

filtered image with optimized non-local means using a small block (suggested:3x3), which corresponds to a “high resolution” filter.

fimao3D double array,

filtered image with optimized non-local means using a small block (suggested:5x5), which corresponds to a “low resolution” filter.

sigmathe estimated standard deviation of the Gaussian random variables

that explain the rician noise. Note: In P. Coupe et al. the rician noise was simulated as sqrt((f+x)^2 + (y)^2) where f is the pixel value and x and y are independent realizations of a random variable with Normal distribution, with mean=0 and standard deviation=h

Returns

fima3D double array

output denoised array which is of the same shape as that of the input

References

[Coupe11]

Pierrick Coupe, Jose Manjon, Montserrat Robles, Louis Collins. “Multiresolution Non-Local Means Filter for 3D MR Image Denoising” IET Image Processing, Institution of Engineering and Technology, 2011

gibbs_removal

dipy.denoise.gibbs.gibbs_removal(vol, slice_axis=2, n_points=3, inplace=True, num_processes=1)

Suppresses Gibbs ringing artefacts of images volumes.

Parameters

volndarray ([X, Y]), ([X, Y, Z]) or ([X, Y, Z, g])

Matrix containing one volume (3D) or multiple (4D) volumes of images.

slice_axisint (0, 1, or 2)

Data axis corresponding to the number of acquired slices. Default is set to the third axis.

n_pointsint, optional

Number of neighbour points to access local TV (see note). Default is set to 3.

inplacebool, optional

If True, the input data is replaced with results. Otherwise, returns a new array. Default is set to True.

num_processesint or None, optional

Split the calculation to a pool of children processes. This only applies to 3D or 4D data arrays. Default is 1. If < 0 the maximal number of cores minus num_processes + 1 is used (enter -1 to use as many cores as possible). 0 raises an error.

Returns

volndarray ([X, Y]), ([X, Y, Z]) or ([X, Y, Z, g])

Matrix containing one volume (3D) or multiple (4D) volumes of corrected images.

Notes

For 4D matrix last element should always correspond to the number of diffusion gradient directions.

References

Please cite the following articles .. [1] Neto Henriques, R., 2018. Advanced Methods for Diffusion MRI Data

Analysis and their Application to the Healthy Ageing Brain (Doctoral thesis). https://doi.org/10.17863/CAM.29356

genpca

dipy.denoise.localpca.genpca(arr, sigma=None, mask=None, patch_radius=2, pca_method='eig', tau_factor=None, return_sigma=False, out_dtype=None, suppress_warning=False)

General function to perform PCA-based denoising of diffusion datasets.

Parameters

arr4D array

Array of data to be denoised. The dimensions are (X, Y, Z, N), where N are the diffusion gradient directions. The first 3 dimension must have size >= 2 * patch_radius + 1 or size = 1.

sigmafloat or 3D array (optional)

Standard deviation of the noise estimated from the data. If no sigma is given, this will be estimated based on random matrix theory [1]_,[2]_

mask3D boolean array (optional)

A mask with voxels that are true inside the brain and false outside of it. The function denoises within the true part and returns zeros outside of those voxels.

patch_radiusint or 1D array (optional)

The radius of the local patch to be taken around each voxel (in voxels). E.g. patch_radius=2 gives 5x5x5 patches.

pca_method‘eig’ or ‘svd’ (optional)

Use either eigenvalue decomposition (eig) or singular value decomposition (svd) for principal component analysis. The default method is ‘eig’ which is faster. However, occasionally ‘svd’ might be more accurate.

tau_factorfloat (optional)

Thresholding of PCA eigenvalues is done by nulling out eigenvalues that are smaller than:

\[\tau = (\tau_{factor} \sigma)^2\]

tau_{factor} can be set to a predefined values (e.g. tau_{factor} = 2.3 [3]), or automatically calculated using random matrix theory (in case that tau_{factor} is set to None).

return_sigmabool (optional)

If true, the Standard deviation of the noise will be returned.

out_dtypestr or dtype (optional)

The dtype for the output array. Default: output has the same dtype as the input.

suppress_warningbool (optional)

If true, suppress warning caused by patch_size < arr.shape[-1].

Returns

denoised_arr4D array

This is the denoised array of the same size as that of the input data, clipped to non-negative values.

References

localpca

dipy.denoise.localpca.localpca(arr, sigma=None, mask=None, patch_radius=2, gtab=None, patch_radius_sigma=1, pca_method='eig', tau_factor=2.3, return_sigma=False, correct_bias=True, out_dtype=None, suppress_warning=False)

Performs local PCA denoising according to Manjon et al. [1]_.

Parameters

arr4D array

Array of data to be denoised. The dimensions are (X, Y, Z, N), where N are the diffusion gradient directions.

sigmafloat or 3D array (optional)

Standard deviation of the noise estimated from the data. If not given, calculate using method in [1]_.

mask3D boolean array (optional)

A mask with voxels that are true inside the brain and false outside of it. The function denoises within the true part and returns zeros outside of those voxels.

patch_radiusint or 1D array (optional)

The radius of the local patch to be taken around each voxel (in voxels). E.g. patch_radius=2 gives 5x5x5 patches.

gtab: gradient table object (optional if sigma is provided)

gradient information for the data gives us the bvals and bvecs of diffusion data, which is needed to calculate noise level if sigma is not provided.

patch_radius_sigmaint (optional)

The radius of the local patch to be taken around each voxel (in voxels) for estimating sigma. E.g. patch_radius_sigma=2 gives 5x5x5 patches.

pca_method‘eig’ or ‘svd’ (optional)

Use either eigenvalue decomposition (eig) or singular value decomposition (svd) for principal component analysis. The default method is ‘eig’ which is faster. However, occasionally ‘svd’ might be more accurate.

tau_factorfloat (optional)

Thresholding of PCA eigenvalues is done by nulling out eigenvalues that are smaller than:

\[\tau = (\tau_{factor} \sigma)^2\]

tau_{factor} can be change to adjust the relationship between the noise standard deviation and the threshold tau. If tau_{factor} is set to None, it will be automatically calculated using the Marcenko-Pastur distribution [2]_. Default: 2.3 according to [1]_.

return_sigmabool (optional)

If true, a noise standard deviation estimate based on the Marcenko-Pastur distribution is returned [2]_.

correct_biasbool (optional)

Whether to correct for bias due to Rician noise. This is an implementation of equation 8 in [1]_.

out_dtypestr or dtype (optional)

The dtype for the output array. Default: output has the same dtype as the input.

suppress_warningbool (optional)

If true, suppress warning caused by patch_size < arr.shape[-1].

Returns

denoised_arr4D array

This is the denoised array of the same size as that of the input data, clipped to non-negative values

References

mppca

dipy.denoise.localpca.mppca(arr, mask=None, patch_radius=2, pca_method='eig', return_sigma=False, out_dtype=None, suppress_warning=False)

Performs PCA-based denoising using the Marcenko-Pastur distribution [1]_.

Parameters

arr4D array

Array of data to be denoised. The dimensions are (X, Y, Z, N), where N are the diffusion gradient directions.

mask3D boolean array (optional)

A mask with voxels that are true inside the brain and false outside of it. The function denoises within the true part and returns zeros outside of those voxels.

patch_radiusint or 1D array (optional)

The radius of the local patch to be taken around each voxel (in voxels). E.g. patch_radius=2 gives 5x5x5 patches.

pca_method‘eig’ or ‘svd’ (optional)

Use either eigenvalue decomposition (eig) or singular value decomposition (svd) for principal component analysis. The default method is ‘eig’ which is faster. However, occasionally ‘svd’ might be more accurate.

return_sigmabool (optional)

If true, a noise standard deviation estimate based on the Marcenko-Pastur distribution is returned [2]_.

out_dtypestr or dtype (optional)

The dtype for the output array. Default: output has the same dtype as the input.

suppress_warningbool (optional)

If true, suppress warning caused by patch_size < arr.shape[-1].

Returns

denoised_arr4D array

This is the denoised array of the same size as that of the input data, clipped to non-negative values

sigma3D array (when return_sigma=True)

Estimate of the spatial varying standard deviation of the noise

References

nlmeans

dipy.denoise.nlmeans.nlmeans(arr, sigma, mask=None, patch_radius=1, block_radius=5, rician=True, num_threads=None)

Non-local means for denoising 3D and 4D images

Parameters

arr3D or 4D ndarray

The array to be denoised

mask : 3D ndarray sigma : float or 3D array

standard deviation of the noise estimated from the data

patch_radiusint

patch size is 2 x patch_radius + 1. Default is 1.

block_radiusint

block size is 2 x block_radius + 1. Default is 5.

ricianboolean

If True the noise is estimated as Rician, otherwise Gaussian noise is assumed.

num_threadsint, optional

Number of threads to be used for OpenMP parallelization. If None (default) the value of OMP_NUM_THREADS environment variable is used if it is set, otherwise all available threads are used. If < 0 the maximal number of threads minus |num_threads + 1| is used (enter -1 to use as many threads as possible). 0 raises an error.

Returns

denoised_arrndarray

the denoised arr which has the same shape as arr.

References

[Descoteaux08]

Descoteaux, Maxime and Wiest-Daesslé, Nicolas and Prima, Sylvain and Barillot, Christian and Deriche, Rachid Impact of Rician Adapted Non-Local Means Filtering on HARDI, MICCAI 2008

piesno

dipy.denoise.noise_estimate.piesno(data, N, alpha=0.01, l=100, itermax=100, eps=1e-05, return_mask=False)

Probabilistic Identification and Estimation of Noise (PIESNO).

Parameters

datandarray

The magnitude signals to analyse. The last dimension must contain the same realisation of the volume, such as dMRI or fMRI data.

Nint

The number of phase array coils of the MRI scanner. If your scanner does a SENSE reconstruction, ALWAYS use N=1, as the noise profile is always Rician. If your scanner does a GRAPPA reconstruction, set N as the number of phase array coils.

alphafloat

Probabilistic estimation threshold for the gamma function.

lint

number of initial estimates for sigma to try.

itermaxint

Maximum number of iterations to execute if convergence is not reached.

epsfloat

Tolerance for the convergence criterion. Convergence is reached if two subsequent estimates are smaller than eps.

return_maskbool

If True, return a mask identifying all the pure noise voxel that were found.

Returns

sigmafloat

The estimated standard deviation of the gaussian noise.

maskndarray (optional)

A boolean mask indicating the voxels identified as pure noise.

Notes

This function assumes two things : 1. The data has a noisy, non-masked background and 2. The data is a repetition of the same measurements along the last axis, i.e. dMRI or fMRI data, not structural data like T1/T2.

This function processes the data slice by slice, as originally designed in the paper. Use it to get a slice by slice estimation of the noise, as in spinal cord imaging for example.

References

“Probabilistic Identification and Estimation of Noise (PIESNO): A self-consistent approach and its applications in MRI.” Journal of Magnetic Resonance 2009; 199: 94-103.

“A signal transformational framework for breaking the noise floor and its applications in MRI.” Journal of Magnetic Resonance 2009; 197: 108-119.

estimate_sigma

dipy.denoise.noise_estimate.estimate_sigma(arr, disable_background_masking=False, N=0)

Standard deviation estimation from local patches

Parameters

arr3D or 4D ndarray

The array to be estimated

disable_background_maskingbool, default False

If True, uses all voxels for the estimation, otherwise, only non-zeros voxels are used. Useful if the background is masked by the scanner.

Nint, default 0

Number of coils of the receiver array. Use N = 1 in case of a SENSE reconstruction (Philips scanners) or the number of coils for a GRAPPA reconstruction (Siemens and GE). Use 0 to disable the correction factor, as for example if the noise is Gaussian distributed. See [1] for more information.

Returns

sigmandarray

standard deviation of the noise, one estimation per volume.

Notes

This function is the same as manually taking the standard deviation of the background and gives one value for the whole 3D array. It also includes the coil-dependent correction factor of Koay 2006 (see [1]_, equation 18) with theta = 0. Since this function was introduced in [2]_ for T1 imaging, it is expected to perform ok on diffusion MRI data, but might oversmooth some regions and leave others un-denoised for spatially varying noise profiles. Consider using piesno() to estimate sigma instead if visual inaccuracies are apparent in the denoised result.

References

scheme for signal extraction from noisy magnitude MR signals. Journal of Magnetic Resonance), 179(2), 317-22.

C., 2008. An optimized blockwise nonlocal means denoising filter for 3-D magnetic resonance images, IEEE Trans. Med. Imaging 27, 425-41.

non_local_means

dipy.denoise.non_local_means.non_local_means(arr, sigma, mask=None, patch_radius=1, block_radius=5, rician=True)
Non-local means for denoising 3D and 4D images, using

blockwise averaging approach

Parameters

arr3D or 4D ndarray

The array to be denoised

mask : 3D ndarray sigma : float

standard deviation of the noise estimated from the data

patch_radiusint

patch size is 2 x patch_radius + 1. Default is 1.

block_radiusint

block size is 2 x block_radius + 1. Default is 5.

ricianboolean

If True the noise is estimated as Rician, otherwise Gaussian noise is assumed.

Returns

denoised_arrndarray

the denoised arr which has the same shape as arr.

References

[Coupe08]

P. Coupe, P. Yger, S. Prima, P. Hellier, C. Kervrann, C. Barillot, An Optimized Blockwise Non Local Means Denoising Filter for 3D Magnetic Resonance Images, IEEE Transactions on Medical Imaging, 27(4):425-441, 2008

[Coupe11]

Pierrick Coupe, Jose Manjon, Montserrat Robles, Louis Collins. Adaptive Multiresolution Non-Local Means Filter for 3D MR Image Denoising IET Image Processing, Institution of Engineering and Technology, 2011

patch2self

dipy.denoise.patch2self.patch2self(data, bvals, patch_radius=(0, 0, 0), model='ols', b0_threshold=50, out_dtype=None, alpha=1.0, verbose=False, b0_denoising=True, clip_negative_vals=False, shift_intensity=True)

Patch2Self Denoiser.

Parameters

datandarray

The 4D noisy DWI data to be denoised.

bvals1D array

Array of the bvals from the DWI acquisition

patch_radiusint or 1D array, optional

The radius of the local patch to be taken around each voxel (in voxels). Default: 0 (denoise in blocks of 1x1x1 voxels).

modelstring, or initialized linear model object.

This will determine the algorithm used to solve the set of linear equations underlying this model. If it is a string it needs to be one of the following: {‘ols’, ‘ridge’, ‘lasso’}. Otherwise, it can be an object that inherits from dipy.optimize.SKLearnLinearSolver or an object with a similar interface from Scikit-Learn: sklearn.linear_model.LinearRegression, sklearn.linear_model.Lasso or sklearn.linear_model.Ridge and other objects that inherit from sklearn.base.RegressorMixin. Default: ‘ols’.

b0_thresholdint, optional

Threshold for considering volumes as b0.

out_dtypestr or dtype, optional

The dtype for the output array. Default: output has the same dtype as the input.

alphafloat, optional

Regularization parameter only for ridge regression model.

verbosebool, optional

Show progress of Patch2Self and time taken.

b0_denoisingbool, optional

Skips denoising b0 volumes if set to False.

clip_negative_valsbool, optional

Sets negative values after denoising to 0 using np.clip.

shift_intensitybool, optional

Shifts the distribution of intensities per volume to give non-negative values

Returns

denoised arrayndarray

This is the denoised array of the same size as that of the input data, clipped to non-negative values.

References

[Fadnavis20] S. Fadnavis, J. Batson, E. Garyfallidis, Patch2Self:

Denoising Diffusion MRI with Self-supervised Learning, Advances in Neural Information Processing Systems 33 (2020)