symjax.tensor.signal

Implementation of various signal processing related techniques such as time-frequency representations convolution/correlation/pooling operations, as well as various apodization windows and filter-banks creations.

Apodization Windows

blackman(*args, **kwargs) Return the Blackman window.
bartlett(*args, **kwargs) Return the Bartlett window.
hamming(*args, **kwargs) Return the Hamming window.
hanning(*args, **kwargs) Return the Hanning window.
kaiser(*args, **kwargs) Return the Kaiser window.
tukey(M[, alpha]) Return a Tukey window, also known as a tapered cosine window.

Time-Frequency Representations

mfcc(signal, window, hop, n_filter, …[, …]) https://librosa.github.io/librosa/_modules/librosa/feature/spectral.html#mfcc
stft(signal, window, hop[, apod, nfft, mode]) Compute the Shoft-Time-Fourier-Transform of a signal given the window length, hop and additional parameters.
dct(signal[, axes]) https://dsp.stackexchange.com/questions/2807/fast-cosine-transform-via-fft
wvd(signal, window, hop, L[, apod, mode])
hilbert_transform(signal) the time should be the last dimension return the analytical signal

Filters

fourier_complex_morlet(bandwidths, centers, N) Complex Morlet wavelet in Fourier
complex_morlet(bandwidths, centers[, time]) Complex Morlet wavelet
sinc_bandpass(time, f0, f1) ensure that f0<f1 and f0>0, f1<1 whenever time is …, -1, 0, 1, …
mel_filterbank(length, n_filter, low, high, …)
hat_1d

Convolution/Correlation/Pooling

convolve(in1, in2[, mode, method, precision]) Convolve two N-dimensional arrays.
batch_convolve(input, filter[, strides, …]) General n-dimensional batch convolution with dilations.
convolve2d(in1, in2[, mode, boundary, …]) Convolve two 2-dimensional arrays.
correlate(in1, in2[, mode, method, precision]) Cross-correlate two N-dimensional arrays.
correlate2d(in1, in2[, mode, boundary, …]) Cross-correlate two 2-dimensional arrays.
batch_pool

Utilities

extract_signal_patches
extract_image_patches(image, window_shape[, …]) extract patches from an input tensor

Detailed Descriptions

symjax.tensor.signal.blackman(*args, **kwargs)

Return the Blackman window.

LAX-backend implementation of blackman(). Original docstring below.

The Blackman window is a taper formed by using the first three terms of a summation of cosines. It was designed to have close to the minimal leakage possible. It is close to optimal, only slightly worse than a Kaiser window.

M : int
Number of points in the output window. If zero or less, an empty array is returned.
out : ndarray
The window, with the maximum value normalized to one (the value one appears only if the number of samples is odd).

bartlett, hamming, hanning, kaiser

The Blackman window is defined as

\[w(n) = 0.42 - 0.5 \cos(2\pi n/M) + 0.08 \cos(4\pi n/M)\]

Most references to the Blackman window come from the signal processing literature, where it is used as one of many windowing functions for smoothing values. It is also known as an apodization (which means “removing the foot”, i.e. smoothing discontinuities at the beginning and end of the sampled signal) or tapering function. It is known as a “near optimal” tapering function, almost as good (by some measures) as the kaiser window.

Blackman, R.B. and Tukey, J.W., (1958) The measurement of power spectra, Dover Publications, New York.

Oppenheim, A.V., and R.W. Schafer. Discrete-Time Signal Processing. Upper Saddle River, NJ: Prentice-Hall, 1999, pp. 468-471.

>>> import matplotlib.pyplot as plt
>>> np.blackman(12)
array([-1.38777878e-17,   3.26064346e-02,   1.59903635e-01, # may vary
        4.14397981e-01,   7.36045180e-01,   9.67046769e-01,
        9.67046769e-01,   7.36045180e-01,   4.14397981e-01,
        1.59903635e-01,   3.26064346e-02,  -1.38777878e-17])

Plot the window and the frequency response:

