Module: restoration

Image restoration module.

References

[R266]

François Orieux, Jean-François Giovannelli, and Thomas Rodet, “Bayesian estimation of regularization and point spread function parameters for Wiener-Hunt deconvolution”, J. Opt. Soc. Am. A 27, 1593-1607 (2010)

http://www.opticsinfobase.org/josaa/abstract.cfm?URI=josaa-27-7-1593

[R267]Richardson, William Hadley, “Bayesian-Based Iterative Method of Image Restoration”. JOSA 62 (1): 55–59. doi:10.1364/JOSA.62.000055, 1972
[R268]B. R. Hunt “A matrix theory proof of the discrete convolution theorem”, IEEE Trans. on Audio and Electroacoustics, vol. au-19, no. 4, pp. 285-288, dec. 1971
skimage.restoration.denoise_bilateral(image) Denoise image using bilateral filter.
skimage.restoration.denoise_tv_bregman(...) Perform total-variation denoising using split-Bregman optimization.
skimage.restoration.denoise_tv_chambolle(im) Perform total-variation denoising on 2D and 3D images.
skimage.restoration.richardson_lucy(image, psf) Richardson-Lucy deconvolution.
skimage.restoration.unsupervised_wiener(...) Unsupervised Wiener-Hunt deconvolution
skimage.restoration.unwrap_phase(image[, ...]) From image, wrapped to lie in the interval [-pi, pi), recover the original, unwrapped image.
skimage.restoration.wiener(image, psf, balance) Wiener-Hunt deconvolution

denoise_bilateral

skimage.restoration.denoise_bilateral(image, win_size=5, sigma_range=None, sigma_spatial=1, bins=10000, mode='constant', cval=0)

Denoise image using bilateral filter.

This is an edge-preserving and noise reducing denoising filter. It averages pixels based on their spatial closeness and radiometric similarity.

Spatial closeness is measured by the gaussian function of the euclidian distance between two pixels and a certain standard deviation (sigma_spatial).

Radiometric similarity is measured by the gaussian function of the euclidian distance between two color values and a certain standard deviation (sigma_range).

Parameters:

image : ndarray

Input image.

win_size : int

Window size for filtering.

sigma_range : float

Standard deviation for grayvalue/color distance (radiometric similarity). A larger value results in averaging of pixels with larger radiometric differences. Note, that the image will be converted using the img_as_float function and thus the standard deviation is in respect to the range [0, 1].

sigma_spatial : float

Standard deviation for range distance. A larger value results in averaging of pixels with larger spatial differences.

bins : int

Number of discrete values for gaussian weights of color filtering. A larger value results in improved accuracy.

mode : string

How to handle values outside the image borders. See scipy.ndimage.map_coordinates for detail.

cval : string

Used in conjunction with mode ‘constant’, the value outside the image boundaries.

Returns:

denoised : ndarray

Denoised image.

References

[R281]http://users.soe.ucsc.edu/~manduchi/Papers/ICCV98.pdf

denoise_tv_bregman

skimage.restoration.denoise_tv_bregman(image, weight, max_iter=100, eps=0.001, isotropic=True)

Perform total-variation denoising using split-Bregman optimization.

Total-variation denoising (also know as total-variation regularization) tries to find an image with less total-variation under the constraint of being similar to the input image, which is controlled by the regularization parameter.

Parameters:

image : ndarray

