denoise#

Module: denoise.adaptive_soft_matching#

adaptive_soft_matching(ima, fimau, fimao, sigma)

Adaptive Soft Coefficient Matching

Module: denoise.denspeed#

add_padding_reflection(arr, padding)

correspond_indices(dim_size, padding)

determine_num_threads(num_threads)

Determine the effective number of threads to be used for OpenMP calls

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

Non-local means for denoising 3D images

remove_padding(arr, padding)

Module: denoise.enhancement_kernel#

EnhancementKernel

HemiSphere(*[, x, y, z, theta, phi, xyz, ...])

Points on the unit sphere.

Sphere(*[, x, y, z, theta, phi, xyz, faces, ...])

Points on the unit sphere.

ceil(x, /)

Return the ceiling of x as an Integral.

disperse_charges(hemi, iters, *[, const])

Models electrostatic repulsion on the unit sphere

get_sphere(*[, name])

provide triangulated spheres

gettempdir()

Returns tempfile.tempdir as str.

Module: denoise.gibbs#

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

Suppresses Gibbs ringing artefacts of images volumes.

Module: denoise.localpca#

dimensionality_problem_message(arr, ...)

Message about the number of samples being smaller than one less the dimensionality of the data to be denoised.

create_patch_radius_arr(arr, pr)

Create the patch radius array from the data to be denoised and the patch radius.

compute_patch_size(patch_radius)

Compute patch size from the patch radius.

compute_num_samples(patch_size)

Compute the number of samples as the dot product of the elements.

compute_suggested_patch_radius(arr, patch_size)

Compute the suggested patch radius.

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

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

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

Performs local PCA denoising.

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

Performs PCA-based denoising using the Marcenko-Pastur distribution.

Module: denoise.nlmeans#

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

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

Module: denoise.nlmeans_block#

firdn(image, h)

Applies the filter given by the convolution kernel 'h' columnwise to 'image', then subsamples by 2.

nlmeans_block(image, mask, patch_radius, ...)

Non-Local Means Denoising Using Blockwise Averaging.

upfir(image, h)

Upsamples the columns of the input image by 2, then applies the convolution kernel 'h' (again, columnwise).

Module: denoise.noise_estimate#

piesno(data, N, *[, alpha, step, itermax, ...])

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 blockwise averaging approach.

Module: denoise.patch2self#

vol_denoise(data_dict, b0_idx, dwi_idx, ...)

Denoise a single 3D volume using train and test phase.

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

Patch2Self Denoiser.

Module: denoise.pca_noise_estimate#

pca_noise_estimate(data, gtab[, ...])

PCA based local noise estimation.

warn()

Issue a warning, or maybe ignore it or raise an exception.

Module: denoise.shift_twist_convolution#

convolve(odfs_sh, kernel, sh_order_max[, ...])

Perform the shift-twist convolution with the ODF data and the lookup-table of the kernel.

convolve_sf(odfs_sf, kernel[, test_mode, ...])

Perform the shift-twist convolution with the ODF data and the lookup-table of the kernel.

deprecated_params(old_name, *[, new_name, ...])

Deprecate a renamed or removed function argument.

determine_num_threads(num_threads)

Determine the effective number of threads to be used for OpenMP calls

sf_to_sh(sf, sphere, *[, sh_order_max, ...])

Spherical function to spherical harmonics (SH).

sh_to_sf(sh, sphere, *[, sh_order_max, ...])

Spherical harmonics (SH) to spherical function (SF).

adaptive_soft_matching#

dipy.denoise.adaptive_soft_matching.adaptive_soft_matching(ima, fimau, fimao, sigma)[source]#

Adaptive Soft Coefficient Matching

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

Parameters:
imandarray

Original (unfiltered) image

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

add_padding_reflection#

dipy.denoise.denspeed.add_padding_reflection(arr, padding)#

correspond_indices#

dipy.denoise.denspeed.correspond_indices(dim_size, padding)#

determine_num_threads#

dipy.denoise.denspeed.determine_num_threads(num_threads)#