>>> from numpy.fft import fft, fftshift
>>> window = np.blackman(51)
>>> plt.plot(window)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Blackman window")
Text(0.5, 1.0, 'Blackman window')
>>> plt.ylabel("Amplitude")
Text(0, 0.5, 'Amplitude')
>>> plt.xlabel("Sample")
Text(0.5, 0, 'Sample')
>>> plt.show()
>>> plt.figure()
<Figure size 640x480 with 0 Axes>
>>> A = fft(window, 2048) / 25.5
>>> mag = np.abs(fftshift(A))
>>> freq = np.linspace(-0.5, 0.5, len(A))
>>> with np.errstate(divide='ignore', invalid='ignore'):
...     response = 20 * np.log10(mag)
...
>>> response = np.clip(response, -100, 100)
>>> plt.plot(freq, response)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Frequency response of Blackman window")
Text(0.5, 1.0, 'Frequency response of Blackman window')
>>> plt.ylabel("Magnitude [dB]")
Text(0, 0.5, 'Magnitude [dB]')
>>> plt.xlabel("Normalized frequency [cycles per sample]")
Text(0.5, 0, 'Normalized frequency [cycles per sample]')
>>> _ = plt.axis('tight')
>>> plt.show()
symjax.tensor.signal.bartlett(*args, **kwargs)

Return the Bartlett window.

LAX-backend implementation of bartlett(). Original docstring below.

The Bartlett window is very similar to a triangular window, except that the end points are at zero. It is often used in signal processing for tapering a signal, without generating too much ripple in the frequency domain.

M : int
Number of points in the output window. If zero or less, an empty array is returned.
out : array
The triangular window, with the maximum value normalized to one (the value one appears only if the number of samples is odd), with the first and last samples equal to zero.

blackman, hamming, hanning, kaiser

The Bartlett window is defined as

\[w(n) = \frac{2}{M-1} \left( \frac{M-1}{2} - \left|n - \frac{M-1}{2}\right| \right)\]

Most references to the Bartlett window come from the signal processing literature, where it is used as one of many windowing functions for smoothing values. Note that convolution with this window produces linear interpolation. It is also known as an apodization (which means”removing the foot”, i.e. smoothing discontinuities at the beginning and end of the sampled signal) or tapering function. The fourier transform of the Bartlett is the product of two sinc functions. Note the excellent discussion in Kanasewich.

[1]M.S. Bartlett, “Periodogram Analysis and Continuous Spectra”, Biometrika 37, 1-16, 1950.
[2]E.R. Kanasewich, “Time Sequence Analysis in Geophysics”, The University of Alberta Press, 1975, pp. 109-110.
[3]A.V. Oppenheim and R.W. Schafer, “Discrete-Time Signal Processing”, Prentice-Hall, 1999, pp. 468-471.
[4]Wikipedia, “Window function”, https://en.wikipedia.org/wiki/Window_function
[5]W.H. Press, B.P. Flannery, S.A. Teukolsky, and W.T. Vetterling, “Numerical Recipes”, Cambridge University Press, 1986, page 429.
>>> import matplotlib.pyplot as plt
>>> np.bartlett(12)
array([ 0.        ,  0.18181818,  0.36363636,  0.54545455,  0.72727273, # may vary
        0.90909091,  0.90909091,  0.72727273,  0.54545455,  0.36363636,
        0.18181818,  0.        ])

Plot the window and its frequency response (requires SciPy and matplotlib):

>>> from numpy.fft import fft, fftshift
>>> window = np.bartlett(51)
>>> plt.plot(window)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Bartlett window")
Text(0.5, 1.0, 'Bartlett window')
>>> plt.ylabel("Amplitude")
Text(0, 0.5, 'Amplitude')
>>> plt.xlabel("Sample")
Text(0.5, 0, 'Sample')
>>> plt.show()
>>> plt.figure()
<Figure size 640x480 with 0 Axes>
>>> A = fft(window, 2048) / 25.5
>>> mag = np.abs(fftshift(A))
>>> freq = np.linspace(-0.5, 0.5, len(A))
>>> with np.errstate(divide='ignore', invalid='ignore'):
...     response = 20 * np.log10(mag)
...
>>> response = np.clip(response, -100, 100)
>>> plt.plot(freq, response)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Frequency response of Bartlett window")
Text(0.5, 1.0, 'Frequency response of Bartlett window')
>>> plt.ylabel("Magnitude [dB]")
Text(0, 0.5, 'Magnitude [dB]')
>>> plt.xlabel("Normalized frequency [cycles per sample]")
Text(0.5, 0, 'Normalized frequency [cycles per sample]')
>>> _ = plt.axis('tight')
>>> plt.show()
symjax.tensor.signal.hamming(*args, **kwargs)

Return the Hamming window.

LAX-backend implementation of hamming(). Original docstring below.