Input data to be denoised (converted using img_as_float`).

weight : float, optional

Denoising weight. The smaller the weight, the more denoising (at the expense of less similarity to the input). The regularization parameter lambda is chosen as 2 * weight.

eps : float, optional

Relative difference of the value of the cost function that determines the stop criterion. The algorithm stops when:

SUM((u(n) - u(n-1))**2) < eps

max_iter : int, optional

Maximal number of iterations used for the optimization.

isotropic : boolean, optional

Switch between isotropic and anisotropic TV denoising.

Returns:

u : ndarray

Denoised image.

References

[R282]http://en.wikipedia.org/wiki/Total_variation_denoising
[R283]Tom Goldstein and Stanley Osher, “The Split Bregman Method For L1 Regularized Problems”, ftp://ftp.math.ucla.edu/pub/camreport/cam08-29.pdf
[R284]Pascal Getreuer, “Rudin–Osher–Fatemi Total Variation Denoising using Split Bregman” in Image Processing On Line on 2012–05–19, http://www.ipol.im/pub/art/2012/g-tvd/article_lr.pdf
[R285]http://www.math.ucsb.edu/~cgarcia/UGProjects/BregmanAlgorithms_JacquelineBush.pdf

denoise_tv_chambolle

skimage.restoration.denoise_tv_chambolle(im, weight=50, eps=0.0002, n_iter_max=200, multichannel=False)

Perform total-variation denoising on 2D and 3D images.

Parameters:

im : ndarray (2d or 3d) of ints, uints or floats

Input data to be denoised. im can be of any numeric type, but it is cast into an ndarray of floats for the computation of the denoised image.

weight : float, optional

Denoising weight. The greater weight, the more denoising (at the expense of fidelity to input).

eps : float, optional

Relative difference of the value of the cost function that determines the stop criterion. The algorithm stops when:

(E_(n-1) - E_n) < eps * E_0

n_iter_max : int, optional

Maximal number of iterations used for the optimization.

multichannel : bool, optional

Apply total-variation denoising separately for each channel. This option should be true for color images, otherwise the denoising is also applied in the 3rd dimension.

Returns:

out : ndarray

Denoised image.

Notes

Make sure to set the multichannel parameter appropriately for color images.

The principle of total variation denoising is explained in http://en.wikipedia.org/wiki/Total_variation_denoising

The principle of total variation denoising is to minimize the total variation of the image, which can be roughly described as the integral of the norm of the image gradient. Total variation denoising tends to produce “cartoon-like” images, that is, piecewise-constant images.

This code is an implementation of the algorithm of Rudin, Fatemi and Osher that was proposed by Chambolle in [R286].

References

[R286](1, 2) A. Chambolle, An algorithm for total variation minimization and applications, Journal of Mathematical Imaging and Vision, Springer, 2004, 20, 89-97.

Examples

2D example on Lena image:

>>> from skimage import color, data
>>> lena = color.rgb2gray(data.lena())[:50, :50]
>>> lena += 0.5 * lena.std() * np.random.randn(*lena.shape)
>>> denoised_lena = denoise_tv_chambolle(lena, weight=60)

3D example on synthetic data:

>>> x, y, z = np.ogrid[0:20, 0:20, 0:20]
>>> mask = (x - 22)**2 + (y - 20)**2 + (z - 17)**2 < 8**2
>>> mask = mask.astype(np.float)
>>> mask += 0.2*np.random.randn(*mask.shape)
>>> res = denoise_tv_chambolle(mask, weight=100)

richardson_lucy

skimage.restoration.richardson_lucy(image, psf, iterations=50, clip=True)

Richardson-Lucy deconvolution.

Parameters:

image : ndarray

Input degraded image

psf : ndarray

The point spread function

iterations : int

Number of iterations. This parameter play to role of regularisation.

clip : boolean, optional

True by default. If true, pixel value of the result above 1 or under -1 are thresholded for skimage pipeline compatibility.

Returns:

im_deconv : ndarray

The deconvolved image

References

[R288]http://en.wikipedia.org/wiki/Richardson%E2%80%93Lucy_deconvolution

Examples

>>> from skimage import color, data, restoration
>>> camera = color.rgb2gray(data.camera())
>>> from scipy.signal import convolve2d
>>> psf = np.ones((5, 5)) / 25
>>> camera = convolve2d(camera, psf, 'same')
>>> camera += 0.1 * camera.std() * np.random.standard_normal(camera.shape)
>>> deconvolved = restoration.richardson_lucy(camera, psf, 5)

unsupervised_wiener

skimage.restoration.unsupervised_wiener(image, psf, reg=None, user_params=None, is_real=True, clip=True)

Unsupervised Wiener-Hunt deconvolution

Return the deconvolution with a Wiener-Hunt approach, where the hyperparameters are automatically estimated. The algorithm is a stochastic iterative process (Gibbs sampler) described in the reference below. See also wiener function.

Parameters:

image : (M, N) ndarray

The input degraded image

psf : ndarray

The impulse response (input image’s space) or the transfer function (Fourier space). Both are accepted. The transfer function is recognize as being complex (np.iscomplexobj(psf)).

reg : ndarray, optional

The regularisation operator. The Laplacian by default. It can be an impulse response or a transfer function, as for the psf.

user_params : dict

dictionary of gibbs parameters. See below.

clip : boolean, optional

True by default. If true, pixel value of the result above 1 or under -1 are thresholded for skimage pipeline compatibility.

Returns:

x_postmean : (M, N) ndarray

The deconvolved image (the posterior mean).

chains : dict

The keys noise and prior contain the chain list of noise and prior precision respectively.

Other Parameters:
 

The keys of ``user_params`` are: :

threshold : float

The stopping criterion: the norm of the difference between to successive approximated solution (empirical mean of object samples, see Notes section). 1e-4 by default.

burnin : int

The number of sample to ignore to start computation of the mean. 100 by default.

min_iter : int

The minimum number of iterations. 30 by default.

max_iter : int

The maximum number of iterations if threshold is not satisfied. 150 by default.

callback : callable (None by default)

A user provided callable to which is passed, if the function exists, the current image sample for whatever purpose. The user can store the sample, or compute other moments than the mean. It has no influence on the algorithm execution and is only for inspection.

Notes

The estimated image is design as the posterior mean of a probability law (from a Bayesian analysis). The mean is defined as a sum over all the possible images weighted by their respective probability. Given the size of the problem, the exact sum is not tractable. This algorithm use of MCMC to draw image under the posterior law. The practical idea is to only draw high probable image since they have the biggest contribution to the mean. At the opposite, the lowest probable image are draw less often since their contribution are low. Finally the empirical mean of these samples give us an estimation of the mean, and an exact computation with an infinite sample set.

References

[R288]

François Orieux, Jean-François Giovannelli, and Thomas Rodet, “Bayesian estimation of regularization and point spread function parameters for Wiener-Hunt deconvolution”, J. Opt. Soc. Am. A 27, 1593-1607 (2010)

http://www.opticsinfobase.org/josaa/abstract.cfm?URI=josaa-27-7-1593

http://research.orieux.fr/files/papers/OGR-JOSA10.pdf

Examples

>>> from skimage import color, data, restoration
>>> lena = color.rgb2gray(data.lena())
>>> from scipy.signal import convolve2d
>>> psf = np.ones((5, 5)) / 25
>>> lena = convolve2d(lena, psf, 'same')
>>> lena += 0.1 * lena.std() * np.random.standard_normal(lena.shape)
>>> deconvolved_lena = restoration.unsupervised_wiener(lena, psf)

unwrap_phase

skimage.restoration.unwrap_phase(image, wrap_around=False)

From image, wrapped to lie in the interval [-pi, pi), recover the original, unwrapped image.

Parameters:

image : 1D, 2D or 3D ndarray of floats, optionally a masked array

The values should be in the range [-pi, pi). If a masked array is provided, the masked entries will not be changed, and their values will not be used to guide the unwrapping of neighboring, unmasked values. Masked 1D arrays are not allowed, and will raise a ValueError.

wrap_around : bool or sequence of bool

When an element of the sequence is True, the unwrapping process will regard the edges along the corresponding axis of the image to be connected and use this connectivity to guide the phase unwrapping process. If only a single boolean is given, it will apply to all axes. Wrap around is not supported for 1D arrays.

Returns:

image_unwrapped : array_like, float32

Unwrapped image of the same shape as the input. If the input image was a masked array, the mask will be preserved.

Raises:

ValueError :

If called with a masked 1D array or called with a 1D array and wrap_around=True.

References

[R289]Miguel Arevallilo Herraez, David R. Burton, Michael J. Lalor, and Munther A. Gdeisat, “Fast two-dimensional phase-unwrapping algorithm based on sorting by reliability following a noncontinuous path”, Journal Applied Optics, Vol. 41, No. 35 (2002) 7437,
[R290]Abdul-Rahman, H., Gdeisat, M., Burton, D., & Lalor, M., “Fast three-dimensional phase-unwrapping algorithm based on sorting by reliability following a non-continuous path. In W. Osten, C. Gorecki, & E. L. Novak (Eds.), Optical Metrology (2005) 32–40, International Society for Optics and Photonics.

Examples

>>> c0, c1 = np.ogrid[-1:1:128j, -1:1:128j]
>>> image = 12 * np.pi * np.exp(-(c0**2 + c1**2))
>>> image_wrapped = np.angle(np.exp(1j * image))
>>> image_unwrapped = unwrap_phase(image_wrapped)
>>> np.std(image_unwrapped - image) < 1e-6   # A constant offset is normal
True

wiener

skimage.restoration.wiener(image, psf, balance, reg=None, is_real=True, clip=True)

Wiener-Hunt deconvolution

Return the deconvolution with a Wiener-Hunt approach (i.e. with Fourier diagonalisation).

Parameters:

image : (M, N) ndarray

Input degraded image

psf : ndarray

Point Spread Function. This is assumed to be the impulse response (input image space) if the data-type is real, or the transfer function (Fourier space) if the data-type is complex. There is no constraints on the shape of the impulse response. The transfer function must be of shape (M, N) if is_real is True, (M, N // 2 + 1) otherwise (see np.fft.rfftn).

balance : float

The regularisation parameter value that tunes the balance between the data adequacy that improve frequency restoration and the prior adequacy that reduce frequency restoration (to avoid noise artifact).

reg : ndarray, optional

The regularisation operator. The Laplacian by default. It can be an impulse response or a transfer function, as for the psf. Shape constraint is the same than for the psf parameter.

is_real : boolean, optional

True by default. Specify if psf and reg are provided with hermitian hypothesis, that is only half of the frequency plane is provided (due to the redundancy of Fourier transform of real signal). It’s apply only if psf and/or reg are provided as transfer function. For the hermitian property see uft module or np.fft.rfftn.

clip : boolean, optional

True by default. If true, pixel value of the result above 1 or under -1 are thresholded for skimage pipeline compatibility.

Returns:

im_deconv : (M, N) ndarray

The deconvolved image

Notes

This function applies the Wiener filter to a noisy and degraded image by an impulse response (or PSF). If the data model is

y = Hx + n

where n is noise, H the PSF and x the unknown original image, the Wiener filter is

\hat x = F^\dag (|\Lambda_H|^2 + \lambda |\Lambda_D|^2) \Lambda_H^\dag F y

where F and F^\dag are the Fourier and inverse Fourier transfroms respectively, \Lambda_H the transfer function (or the Fourier transfrom of the PSF, see [Hunt] below) and \Lambda_D the filter to penalize the restored image frequencies (Laplacian by default, that is penalization of high frequency). The parameter \lambda tunes the balance between the data (that tends to increase high frequency, even those coming from noise), and the regularization.

These methods are then specific to a prior model. Consequently, the application or the true image nature must corresponds to the prior model. By default, the prior model (Laplacian) introduce image smoothness or pixel correlation. It can also be interpreted as high-frequency penalization to compensate the instability of the solution wrt. data (sometimes called noise amplification or “explosive” solution).

Finally, the use of Fourier space implies a circulant property of H, see [Hunt].

References

[R291]

François Orieux, Jean-François Giovannelli, and Thomas Rodet, “Bayesian estimation of regularization and point spread function parameters for Wiener-Hunt deconvolution”, J. Opt. Soc. Am. A 27, 1593-1607 (2010)

http://www.opticsinfobase.org/josaa/abstract.cfm?URI=josaa-27-7-1593

http://research.orieux.fr/files/papers/OGR-JOSA10.pdf

[R292]B. R. Hunt “A matrix theory proof of the discrete convolution theorem”, IEEE Trans. on Audio and Electroacoustics, vol. au-19, no. 4, pp. 285-288, dec. 1971

Examples

>>> from skimage import color, data, restoration
>>> lena = color.rgb2gray(data.lena())
>>> from scipy.signal import convolve2d
>>> psf = np.ones((5, 5)) / 25
>>> lena = convolve2d(lena, psf, 'same')
>>> lena += 0.1 * lena.std() * np.random.standard_normal(lena.shape)
>>> deconvolved_lena = restoration.wiener(lena, psf, 1100)