Determine the effective number of threads to be used for OpenMP calls

  • For num_threads = None, - if the OMP_NUM_THREADS environment variable is set, return that value - otherwise, return the maximum number of cores retrieved by openmp.opm_get_num_procs().

  • For num_threads > 0, return this value.

  • For num_threads < 0, return the maximal number of threads minus |num_threads + 1|. In particular num_threads = -1 will use as many threads as there are available cores on the machine.

  • For num_threads = 0 a ValueError is raised.

Parameters:
num_threadsint or None

Desired number of threads to be used.

nlmeans_3d#

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

Non-local means for denoising 3D images

Parameters:
arr3D ndarray

The array to be denoised

mask3D ndarray
sigmafloat 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.

remove_padding#

dipy.denoise.denspeed.remove_padding(arr, padding)#

EnhancementKernel#

class dipy.denoise.enhancement_kernel.EnhancementKernel#

Bases: object

Methods

evaluate_kernel(x, y, r, v)

Evaluate the kernel at position x relative to position y, with orientation r relative to orientation v.

get_lookup_table()

Return the computed look-up table.

get_orientations()

Return the orientations.

get_sphere()

Get the sphere corresponding with the orientations

evaluate_kernel(x, y, r, v)#

Evaluate the kernel at position x relative to position y, with orientation r relative to orientation v.

Parameters:
x1D ndarray

Position x

y1D ndarray

Position y

r1D ndarray

Orientation r

v1D ndarray

Orientation v

Returns:
kernel_valuedouble
get_lookup_table()#

Return the computed look-up table.

get_orientations()#

Return the orientations.

get_sphere()#

Get the sphere corresponding with the orientations

HemiSphere#

class dipy.denoise.enhancement_kernel.HemiSphere(*, x=None, y=None, z=None, theta=None, phi=None, xyz=None, faces=None, edges=None, tol=1e-05)[source]#

Bases: Sphere

Points on the unit sphere.

A HemiSphere is similar to a Sphere but it takes antipodal symmetry into account. Antipodal symmetry means that point v on a HemiSphere is the same as the point -v. Duplicate points are discarded when constructing a HemiSphere (including antipodal duplicates). edges and faces are remapped to the remaining points as closely as possible.

The HemiSphere can be constructed using one of three conventions:

HemiSphere(x, y, z)
HemiSphere(xyz=xyz)
HemiSphere(theta=theta, phi=phi)
Parameters:
x, y, z1-D array_like

Vertices as x-y-z coordinates.

theta, phi1-D array_like

Vertices as spherical coordinates. Theta and phi are the inclination and azimuth angles respectively.

xyz(N, 3) ndarray

Vertices as x-y-z coordinates.

faces(N, 3) ndarray

Indices into vertices that form triangular faces. If unspecified, the faces are computed using a Delaunay triangulation.

edges(N, 2) ndarray

Edges between vertices. If unspecified, the edges are derived from the faces.

tolfloat

Angle in degrees. Vertices that are less than tol degrees apart are treated as duplicates.

Attributes:
x
y
z

Methods

find_closest(xyz)

Find the index of the vertex in the Sphere closest to the input vector, taking into account antipodal symmetry

from_sphere(sphere, *[, tol])

Create instance from a Sphere

mirror()

Create a full Sphere from a HemiSphere

subdivide(*[, n])

Create a more subdivided HemiSphere

edges

faces

vertices

See also

Sphere
faces()[source]#
find_closest(xyz)[source]#

Find the index of the vertex in the Sphere closest to the input vector, taking into account antipodal symmetry

Parameters:
xyzarray-like, 3 elements

A unit vector

Returns:
idxint

The index into the Sphere.vertices array that gives the closest vertex (in angle).

classmethod from_sphere(sphere, *, tol=1e-05)[source]#

Create instance from a Sphere

mirror()[source]#

Create a full Sphere from a HemiSphere

subdivide(*, n=1)[source]#

Create a more subdivided HemiSphere

See Sphere.subdivide for full documentation.

Sphere#

class dipy.denoise.enhancement_kernel.Sphere(*, x=None, y=None, z=None, theta=None, phi=None, xyz=None, faces=None, edges=None)[source]#