The Hamming window is a taper formed by using a weighted cosine.

M : int
Number of points in the output window. If zero or less, an empty array is returned.
out : ndarray
The window, with the maximum value normalized to one (the value one appears only if the number of samples is odd).

bartlett, blackman, hanning, kaiser

The Hamming window is defined as

\[w(n) = 0.54 - 0.46cos\left(\frac{2\pi{n}}{M-1}\right) \qquad 0 \leq n \leq M-1\]

The Hamming was named for R. W. Hamming, an associate of J. W. Tukey and is described in Blackman and Tukey. It was recommended for smoothing the truncated autocovariance function in the time domain. Most references to the Hamming window come from the signal processing literature, where it is used as one of many windowing functions for smoothing values. It is also known as an apodization (which means “removing the foot”, i.e. smoothing discontinuities at the beginning and end of the sampled signal) or tapering function.

[1]Blackman, R.B. and Tukey, J.W., (1958) The measurement of power spectra, Dover Publications, New York.
[2]E.R. Kanasewich, “Time Sequence Analysis in Geophysics”, The University of Alberta Press, 1975, pp. 109-110.
[3]Wikipedia, “Window function”, https://en.wikipedia.org/wiki/Window_function
[4]W.H. Press, B.P. Flannery, S.A. Teukolsky, and W.T. Vetterling, “Numerical Recipes”, Cambridge University Press, 1986, page 425.
>>> np.hamming(12)
array([ 0.08      ,  0.15302337,  0.34890909,  0.60546483,  0.84123594, # may vary
        0.98136677,  0.98136677,  0.84123594,  0.60546483,  0.34890909,
        0.15302337,  0.08      ])

Plot the window and the frequency response:

>>> import matplotlib.pyplot as plt
>>> from numpy.fft import fft, fftshift
>>> window = np.hamming(51)
>>> plt.plot(window)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Hamming window")
Text(0.5, 1.0, 'Hamming window')
>>> plt.ylabel("Amplitude")
Text(0, 0.5, 'Amplitude')
>>> plt.xlabel("Sample")
Text(0.5, 0, 'Sample')
>>> plt.show()
>>> plt.figure()
<Figure size 640x480 with 0 Axes>
>>> A = fft(window, 2048) / 25.5
>>> mag = np.abs(fftshift(A))
>>> freq = np.linspace(-0.5, 0.5, len(A))
>>> response = 20 * np.log10(mag)
>>> response = np.clip(response, -100, 100)
>>> plt.plot(freq, response)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Frequency response of Hamming window")
Text(0.5, 1.0, 'Frequency response of Hamming window')
>>> plt.ylabel("Magnitude [dB]")
Text(0, 0.5, 'Magnitude [dB]')
>>> plt.xlabel("Normalized frequency [cycles per sample]")
Text(0.5, 0, 'Normalized frequency [cycles per sample]')
>>> plt.axis('tight')
...
>>> plt.show()
symjax.tensor.signal.hanning(*args, **kwargs)

Return the Hanning window.

LAX-backend implementation of hanning(). Original docstring below.

The Hanning window is a taper formed by using a weighted cosine.

M : int
Number of points in the output window. If zero or less, an empty array is returned.
out : ndarray, shape(M,)
The window, with the maximum value normalized to one (the value one appears only if M is odd).

bartlett, blackman, hamming, kaiser

The Hanning window is defined as

\[w(n) = 0.5 - 0.5cos\left(\frac{2\pi{n}}{M-1}\right) \qquad 0 \leq n \leq M-1\]

The Hanning was named for Julius von Hann, an Austrian meteorologist. It is also known as the Cosine Bell. Some authors prefer that it be called a Hann window, to help avoid confusion with the very similar Hamming window.

Most references to the Hanning window come from the signal processing literature, where it is used as one of many windowing functions for smoothing values. It is also known as an apodization (which means “removing the foot”, i.e. smoothing discontinuities at the beginning and end of the sampled signal) or tapering function.

[1]Blackman, R.B. and Tukey, J.W., (1958) The measurement of power spectra, Dover Publications, New York.
[2]E.R. Kanasewich, “Time Sequence Analysis in Geophysics”, The University of Alberta Press, 1975, pp. 106-108.
[3]Wikipedia, “Window function”, https://en.wikipedia.org/wiki/Window_function
[4]W.H. Press, B.P. Flannery, S.A. Teukolsky, and W.T. Vetterling, “Numerical Recipes”, Cambridge University Press, 1986, page 425.
>>> np.hanning(12)
array([0.        , 0.07937323, 0.29229249, 0.57115742, 0.82743037,
       0.97974649, 0.97974649, 0.82743037, 0.57115742, 0.29229249,
       0.07937323, 0.        ])

