Denoising

This module supports these filter algorithms:

See also

See the file test_dnoise.py for an example.

Whiten

Function

pysptools.noise.dnoise.whiten(M)[source]

Whitens a HSI cube. Use the noise covariance matrix to decorrelate and rescale the noise in the data (noise whitening). Results in transformed data in which the noise has unit variance and no band-to-band correlations.

Parameters:Mnumpy array 2d matrix of HSI data (N x p).
Returns: numpy array
Whitened HSI data (N x p).
Reference:
Krizhevsky, Alex, Learning Multiple Layers of Features from Tiny Images, MSc thesis, University of Toronto, 2009. See Appendix A.

Class

class pysptools.noise.Whiten[source]

Whiten the cube.

apply(M)

Whitens a HSI cube. Use the noise covariance matrix to decorrelate and rescale the noise in the data (noise whitening). Results in transformed data in which the noise has unit variance and no band-to-band correlations.

Parameters:Mnumpy array A HSI cube (m x n x p).
Returns: numpy array
A whitened HSI cube (m x n x p).
get()[source]
Returns: numpy array
The whitened HSI cube (m x n x p).

Minimum Noise Fraction (MNF) Transformation

class pysptools.noise.MNF[source]

Transform a HSI cube.

apply(M)

A linear transformation that consists of a noise whitening step and one PCA rotation.

This process is designed to
  • determine the inherent dimensionality of image data,
  • segregate noise in the data,
  • allow efficient elimination and/or reduction of noise, and
  • reduce the computational requirements for subsequent processing.
Parameters:Mnumpy array A HSI cube (m x n x p).
Returns: numpy array
A MNF transformed cube (m x n x p).

References

C-I Change and Q Du, “Interference and Noise-Adjusted Principal Components Analysis,” IEEE TGRS, Vol 36, No 5, September 1999.

display_components(n_first=None, colorMap='jet', suffix=None)[source]

Display some bands.

Parameters:
  • n_firstint [default None] Display the first n components.
  • colorMapstring [default jet] A matplotlib color map.
  • suffixstring [default None] Suffix to add to the title.
get_components(n)[source]
Return: numpy array
Return the first n bands (maximum variance bands).
inverse_transform(X)

Inverse the PCA rotation step. The cube stay whitened. Usefull if you want to denoise noisy bands before the rotation.

X: numpy array
A transformed (MNF) cube (m x n x p).
Return: numpy array
A inverted cube (m x n x p).
plot_components(path, n_first=None, colorMap='jet', suffix=None)[source]

Plot some bands.

Parameters:
  • pathstring The path where to put the plot.
  • n_firstint [default None] Print the first n components.
  • colorMapstring [default jet] A matplotlib color map.
  • suffixstring [default None] Suffix to add to the title.

Savitzky Golay

class pysptools.noise.SavitzkyGolay[source]

Apply a Savitzky Golay low pass filter.

denoise_bands(M, window_size, order, derivative=None)

Apply the Savitzky Golay filter on each band.

Parameters:
  • Mnumpy array A HSI cube (m x n x p).
  • window_sizeint the length of the window. Must be an odd integer number.
  • orderint the order of the polynomial used in the filtering. Must be less then window_size - 1.
  • derivativestring [default None] direction of the derivative to compute, can be None, ‘col’, ‘row’, ‘both’.
Returns: numpy array
the smoothed signal (or it’s n-th derivative) (m x n x p).
Code source:
The scipy Cookbook, SavitzkyGolay section. This class is not under the copyright of this file.
denoise_spectra(M, window_size, order, deriv=0, rate=1)

Apply the Savitzky Golay filter on each spectrum. Smooth (and optionally differentiate) data with a Savitzky-Golay filter. The Savitzky-Golay filter removes high frequency noise from data. It has the advantage of preserving the original shape and features of the signal better than other types of filtering approaches, such as moving averages techniques.

Parameters:
  • Mnumpy array A HSI cube (m x n x p).
  • window_sizeint the length of the window. Must be an odd integer number.
  • orderint the order of the polynomial used in the filtering. Must be less then window_size - 1.
  • derivint [default 0] the order of the derivative to compute (default = 0 means only smoothing).
Returns: numpy array
the smoothed signal (or it’s n-th derivative) (m x n x p).
Code source:
The scipy Cookbook, SavitzkyGolay section. This class is not under the copyright of this file.
display_bands_sample(band_no, suffix=None)[source]

Display a filtered band.

Parameters:
  • band_noint or string The band index. If band_no == ‘all’, plot all the bands.
  • suffixstring [default None] Add a suffix to the title.
display_spectrum_sample(M, x, y, suffix=None)[source]

Display a spectrum sample with the original and the filtered signal.

Parameters:
  • M – ‘numpy array’ The original cube (m x n x p).
  • xint The x coordinate.
  • yint The y coordinate.
  • suffixstring [default None] Add a suffix to the title.
plot_bands_sample(path, band_no, suffix=None)[source]

Plot a filtered band.

Parameters:
  • pathstring The path where to put the plot.
  • band_noint or string The band index. If band_no == ‘all’, plot all the bands.
  • suffixstring [default None] Add a suffix to the file name.
plot_spectrum_sample(M, path, x, y, suffix=None)[source]

Plot a spectrum sample with the original and the filtered signal.

Parameters:
  • M – ‘numpy array’ The original cube (m x n x p).
  • pathstring The path where to put the plot.
  • xint The x coordinate.
  • yint The y coordinate.
  • suffixstring [default None] Add a suffix to the file name.