Bases: object

Points on the unit sphere.

The sphere can be constructed using one of three conventions:

Sphere(x, y, z)
Sphere(xyz=xyz)
Sphere(theta=theta, phi=phi)
Parameters:
x, y, z1-D array_like

Vertices as x-y-z coordinates.

theta, phi1-D array_like

Vertices as spherical coordinates. Theta and phi are the inclination and azimuth angles respectively.

xyz(N, 3) ndarray

Vertices as x-y-z coordinates.

faces(N, 3) ndarray

Indices into vertices that form triangular faces. If unspecified, the faces are computed using a Delaunay triangulation.

edges(N, 2) ndarray

Edges between vertices. If unspecified, the edges are derived from the faces.

Attributes:
x
y
z

Methods

find_closest(xyz)

Find the index of the vertex in the Sphere closest to the input vector

subdivide(*[, n])

Subdivides each face of the sphere into four new faces.

edges

faces

vertices

edges()[source]#
faces()[source]#
find_closest(xyz)[source]#

Find the index of the vertex in the Sphere closest to the input vector

Parameters:
xyzarray-like, 3 elements

A unit vector

Returns:
idxint

The index into the Sphere.vertices array that gives the closest vertex (in angle).

subdivide(*, n=1)[source]#

Subdivides each face of the sphere into four new faces.

New vertices are created at a, b, and c. Then each face [x, y, z] is divided into faces [x, a, c], [y, a, b], [z, b, c], and [a, b, c].

      y
      /\
     /  \
   a/____\b
   /\    /\
  /  \  /  \
 /____\/____\
x      c     z
Parameters:
nint, optional

The number of subdivisions to perform.

Returns:
new_sphereSphere

The subdivided sphere.

vertices()[source]#
property x#
property y#
property z#

ceil#

dipy.denoise.enhancement_kernel.ceil(x, /)#

Return the ceiling of x as an Integral.

This is the smallest integer >= x.

disperse_charges#

dipy.denoise.enhancement_kernel.disperse_charges(hemi, iters, *, const=0.2)[source]#

Models electrostatic repulsion on the unit sphere

Places charges on a sphere and simulates the repulsive forces felt by each one. Allows the charges to move for some number of iterations and returns their final location as well as the total potential of the system at each step.

Parameters:
hemiHemiSphere

Points on a unit sphere.

itersint

Number of iterations to run.

constfloat

Using a smaller const could provide a more accurate result, but will need more iterations to converge.

Returns:
hemiHemiSphere

Distributed points on a unit sphere.

potentialndarray

The electrostatic potential at each iteration. This can be useful to check if the repulsion converged to a minimum.

Notes

This function is meant to be used with diffusion imaging so antipodal symmetry is assumed. Therefore, each charge must not only be unique, but if there is a charge at +x, there cannot be a charge at -x. These are treated as the same location and because the distance between the two charges will be zero, the result will be unstable.

get_sphere#

dipy.denoise.enhancement_kernel.get_sphere(*, name='symmetric362')[source]#

provide triangulated spheres

Parameters:
namestr

which sphere - one of: * ‘symmetric362’ * ‘symmetric642’ * ‘symmetric724’ * ‘repulsion724’ * ‘repulsion100’ * ‘repulsion200’

Returns:
spherea dipy.core.sphere.Sphere class instance

Examples

>>> import numpy as np
>>> from dipy.data import get_sphere
>>> sphere = get_sphere(name="symmetric362")
>>> verts, faces = sphere.vertices, sphere.faces
>>> verts.shape == (362, 3)
True
>>> faces.shape == (720, 3)
True
>>> verts, faces = get_sphere(name="not a sphere name") 
Traceback (most recent call last):
    ...
DataError: No sphere called "not a sphere name"

gettempdir#

dipy.denoise.enhancement_kernel.gettempdir()[source]#

Returns tempfile.tempdir as str.

gibbs_removal#

dipy.denoise.gibbs.gibbs_removal(vol, *, slice_axis=2, n_points=3, inplace=True, num_processes=1)[source]#