Plot the window and its frequency response:

>>> import matplotlib.pyplot as plt
>>> from numpy.fft import fft, fftshift
>>> window = np.hanning(51)
>>> plt.plot(window)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Hann window")
Text(0.5, 1.0, 'Hann window')
>>> plt.ylabel("Amplitude")
Text(0, 0.5, 'Amplitude')
>>> plt.xlabel("Sample")
Text(0.5, 0, 'Sample')
>>> plt.show()
>>> plt.figure()
<Figure size 640x480 with 0 Axes>
>>> A = fft(window, 2048) / 25.5
>>> mag = np.abs(fftshift(A))
>>> freq = np.linspace(-0.5, 0.5, len(A))
>>> with np.errstate(divide='ignore', invalid='ignore'):
...     response = 20 * np.log10(mag)
...
>>> response = np.clip(response, -100, 100)
>>> plt.plot(freq, response)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Frequency response of the Hann window")
Text(0.5, 1.0, 'Frequency response of the Hann window')
>>> plt.ylabel("Magnitude [dB]")
Text(0, 0.5, 'Magnitude [dB]')
>>> plt.xlabel("Normalized frequency [cycles per sample]")
Text(0.5, 0, 'Normalized frequency [cycles per sample]')
>>> plt.axis('tight')
...
>>> plt.show()
symjax.tensor.signal.kaiser(*args, **kwargs)

Return the Kaiser window.

LAX-backend implementation of kaiser(). Original docstring below.

The Kaiser window is a taper formed by using a Bessel function.

M : int
Number of points in the output window. If zero or less, an empty array is returned.
beta : float
Shape parameter for window.
out : array
The window, with the maximum value normalized to one (the value one appears only if the number of samples is odd).

bartlett, blackman, hamming, hanning

The Kaiser window is defined as

\[w(n) = I_0\left( \beta \sqrt{1-\frac{4n^2}{(M-1)^2}} \right)/I_0(\beta)\]

with

\[\quad -\frac{M-1}{2} \leq n \leq \frac{M-1}{2},\]

where \(I_0\) is the modified zeroth-order Bessel function.

The Kaiser was named for Jim Kaiser, who discovered a simple approximation to the DPSS window based on Bessel functions. The Kaiser window is a very good approximation to the Digital Prolate Spheroidal Sequence, or Slepian window, which is the transform which maximizes the energy in the main lobe of the window relative to total energy.

The Kaiser can approximate many other windows by varying the beta parameter.

beta Window shape
0 Rectangular
5 Similar to a Hamming
6 Similar to a Hanning
8.6 Similar to a Blackman

A beta value of 14 is probably a good starting point. Note that as beta gets large, the window narrows, and so the number of samples needs to be large enough to sample the increasingly narrow spike, otherwise NaNs will get returned.

Most references to the Kaiser window come from the signal processing literature, where it is used as one of many windowing functions for smoothing values. It is also known as an apodization (which means “removing the foot”, i.e. smoothing discontinuities at the beginning and end of the sampled signal) or tapering function.

[1]J. F. Kaiser, “Digital Filters” - Ch 7 in “Systems analysis by digital computer”, Editors: F.F. Kuo and J.F. Kaiser, p 218-285. John Wiley and Sons, New York, (1966).
[2]E.R. Kanasewich, “Time Sequence Analysis in Geophysics”, The University of Alberta Press, 1975, pp. 177-178.
[3]Wikipedia, “Window function”, https://en.wikipedia.org/wiki/Window_function
>>> import matplotlib.pyplot as plt
>>> np.kaiser(12, 14)
 array([7.72686684e-06, 3.46009194e-03, 4.65200189e-02, # may vary
        2.29737120e-01, 5.99885316e-01, 9.45674898e-01,
        9.45674898e-01, 5.99885316e-01, 2.29737120e-01,
        4.65200189e-02, 3.46009194e-03, 7.72686684e-06])

Plot the window and the frequency response:

