denoise

bench Run benchmarks for module using nose.
test Run tests for module using nose.

Module: denoise.adaptive_soft_matching

adaptive_soft_matching(ima, fimau, fimao, sigma) Adaptive Soft Coefficient Matching

Module: denoise.localpca

eigh(a[, b, lower, eigvals_only, …]) Solve an ordinary or generalized eigenvalue problem for a complex Hermitian or real symmetric matrix.
localpca(arr, sigma[, mask, pca_method, …]) Local PCA-based denoising of diffusion datasets.

Module: denoise.nlmeans

nlmeans(arr, sigma[, mask, patch_radius, …]) Non-local means for denoising 3D and 4D images
nlmeans_3d Non-local means for denoising 3D images

Module: denoise.noise_estimate

convolve(input, weights[, output, mode, …]) Multidimensional convolution.
estimate_sigma(arr[, …]) Standard deviation estimation from local patches
piesno(data, N[, alpha, l, itermax, eps, …]) Probabilistic Identification and Estimation of Noise (PIESNO).

Module: denoise.non_local_means

nlmeans_block Non-Local Means Denoising Using Blockwise Averaging
non_local_means(arr, sigma[, mask, …]) Non-local means for denoising 3D and 4D images, using

bench

dipy.denoise.bench(label='fast', verbose=1, extra_argv=None)

Run benchmarks for module using nose.

Parameters:
label : {‘fast’, ‘full’, ‘’, attribute identifier}, optional

Identifies the benchmarks to run. This can be a string to pass to the nosetests executable with the ‘-A’ option, or one of several special values. Special values are: * ‘fast’ - the default - which corresponds to the nosetests -A

option of ‘not slow’.

  • ‘full’ - fast (as above) and slow benchmarks as in the ‘no -A’ option to nosetests - this is the same as ‘’.
  • None or ‘’ - run all tests.

attribute_identifier - string passed directly to nosetests as ‘-A’.

verbose : int, optional

Verbosity value for benchmark outputs, in the range 1-10. Default is 1.

extra_argv : list, optional

List with any extra arguments to pass to nosetests.

Returns:
success : bool

Returns True if running the benchmarks works, False if an error occurred.

Notes

Benchmarks are like tests, but have names starting with “bench” instead of “test”, and can be found under the “benchmarks” sub-directory of the module.

Each NumPy module exposes bench in its namespace to run all benchmarks for it.

Examples

>>> success = np.lib.bench() 
Running benchmarks for numpy.lib
...
using 562341 items:
unique:
0.11
unique1d:
0.11
ratio: 1.0
nUnique: 56230 == 56230
...
OK
>>> success 
True

test

dipy.denoise.test(label='fast', verbose=1, extra_argv=None, doctests=False, coverage=False, raise_warnings=None, timer=False)

Run tests for module using nose.

Parameters:
label : {‘fast’, ‘full’, ‘’, attribute identifier}, optional

Identifies the tests to run. This can be a string to pass to the nosetests executable with the ‘-A’ option, or one of several special values. Special values are: * ‘fast’ - the default - which corresponds to the nosetests -A

option of ‘not slow’.

  • ‘full’ - fast (as above) and slow tests as in the ‘no -A’ option to nosetests - this is the same as ‘’.
  • None or ‘’ - run all tests.

attribute_identifier - string passed directly to nosetests as ‘-A’.

verbose : int, optional

Verbosity value for test outputs, in the range 1-10. Default is 1.

extra_argv : list, optional

List with any extra arguments to pass to nosetests.

doctests : bool, optional

If True, run doctests in module. Default is False.

coverage : bool, optional