Suppresses Gibbs ringing artefacts of images volumes.

See [2] and [3] for further details about the method.

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.

n_pointsint, optional

Number of neighbour points to access local TV (see note).

inplacebool, optional

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

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

dimensionality_problem_message#

dipy.denoise.localpca.dimensionality_problem_message(arr, num_samples, spr)[source]#

Message about the number of samples being smaller than one less the dimensionality of the data to be denoised.

Parameters:
arrndarray

Data to be denoised.

num_samplesint

Number of samples.

sprint

Suggested patch radius.

Returns:
str

Message.

create_patch_radius_arr#

dipy.denoise.localpca.create_patch_radius_arr(arr, pr)[source]#

Create the patch radius array from the data to be denoised and the patch radius.

Parameters:
arrndarray

Data to be denoised.

print or ndarray

Patch radius.

Returns:
patch_radiusndarray

Number of samples.

compute_patch_size#

dipy.denoise.localpca.compute_patch_size(patch_radius)[source]#

Compute patch size from the patch radius.

it is twice the radius plus one.

Parameters:
patch_radiusndarray

Patch radius.

Returns:
int

Number of samples.

compute_num_samples#

dipy.denoise.localpca.compute_num_samples(patch_size)[source]#

Compute the number of samples as the dot product of the elements.

Parameters:
patch_sizendarray

Patch size.

Returns:
int

Number of samples.

compute_suggested_patch_radius#

dipy.denoise.localpca.compute_suggested_patch_radius(arr, patch_size)[source]#

Compute the suggested patch radius.

Parameters:
arrndarray

Data to be denoised.

patch_sizendarray

Patch size.

Returns:
int

Suggested patch radius.

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)[source]#

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 [4], [5].

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\) [6]), 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)[source]#

Performs local PCA denoising.

The method follows the algorithm in Manjón et al.[6].

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 Manjón et al.[6].

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 [5]. Default: 2.3 according to Manjón et al.[6].

return_sigmabool, optional

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

correct_biasbool, optional

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

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)[source]#

Performs PCA-based denoising using the Marcenko-Pastur distribution.

See [5] for further details about the method.

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 [4].

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)[source]#

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