>>> from numpy.fft import fft, fftshift
>>> window = np.kaiser(51, 14)
>>> plt.plot(window)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Kaiser window")
Text(0.5, 1.0, 'Kaiser window')
>>> plt.ylabel("Amplitude")
Text(0, 0.5, 'Amplitude')
>>> plt.xlabel("Sample")
Text(0.5, 0, 'Sample')
>>> plt.show()
>>> plt.figure()
<Figure size 640x480 with 0 Axes>
>>> A = fft(window, 2048) / 25.5
>>> mag = np.abs(fftshift(A))
>>> freq = np.linspace(-0.5, 0.5, len(A))
>>> response = 20 * np.log10(mag)
>>> response = np.clip(response, -100, 100)
>>> plt.plot(freq, response)
[<matplotlib.lines.Line2D object at 0x...>]
>>> plt.title("Frequency response of Kaiser window")
Text(0.5, 1.0, 'Frequency response of Kaiser window')
>>> plt.ylabel("Magnitude [dB]")
Text(0, 0.5, 'Magnitude [dB]')
>>> plt.xlabel("Normalized frequency [cycles per sample]")
Text(0.5, 0, 'Normalized frequency [cycles per sample]')
>>> plt.axis('tight')
(-0.5, 0.5, -100.0, ...) # may vary
>>> plt.show()
symjax.tensor.signal.tukey(M, alpha=0.5)[source]

Return a Tukey window, also known as a tapered cosine window. :param M: Number of points in the output window. If zero or less, an empty

array is returned.
Parameters:alpha (float, optional) – Shape parameter of the Tukey window, representing the fraction of the window inside the cosine tapered region. If zero, the Tukey window is equivalent to a rectangular window. If one, the Tukey window is equivalent to a Hann window.
Returns:w – The window, with the maximum value normalized to 1 (though the value 1 does not appear if M is even and sym is True).
Return type:ndarray

References

[1]Harris, Fredric J. (Jan 1978). “On the use of Windows for Harmonic Analysis with the Discrete Fourier Transform”. Proceedings of the IEEE 66 (1): 51-83. :doi:`10.1109/PROC.1978.10837`
[2]Wikipedia, “Window function”, https://en.wikipedia.org/wiki/Window_function#Tukey_window
symjax.tensor.signal.mfcc(signal, window, hop, n_filter, low_freq, high_freq, nyquist, n_mfcc, nfft=None, mode='valid', apod=<function _wrap_numpy_nullary_function.<locals>.wrapper>)[source]

https://librosa.github.io/librosa/_modules/librosa/feature/spectral.html#mfcc

symjax.tensor.signal.stft(signal, window, hop, apod=<function ones>, nfft=None, mode='valid')[source]

Compute the Shoft-Time-Fourier-Transform of a signal given the window length, hop and additional parameters.

Parameters:
  • signal (array) – the signal (possibly stacked of signals)
  • window (int) – the window length to be considered for the fft
  • hop (int) – the amount by which the window is moved
  • apod (func) – a function that takes an integer as inumpy.t and return the apodization window of the same length
  • nfft (int (optional)) – the number of bin that the fft on the window will use. If not given it is set the same as window.
  • mode ('valid', 'same' or 'full') – the padding of the inumpy.t signals
Returns:

output – the complex stft

Return type:

complex array

symjax.tensor.signal.dct(signal, axes=(-1, ))[source]

https://dsp.stackexchange.com/questions/2807/fast-cosine-transform-via-fft

symjax.tensor.signal.wvd(signal, window, hop, L, apod=<function _wrap_numpy_nullary_function.<locals>.wrapper>, mode='valid')[source]
symjax.tensor.signal.hilbert_transform(signal)[source]

the time should be the last dimension return the analytical signal

symjax.tensor.signal.fourier_complex_morlet(bandwidths, centers, N)[source]

Complex Morlet wavelet in Fourier

Parameters:
  • bandwidths (array) – the bandwidth of the wavelet
  • centers (array) – the centers of the wavelet
  • freqs (array (optional)) – the frequency sampling in radion going from 0 to pi and back to 0 :param N:
symjax.tensor.signal.complex_morlet(bandwidths, centers, time=None)[source]

Complex Morlet wavelet

It corresponds to with (B, C):

\phi(t) =

rac{1}{pi B} e^{- rac{t^2}{B}}e^{j2pi C t}

For a filter bank do

J = 8 Q = 1 scales = T.power(2,T.linspace(0, J, J*Q)) scales = scales[:, None] complex_morlet(scales, 1/scales)