If True, report coverage of NumPy code. Default is False. (This requires the `coverage module:

raise_warnings : None, str or sequence of warnings, optional

This specifies which warnings to configure as ‘raise’ instead of being shown once during the test execution. Valid strings are:

  • “develop” : equals (Warning,)
  • “release” : equals (), don’t raise on any warnings.

The default is to use the class initialization value.

timer : bool or int, optional

Timing of individual tests with nose-timer (which needs to be installed). If True, time tests and report on all of them. If an integer (say N), report timing results for N slowest tests.

Returns:
result : object

Returns the result of running the tests as a nose.result.TextTestResult object.

Notes

Each NumPy module exposes test in its namespace to run all tests for it. For example, to run all tests for numpy.lib:

>>> np.lib.test() 

Examples

>>> result = np.lib.test() 
Running unit tests for numpy.lib
...
Ran 976 tests in 3.933s

OK

>>> result.errors 
[]
>>> result.knownfail 
[]

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 orginal 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.

fimao : 3D double array,

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

sigma : the 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:
fima : 3D 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

eigh

dipy.denoise.localpca.eigh(a, b=None, lower=True, eigvals_only=False, overwrite_a=False, overwrite_b=False, turbo=True, eigvals=None, type=1, check_finite=True)

Solve an ordinary or generalized eigenvalue problem for a complex Hermitian or real symmetric matrix.

Find eigenvalues w and optionally eigenvectors v of matrix a, where b is positive definite:

              a v[:,i] = w[i] b v[:,i]
v[i,:].conj() a v[:,i] = w[i]
v[i,:].conj() b v[:,i] = 1
Parameters:
a : (M, M) array_like

A complex Hermitian or real symmetric matrix whose eigenvalues and eigenvectors will be computed.

b : (M, M) array_like, optional

A complex Hermitian or real symmetric definite positive matrix in. If omitted, identity matrix is assumed.

lower : bool, optional

Whether the pertinent array data is taken from the lower or upper triangle of a. (Default: lower)

eigvals_only : bool, optional

Whether to calculate only eigenvalues and no eigenvectors. (Default: both are calculated)

turbo : bool, optional

Use divide and conquer algorithm (faster but expensive in memory, only for generalized eigenvalue problem and if eigvals=None)

eigvals : tuple (lo, hi), optional

Indexes of the smallest and largest (in ascending order) eigenvalues and corresponding eigenvectors to be returned: 0 <= lo <= hi <= M-1. If omitted, all eigenvalues and eigenvectors are returned.

type : int, optional

Specifies the problem type to be solved:

type = 1: a v[:,i] = w[i] b v[:,i]

type = 2: a b v[:,i] = w[i] v[:,i]

type = 3: b a v[:,i] = w[i] v[:,i]

overwrite_a : bool, optional

Whether to overwrite data in a (may improve performance)

overwrite_b : bool, optional

Whether to overwrite data in b (may improve performance)

check_finite : bool, optional

Whether to check that the input matrices contain only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination) if the inputs do contain infinities or NaNs.

Returns:
w : (N,) float ndarray

The N (1<=N<=M) selected eigenvalues, in ascending order, each repeated according to its multiplicity.

v : (M, N) complex ndarray

(if eigvals_only == False)

The normalized selected eigenvector corresponding to the eigenvalue w[i] is the column v[:,i].

Normalization:

type 1 and 3: v.conj() a v = w

type 2: inv(v).conj() a inv(v) = w

type = 1 or 2: v.conj() b v = I

type = 3: v.conj() inv(b) v = I

Raises:
LinAlgError

If eigenvalue computation does not converge, an error occurred, or b matrix is not definite positive. Note that if input matrices are not symmetric or hermitian, no error is reported but results will be wrong.

See also

eigvalsh
eigenvalues of symmetric or Hermitian arrays
eig
eigenvalues and right eigenvectors for non-symmetric arrays
eigh
eigenvalues and right eigenvectors for symmetric/Hermitian arrays
eigh_tridiagonal
eigenvalues and right eiegenvectors for symmetric/Hermitian tridiagonal matrices

Notes

This function does not check the input array for being hermitian/symmetric in order to allow for representing arrays with only their upper/lower triangular parts.

Examples

>>> from scipy.linalg import eigh
>>> A = np.array([[6, 3, 1, 5], [3, 0, 5, 1], [1, 5, 6, 2], [5, 1, 2, 2]])
>>> w, v = eigh(A)
>>> np.allclose(A @ v - v @ np.diag(w), np.zeros((4, 4)))
True

localpca

dipy.denoise.localpca.localpca(arr, sigma, mask=None, pca_method='eig', patch_radius=2, tau_factor=2.3, out_dtype=None)

Local PCA-based denoising of diffusion datasets.

Parameters:
arr : 4D array

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

mask : 3D boolean array

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.

sigma : float or 3D array

Standard deviation of the noise estimated from the data.

pca_method : ‘eig’ or ‘svd’

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.

patch_radius : int, optional

The radius of the local patch to be taken around each voxel (in voxels). Default: 2 (denoise in blocks of 5x5x5 voxels).

tau_factor : float, optional

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

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

Default: 2.3, based on the results described in [Manjon13].

out_dtype : str or dtype, optional

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

Returns:
denoised_arr : 4D array

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

References

[Manjon13](1, 2) Manjon JV, Coupe P, Concha L, Buades A, Collins DL (2013) Diffusion Weighted Image Denoising Using Overcomplete Local PCA. PLoS ONE 8(9): e73021. https://doi.org/10.1371/journal.pone.0073021

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:
arr : 3D 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_radius : int

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

block_radius : int

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

rician : boolean

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

num_threads : int

Number of threads. If None (default) then all available threads will be used (all CPU cores).

Returns:
denoised_arr : ndarray

the denoised arr which has the same shape as arr.

References

[Descoteaux08]Descoteaux, Maxim and Wiest-Daessle`, Nicolas and Prima, Sylvain and Barillot, Christian and Deriche, Rachid Impact of Rician Adapted Non-Local Means Filtering on HARDI, MICCAI 2008

nlmeans_3d

dipy.denoise.nlmeans.nlmeans_3d()

Non-local means for denoising 3D images

Parameters:
arr : 3D ndarray

The array to be denoised

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

standard deviation of the noise estimated from the data

patch_radius : int

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

block_radius : int

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

rician : boolean

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

num_threads : int

Number of threads. If None (default) then all available threads will be used.

Returns:
denoised_arr : ndarray

the denoised arr which has the same shape as arr.

convolve

dipy.denoise.noise_estimate.convolve(input, weights, output=None, mode='reflect', cval=0.0, origin=0)

Multidimensional convolution.

The array is convolved with the given kernel.

Parameters:
input : array_like

The input array.

weights : array_like

Array of weights, same number of dimensions as input

output : array or dtype, optional

The array in which to place the output, or the dtype of the returned array. By default an array of the same dtype as input will be created.

mode : str or sequence, optional

The mode parameter determines how the input array is extended when the filter overlaps a border. By passing a sequence of modes with length equal to the number of dimensions of the input array, different modes can be specified along each axis. Default value is ‘reflect’. The valid values and their behavior is as follows:

‘reflect’ (d c b a | a b c d | d c b a)

The input is extended by reflecting about the edge of the last pixel.

‘constant’ (k k k k | a b c d | k k k k)

The input is extended by filling all values beyond the edge with the same constant value, defined by the cval parameter.

‘nearest’ (a a a a | a b c d | d d d d)

The input is extended by replicating the last pixel.

‘mirror’ (d c b | a b c d | c b a)

The input is extended by reflecting about the center of the last pixel.

‘wrap’ (a b c d | a b c d | a b c d)

The input is extended by wrapping around to the opposite edge.

cval : scalar, optional

Value to fill past edges of input if mode is ‘constant’. Default is 0.0

origin : int or sequence, optional

Controls the placement of the filter on the input array’s pixels. A value of 0 (the default) centers the filter over the pixel, with positive values shifting the filter to the left, and negative ones to the right. By passing a sequence of origins with length equal to the number of dimensions of the input array, different shifts can be specified along each axis.

Returns:
result : ndarray

The result of convolution of input with weights.

See also

correlate
Correlate an image with a kernel.

Notes

Each value in result is \(C_i = \sum_j{I_{i+k-j} W_j}\), where W is the weights kernel, j is the n-D spatial index over \(W\), I is the input and k is the coordinate of the center of W, specified by origin in the input parameters.

Examples

Perhaps the simplest case to understand is mode='constant', cval=0.0, because in this case borders (i.e. where the weights kernel, centered on any one value, extends beyond an edge of input.

>>> a = np.array([[1, 2, 0, 0],
...               [5, 3, 0, 4],
...               [0, 0, 0, 7],
...               [9, 3, 0, 0]])
>>> k = np.array([[1,1,1],[1,1,0],[1,0,0]])
>>> from scipy import ndimage
>>> ndimage.convolve(a, k, mode='constant', cval=0.0)
array([[11, 10,  7,  4],
       [10,  3, 11, 11],
       [15, 12, 14,  7],
       [12,  3,  7,  0]])

Setting cval=1.0 is equivalent to padding the outer edge of input with 1.0’s (and then extracting only the original region of the result).

>>> ndimage.convolve(a, k, mode='constant', cval=1.0)
array([[13, 11,  8,  7],
       [11,  3, 11, 14],
       [16, 12, 14, 10],
       [15,  6, 10,  5]])

With mode='reflect' (the default), outer values are reflected at the edge of input to fill in missing values.

>>> b = np.array([[2, 0, 0],
...               [1, 0, 0],
...               [0, 0, 0]])
>>> k = np.array([[0,1,0], [0,1,0], [0,1,0]])
>>> ndimage.convolve(b, k, mode='reflect')
array([[5, 0, 0],
       [3, 0, 0],
       [1, 0, 0]])

This includes diagonally at the corners.

>>> k = np.array([[1,0,0],[0,1,0],[0,0,1]])
>>> ndimage.convolve(b, k)
array([[4, 2, 0],
       [3, 2, 0],
       [1, 1, 0]])

With mode='nearest', the single nearest value in to an edge in input is repeated as many times as needed to match the overlapping weights.

>>> c = np.array([[2, 0, 1],
...               [1, 0, 0],
...               [0, 0, 0]])
>>> k = np.array([[0, 1, 0],
...               [0, 1, 0],
...               [0, 1, 0],
...               [0, 1, 0],
...               [0, 1, 0]])
>>> ndimage.convolve(c, k, mode='nearest')
array([[7, 0, 3],
       [5, 0, 2],
       [3, 0, 1]])

estimate_sigma

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

Standard deviation estimation from local patches

Parameters:
arr : 3D or 4D ndarray

The array to be estimated

disable_background_masking : bool, 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.

N : int, 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:
sigma : ndarray

standard deviation of the noise, one estimation per volume.

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:
data : ndarray

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

N : int

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.

alpha : float

Probabilistic estimation threshold for the gamma function.

l : int

number of initial estimates for sigma to try.

itermax : int

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

eps : float

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

return_mask : bool

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

Returns:
sigma : float

The estimated standard deviation of the gaussian noise.

mask : ndarray (optional)

A boolean mask indicating the voxels identified as pure noise.

References

[1]Koay CG, Ozarslan E and Pierpaoli C.

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

[2]Koay CG, Ozarslan E and Basser PJ.

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

nlmeans_block

dipy.denoise.non_local_means.nlmeans_block()

Non-Local Means Denoising Using Blockwise Averaging

Parameters:
image : 3D array of doubles

the input image, corrupted with rician noise

mask : 3D array of doubles

the input mask

patch_radius : int

similar patches in the non-local means are searched for locally, inside a cube of side 2*v+1 centered at each voxel of interest.

block_radius : int

the size of the block to be used (2*f+1)x(2*f+1)x(2*f+1) in the blockwise non-local means implementation (the Coupe’s proposal).

h : double

the estimated amount of rician noise in the input image: 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

rician : boolean

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

Returns:
fima: 3D double array

the denoised output which has the same shape as input image.

References

[1] 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
[2] 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

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:
arr : 3D or 4D ndarray

The array to be denoised

mask : 3D ndarray
sigma : float

standard deviation of the noise estimated from the data

patch_radius : int

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

block_radius : int

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

rician : boolean

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

Returns:
denoised_arr : ndarray

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