See :footcite:p:Descoteaux2008a` for further details about the method.

Parameters:
arr3D or 4D ndarray

The array to be denoised

mask3D ndarray
sigmafloat 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

firdn#

dipy.denoise.nlmeans_block.firdn(image, h)#

Applies the filter given by the convolution kernel ‘h’ columnwise to ‘image’, then subsamples by 2. This is a special case of the matlab’s ‘upfirdn’ function, ported to python. Returns the filtered image.

Parameters:
image: 2D array of doubles

the input image to be filtered

h: double array

the convolution kernel

nlmeans_block#

dipy.denoise.nlmeans_block.nlmeans_block(image, mask, patch_radius, block_radius, h, rician)#

Non-Local Means Denoising Using Blockwise Averaging.

See [7] and [1] for further details about the method.

Parameters:
image3D array of doubles

the input image, corrupted with rician noise

mask3D array of doubles

the input mask

patch_radiusint

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_radiusint

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

hdouble

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

ricianboolean

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

upfir#

dipy.denoise.nlmeans_block.upfir(image, h)#

Upsamples the columns of the input image by 2, then applies the convolution kernel ‘h’ (again, columnwise). This is a special case of the matlab’s ‘upfirdn’ function, ported to python. Returns the filtered image.

Parameters:
image: 2D array of doubles

the input image to be filtered

h: double array

the convolution kernel

piesno#

dipy.denoise.noise_estimate.piesno(data, N, *, alpha=0.01, step=100, itermax=100, eps=1e-05, return_mask=False)[source]#

Probabilistic Identification and Estimation of Noise (PIESNO).

See [8] and [9] for further details about the method.

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.

stepint

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

estimate_sigma#

dipy.denoise.noise_estimate.estimate_sigma(arr, *, disable_background_masking=False, N=0)[source]#

Standard deviation estimation from local patches

Parameters:
arr3D or 4D ndarray

The array to be estimated

disable_background_maskingbool, optional

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, optional

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 [10] 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 from Koay and Basser[10] (see equation 18 in [10]) with theta = 0. Since this function was introduced in [7] 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

non_local_means#

dipy.denoise.non_local_means.non_local_means(arr, sigma, *, mask=None, patch_radius=1, block_radius=5, rician=True)[source]#

Non-local means for denoising 3D and 4D images, using blockwise averaging approach.

See [7] and [1] for further details about the method.

Parameters:
arr3D or 4D ndarray

The array to be denoised

mask3D ndarray

Mask on data where the non-local means will be applied.

sigmafloat or ndarray

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

vol_denoise#

dipy.denoise.patch2self.vol_denoise(data_dict, b0_idx, dwi_idx, model, alpha, b0_denoising, verbose, tmp_dir)[source]#

Denoise a single 3D volume using train and test phase.

Parameters:
data_dictdict
Dictionary containing the following:
data_namestr

The name of the memmap file containing the memmaped data.

data_dtypedtype

The dtype of the data.

data_shapetuple

The shape of the data.

data_b0sndarray

Array of all 3D patches flattened out to be 2D for b0 volumes.

data_dwindarray

Array of all 3D patches flattened out to be 2D for dwi volumes.

b0_idxndarray

The indices of the b0 volumes.

dwi_idxndarray

The indices of the dwi volumes.

modelsklearn.base.RegressorMixin

This is the model that is initialized from the _fit_denoising_model function.

alphafloat

Regularization parameter only for ridge and lasso regression models.

b0_denoisingbool

Skips denoising b0 volumes if set to False.

verbosebool

Show progress of Patch2Self and time taken.

tmp_dirstr

The directory to save the temporary files.

Returns:
denoised_arr.namestr

The name of the memmap file containing the denoised array.

denoised_arr.dtypedtype

The dtype of the denoised array.

denoised_arr.shapetuple

The shape of the denoised array.

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, tmp_dir=None, version=3)[source]#

Patch2Self Denoiser.

See [11] for further details about the method. See [12] for further details about the new method.

Parameters:
datandarray

The 4D noisy DWI data to be denoised.

bvalsarray of shape (N,)

Array of the bvals from the DWI acquisition

patch_radiusint or array of shape (3,), 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 sklearn.base.RegressorMixin

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.

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.

tmp_dirstr, optional

The directory to save the temporary files. If None, the temporary files are saved in the system’s default temporary directory. Default: None.

versionint, optional

Version 1 or 3 of Patch2Self to use. Default: 3

Returns:
denoised arrayndarray

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

References

pca_noise_estimate#

dipy.denoise.pca_noise_estimate.pca_noise_estimate(data, gtab, patch_radius=1, correct_bias=True, smooth=2, images_as_samples=False)#

PCA based local noise estimation.

Parameters:
data4D array

the input dMRI data. The first 3 dimension must have size >= 2 * patch_radius + 1 or size = 1.

gtabgradient table object

gradient information for the data gives us the bvals and bvecs of diffusion data, which is needed here to select between the noise estimation methods.

patch_radiusint

The radius of the local patch to be taken around each voxel (in voxels). Default: 1 (estimate noise in blocks of 3x3x3 voxels).

correct_biasbool

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

smoothint

Radius of a Gaussian smoothing filter to apply to the noise estimate before returning. Default: 2.

images_as_samplesbool, optional

Whether to use images as rows (samples) for PCA (algorithm in [6]) or to use images as columns (features).

Returns:
sigma_corr: 3D array

The local noise standard deviation estimate.

Notes

In [6], images are used as samples, so voxels are features, therefore eigenvectors are image-shaped. However, Manjón et al.[6] is not clear on how to use these eigenvectors to determine the noise level, so here eigenvalues (variance over samples explained by eigenvectors) are used to scale the eigenvectors. Use images_as_samples=True to use this algorithm. Alternatively, voxels can be used as samples using images_as_samples=False. This is not the canonical algorithm of Manjón et al.[6].

References

warn#

warn(/, message, category=None, stacklevel=1, source=None, *,
warn     skip_file_prefixes=<unrepresentable>)

Issue a warning, or maybe ignore it or raise an exception.

message

Text of the warning message.

category

The Warning category subclass. Defaults to UserWarning.

stacklevel

How far up the call stack to make this warning appear. A value of 2 for example attributes the warning to the caller of the code calling warn().

source

If supplied, the destroyed object which emitted a ResourceWarning

skip_file_prefixes

An optional tuple of module filename prefixes indicating frames to skip during stacklevel computations for stack frame attribution.

convolve#

dipy.denoise.shift_twist_convolution.convolve(odfs_sh, kernel, sh_order_max, test_mode=False, num_threads=None, normalize=True)#

Perform the shift-twist convolution with the ODF data and the lookup-table of the kernel.

See [13], [14], [15] and [16] for further details about the method.

Parameters:
odfsarray of double

The ODF data in spherical harmonics format

kernelarray of double

The 5D lookup table

sh_order_maxinteger

Maximal spherical harmonics order (l)

test_modeboolean

Reduced convolution in one direction only for testing

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.

normalizeboolean

Apply max-normalization to the output such that its value range matches the input ODF data.

Returns:
outputarray of double

The ODF data after convolution enhancement in spherical harmonics format

References

convolve_sf#

dipy.denoise.shift_twist_convolution.convolve_sf(odfs_sf, kernel, test_mode=False, num_threads=None, normalize=True)#

Perform the shift-twist convolution with the ODF data and the lookup-table of the kernel.

Parameters:
odfsarray of double

The ODF data sampled on a sphere

kernelarray of double

The 5D lookup table

test_modeboolean

Reduced convolution in one direction only for testing

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.

normalizeboolean

Apply max-normalization to the output such that its value range matches the input ODF data.

Returns:
outputarray of double

The ODF data after convolution enhancement, sampled on a sphere

deprecated_params#

dipy.denoise.shift_twist_convolution.deprecated_params(old_name, *, new_name=None, since='', until='', version_comparator=<function cmp_pkg_version>, arg_in_kwargs=False, warn_class=<class 'dipy.utils.deprecator.ArgsDeprecationWarning'>, error_class=<class 'dipy.utils.deprecator.ExpiredDeprecationError'>, alternative='')[source]#

Deprecate a renamed or removed function argument.

The decorator assumes that the argument with the old_name was removed from the function signature and the new_name replaced it at the same position in the signature. If the old_name argument is given when calling the decorated function the decorator will catch it and issue a deprecation warning and pass it on as new_name argument.

Parameters:
old_namestr or list/tuple thereof

The old name of the argument.

new_namestr or list/tuple thereof or None, optional

The new name of the argument. Set this to None to remove the argument old_name instead of renaming it.

sincestr or number or list/tuple thereof, optional

The release at which the old argument became deprecated.

untilstr or number or list/tuple thereof, optional

Last released version at which this function will still raise a deprecation warning. Versions higher than this will raise an error.

version_comparatorcallable

Callable accepting string as argument, and return 1 if string represents a higher version than encoded in the version_comparator, 0 if the version is equal, and -1 if the version is lower. For example, the version_comparator may compare the input version string to the current package version string.

arg_in_kwargsbool or list/tuple thereof, optional

If the argument is not a named argument (for example it was meant to be consumed by **kwargs) set this to True. Otherwise the decorator will throw an Exception if the new_name cannot be found in the signature of the decorated function. Default is False.

warn_classwarning, optional

Warning to be issued.

error_classException, optional

Error to be issued

alternativestr, optional

An alternative function or class name that the user may use in place of the deprecated object if new_name is None. The deprecation warning will tell the user about this alternative if provided.

Raises:
TypeError

If the new argument name cannot be found in the function signature and arg_in_kwargs was False or if it is used to deprecate the name of the *args-, **kwargs-like arguments. At runtime such an Error is raised if both the new_name and old_name were specified when calling the function and “relax=False”.

Notes

This function is based on the Astropy (major modification). astropy/astropy. See COPYING file distributed along with the astropy package for the copyright and license terms.

Examples

The deprecation warnings are not shown in the following examples. To deprecate a positional or keyword argument:: >>> from dipy.utils.deprecator import deprecated_params >>> @deprecated_params(‘sig’, new_name=’sigma’, since=’0.3’) … def test(sigma): … return sigma >>> test(2) 2 >>> test(sigma=2) 2 >>> test(sig=2) # doctest: +SKIP 2

It is also possible to replace multiple arguments. The old_name, new_name and since have to be tuple or list and contain the same number of entries:: >>> @deprecated_params([‘a’, ‘b’], new_name=[‘alpha’, ‘beta’], … since=[‘0.2’, 0.4]) … def test(alpha, beta): … return alpha, beta >>> test(a=2, b=3) # doctest: +SKIP (2, 3)

determine_num_threads#

dipy.denoise.shift_twist_convolution.determine_num_threads(num_threads)#

Determine the effective number of threads to be used for OpenMP calls

  • For num_threads = None, - if the OMP_NUM_THREADS environment variable is set, return that value - otherwise, return the maximum number of cores retrieved by openmp.opm_get_num_procs().

  • For num_threads > 0, return this value.

  • For num_threads < 0, return the maximal number of threads minus |num_threads + 1|. In particular num_threads = -1 will use as many threads as there are available cores on the machine.

  • For num_threads = 0 a ValueError is raised.

Parameters:
num_threadsint or None

Desired number of threads to be used.

sf_to_sh#

dipy.denoise.shift_twist_convolution.sf_to_sh(sf, sphere, *, sh_order_max=4, basis_type=None, full_basis=False, legacy=True, smooth=0.0)[source]#

Spherical function to spherical harmonics (SH).

Parameters:
sfndarray

Values of a function on the given sphere.

sphereSphere

The points on which the sf is defined.

sh_order_maxint, optional

Maximum SH order (l) in the SH fit. For sh_order_max, there will be (sh_order_max + 1) * (sh_order_max + 2) / 2 SH coefficients for a symmetric basis and (sh_order_max + 1) * (sh_order_max + 1) coefficients for a full SH basis.

basis_type{None, ‘tournier07’, ‘descoteaux07’}, optional

None for the default DIPY basis, tournier07 for the Tournier 2007 [17], [18] basis, descoteaux07 for the Descoteaux 2007 [19] basis, (None defaults to descoteaux07).

full_basis: bool, optional

True for using a SH basis containing even and odd order SH functions. False for using a SH basis consisting only of even order SH functions.

legacy: bool, optional

True to use a legacy basis definition for backward compatibility with previous tournier07 and descoteaux07 implementations.

smoothfloat, optional

Lambda-regularization in the SH fit.

Returns:
shndarray

SH coefficients representing the input function.

References

sh_to_sf#

dipy.denoise.shift_twist_convolution.sh_to_sf(sh, sphere, *, sh_order_max=4, basis_type=None, full_basis=False, legacy=True)[source]#

Spherical harmonics (SH) to spherical function (SF).

Parameters:
shndarray

SH coefficients representing a spherical function.

sphereSphere

The points on which to sample the spherical function.

sh_order_maxint, optional

Maximum SH order (l) in the SH fit. For sh_order_max, there will be (sh_order_max + 1) * (sh_order_max + 2) / 2 SH coefficients for a symmetric basis and (sh_order_max + 1) * (sh_order_max + 1) coefficients for a full SH basis.

basis_type{None, ‘tournier07’, ‘descoteaux07’}, optional

None for the default DIPY basis, tournier07 for the Tournier 2007 [17], [18] basis, descoteaux07 for the Descoteaux 2007 [19] basis, (None defaults to descoteaux07).

full_basis: bool, optional

True to use a SH basis containing even and odd order SH functions. Else, use a SH basis consisting only of even order SH functions.

legacy: bool, optional

True to use a legacy basis definition for backward compatibility with previous tournier07 and descoteaux07 implementations.

Returns:
sfndarray

Spherical function values on the sphere.

References