bandwidths: array
the bandwidth of the wavelet
centers: array
the centers of the wavelet
time: array (optional)
the time sampling
wavelet: array like
the wavelet centered at 0
symjax.tensor.signal.sinc_bandpass(time, f0, f1)[source]

ensure that f0<f1 and f0>0, f1<1 whenever time is …, -1, 0, 1, …

symjax.tensor.signal.mel_filterbank(length, n_filter, low, high, nyquist)[source]
symjax.tensor.signal.convolve(in1, in2, mode='full', method='auto', precision=None)[source]

Convolve two N-dimensional arrays.

LAX-backend implementation of convolve(). Original docstring below.

Convolve in1 and in2, with the output size determined by the mode argument.

Parameters:
  • in1 (array_like) – First input.
  • in2 (array_like) – Second input. Should have the same number of dimensions as in1.
  • mode (str {'full', 'valid', 'same'}, optional) – A string indicating the size of the output:
  • method (str {'auto', 'direct', 'fft'}, optional) – A string indicating which method to use to calculate the convolution.
Returns:

convolve – An N-dimensional array containing a subset of the discrete linear convolution of in1 with in2.

Return type:

array

See also

numpy.polymul()
performs polynomial multiplication (same operation, but also accepts poly1d objects)
choose_conv_method()
chooses the fastest appropriate convolution method
fftconvolve()
Always uses the FFT method.
oaconvolve()
Uses the overlap-add method to do convolution, which is generally faster when the input arrays are large and significantly different in size.

Notes

By default, convolve and correlate use method='auto', which calls choose_conv_method to choose the fastest method using pre-computed values (choose_conv_method can also measure real-world timing with a keyword argument). Because fftconvolve relies on floating point numbers, there are certain constraints that may force method=direct (more detail in choose_conv_method docstring).

Examples

Smooth a square pulse using a Hann window:

>>> from scipy import signal
>>> sig = np.repeat([0., 1., 0.], 100)
>>> win = signal.hann(50)
>>> filtered = signal.convolve(sig, win, mode='same') / sum(win)
>>> import matplotlib.pyplot as plt
>>> fig, (ax_orig, ax_win, ax_filt) = plt.subplots(3, 1, sharex=True)
>>> ax_orig.plot(sig)
>>> ax_orig.set_title('Original pulse')
>>> ax_orig.margins(0, 0.1)
>>> ax_win.plot(win)
>>> ax_win.set_title('Filter impulse response')
>>> ax_win.margins(0, 0.1)
>>> ax_filt.plot(filtered)
>>> ax_filt.set_title('Filtered signal')
>>> ax_filt.margins(0, 0.1)
>>> fig.tight_layout()
>>> fig.show()
symjax.tensor.signal.batch_convolve(input, filter, strides=1, padding='VALID', input_format=None, filter_format=None, output_format=None, input_dilation=None, filter_dilation=None)[source]

General n-dimensional batch convolution with dilations.

Wraps Jax’s conv_general_dilated functin, and thus also the XLA’s Conv operator.

Parameters:
  • input (Tensor) – a rank n+2 dimensional input array.
  • filter (Tensor) – a rank n+2 dimensional array of kernel weights.
  • strides (int, sequence of int, optional) – a (sequence) of n integers, representing the inter-window strides. If a scalar is given, it is used n times. Defaults to 1.
  • padding (sequence of couple, ‘SAME’, ‘VALID’, optional) – a sequence of n (low, high) integer pairs that give the padding to apply before and after each spatial dimension. For ‘VALID’, those are 0. For ‘SAME’, they are the input length - filter length + 1 for each dim. Defaults to ‘Valid’.
  • input_format (None or str, optional) –

    a string of same length as the number of dimensions in input which specify their role (see below). Defaults to ‘NCW’ for 1d conv, ‘NCHW’ for 2d conv,

    and ‘NDCHW’ for 3d conv.
  • input_dilation (None, int or sequence of int, optional) – giving the dilation factor to apply in each spatial dimension of input. Inumpy.t dilation is also known as transposed convolution as it allows to increase the output spatial dimension by inserting in the input any number of `0`s between each spatial value.
  • filter_dilation (None, int or sequence of int) – giving the dilation factor to apply in each spatial dimension of filter. Filter dilation is also known as atrous convolution as it corresponds to inserting any number of `0`s in between the filter values, similar to performing the non-dilated filter convolution with a subsample version of the input across the spatial dimensions.
Returns:

An array containing the convolution result.

Return type:

Tensor

Format of input, filter and output: For example, to indicate dimension numbers consistent with the conv function with two spatial dimensions, one could use (‘NCHW’, ‘OIHW’, ‘NCHW’). As another example, to indicate dimension numbers consistent with the TensorFlow Conv2D operation, one could use (‘NHWC’, ‘HWIO’, ‘NHWC’). When using the latter form of convolution dimension specification, window strides are associated with spatial dimension character labels according to the order in which the labels appear in the rhs_spec string, so that window_strides[0] is matched with the dimension corresponding to the first character appearing in rhs_spec that is not ‘I’ or ‘O’. :param filter_format: :param output_format:

symjax.tensor.signal.convolve2d(in1, in2, mode='full', boundary='fill', fillvalue=0, precision=None)[source]

Convolve two 2-dimensional arrays.

LAX-backend implementation of convolve2d(). Original docstring below.

Convolve in1 and in2 with output size determined by mode, and boundary conditions determined by boundary and fillvalue.

Parameters:
  • in1 (array_like) – First input.
  • in2 (array_like) – Second input. Should have the same number of dimensions as in1.
  • mode (str {'full', 'valid', 'same'}, optional) – A string indicating the size of the output:
  • boundary (str {'fill', 'wrap', 'symm'}, optional) – A flag indicating how to handle boundaries:
  • fillvalue (scalar, optional) – Value to fill pad input arrays with. Default is 0.
Returns:

out – A 2-dimensional array containing a subset of the discrete linear convolution of in1 with in2.

Return type:

ndarray

Examples

Compute the gradient of an image by 2D convolution with a complex Scharr operator. (Horizontal operator is real, vertical is imaginary.) Use symmetric boundary condition to avoid creating edges at the image boundaries.

>>> from scipy import signal
>>> from scipy import misc
>>> ascent = misc.ascent()
>>> scharr = np.array([[ -3-3j, 0-10j,  +3 -3j],
...                    [-10+0j, 0+ 0j, +10 +0j],
...                    [ -3+3j, 0+10j,  +3 +3j]]) # Gx + j*Gy
>>> grad = signal.convolve2d(ascent, scharr, boundary='symm', mode='same')
>>> import matplotlib.pyplot as plt
>>> fig, (ax_orig, ax_mag, ax_ang) = plt.subplots(3, 1, figsize=(6, 15))
>>> ax_orig.imshow(ascent, cmap='gray')
>>> ax_orig.set_title('Original')
>>> ax_orig.set_axis_off()
>>> ax_mag.imshow(np.absolute(grad), cmap='gray')
>>> ax_mag.set_title('Gradient magnitude')
>>> ax_mag.set_axis_off()
>>> ax_ang.imshow(np.angle(grad), cmap='hsv') # hsv is cyclic, like angles
>>> ax_ang.set_title('Gradient orientation')
>>> ax_ang.set_axis_off()
>>> fig.show()
symjax.tensor.signal.correlate(in1, in2, mode='full', method='auto', precision=None)[source]

Cross-correlate two N-dimensional arrays.

LAX-backend implementation of correlate(). Original docstring below.

Cross-correlate in1 and in2, with the output size determined by the mode argument.

Parameters:
  • in1 (array_like) – First input.
  • in2 (array_like) – Second input. Should have the same number of dimensions as in1.
  • mode (str {'full', 'valid', 'same'}, optional) – A string indicating the size of the output:
  • method (str {'auto', 'direct', 'fft'}, optional) – A string indicating which method to use to calculate the correlation.
Returns:

correlate – An N-dimensional array containing a subset of the discrete linear cross-correlation of in1 with in2.

Return type:

array

See also

choose_conv_method()
contains more documentation on method.

Notes

The correlation z of two d-dimensional arrays x and y is defined as:

z[...,k,...] = sum[..., i_l, ...] x[..., i_l,...] * conj(y[..., i_l - k,...])

This way, if x and y are 1-D arrays and z = correlate(x, y, 'full') then

\[z[k] = (x * y)(k - N + 1) = \sum_{l=0}^{||x||-1}x_l y_{l-k+N-1}^{*}\]

for \(k = 0, 1, ..., ||x|| + ||y|| - 2\)

where \(||x||\) is the length of x, \(N = \max(||x||,||y||)\), and \(y_m\) is 0 when m is outside the range of y.

method='fft' only works for numerical arrays as it relies on fftconvolve. In certain cases (i.e., arrays of objects or when rounding integers can lose precision), method='direct' is always used.

When using “same” mode with even-length inputs, the outputs of correlate and correlate2d differ: There is a 1-index offset between them.

Examples

Implement a matched filter using cross-correlation, to recover a signal that has passed through a noisy channel.

>>> from scipy import signal
>>> sig = np.repeat([0., 1., 1., 0., 1., 0., 0., 1.], 128)
>>> sig_noise = sig + np.random.randn(len(sig))
>>> corr = signal.correlate(sig_noise, np.ones(128), mode='same') / 128
>>> import matplotlib.pyplot as plt
>>> clock = np.arange(64, len(sig), 128)
>>> fig, (ax_orig, ax_noise, ax_corr) = plt.subplots(3, 1, sharex=True)
>>> ax_orig.plot(sig)
>>> ax_orig.plot(clock, sig[clock], 'ro')
>>> ax_orig.set_title('Original signal')
>>> ax_noise.plot(sig_noise)
>>> ax_noise.set_title('Signal with noise')
>>> ax_corr.plot(corr)
>>> ax_corr.plot(clock, corr[clock], 'ro')
>>> ax_corr.axhline(0.5, ls=':')
>>> ax_corr.set_title('Cross-correlated with rectangular pulse')
>>> ax_orig.margins(0, 0.1)
>>> fig.tight_layout()
>>> fig.show()
symjax.tensor.signal.correlate2d(in1, in2, mode='full', boundary='fill', fillvalue=0, precision=None)[source]

Cross-correlate two 2-dimensional arrays.

LAX-backend implementation of correlate2d(). Original docstring below.

Cross correlate in1 and in2 with output size determined by mode, and boundary conditions determined by boundary and fillvalue.

Parameters:
  • in1 (array_like) – First input.
  • in2 (array_like) – Second input. Should have the same number of dimensions as in1.
  • mode (str {'full', 'valid', 'same'}, optional) – A string indicating the size of the output:
  • boundary (str {'fill', 'wrap', 'symm'}, optional) – A flag indicating how to handle boundaries:
  • fillvalue (scalar, optional) – Value to fill pad input arrays with. Default is 0.
Returns:

correlate2d – A 2-dimensional array containing a subset of the discrete linear cross-correlation of in1 with in2.

Return type:

ndarray

Notes

When using “same” mode with even-length inputs, the outputs of correlate and correlate2d differ: There is a 1-index offset between them.

Examples

Use 2D cross-correlation to find the location of a template in a noisy image:

>>> from scipy import signal
>>> from scipy import misc
>>> face = misc.face(gray=True) - misc.face(gray=True).mean()
>>> template = np.copy(face[300:365, 670:750])  # right eye
>>> template -= template.mean()
>>> face = face + np.random.randn(*face.shape) * 50  # add noise
>>> corr = signal.correlate2d(face, template, boundary='symm', mode='same')
>>> y, x = np.unravel_index(np.argmax(corr), corr.shape)  # find the match
>>> import matplotlib.pyplot as plt
>>> fig, (ax_orig, ax_template, ax_corr) = plt.subplots(3, 1,
...                                                     figsize=(6, 15))
>>> ax_orig.imshow(face, cmap='gray')
>>> ax_orig.set_title('Original')
>>> ax_orig.set_axis_off()
>>> ax_template.imshow(template, cmap='gray')
>>> ax_template.set_title('Template')
>>> ax_template.set_axis_off()
>>> ax_corr.imshow(corr, cmap='gray')
>>> ax_corr.set_title('Cross-correlation')
>>> ax_corr.set_axis_off()
>>> ax_orig.plot(x, y, 'ro')
>>> fig.show()
symjax.tensor.signal.extract_signal_patches()
symjax.tensor.signal.extract_image_patches(image, window_shape, strides=1, data_format='channel_first', padding='valid', flatten_patch=False)

extract patches from an input tensor

image: Tensor-like
the input to extract patches from
window_shape: int or tuple of ints
the spatial shape of the patch to extract
hop: int or tuple of ints
the spatial hop of the patch to extract
data_format: str
either channel_first or channel_last
padding: str
either same or valid
flatten_patch: bool
whether to return patches as flattened or not