Signal

The signal module constains all kinds of signal processing related functions.

Filtering

class acoustics.signal.Filterbank(frequencies, sample_frequency=44100, order=8)[source]

Fractional-Octave filter bank.

Warning

For high frequencies the filter coefficients are wrong for low frequencies. Therefore, to improve the response for lower frequencies the signal should be downsampled. Currently, there is no easy way to do so within the Filterbank.

__init__(frequencies, sample_frequency=44100, order=8)[source]

Initialize self. See help(type(self)) for accurate signature.

filters

Filters this filterbank consists of.

filtfilt(signal)[source]

Filter signal with filterbank. Returns a list consisting of a filtered signal per filter.

Note

This function uses scipy.signal.filtfilt() and therefore has a zero-phase response.

frequencies = None

Frequencies object.

See also Frequencies and subclasses.

Note

A frequencies object should have the attributes center, lower and upper.

lfilter(signal)[source]

Filter signal with filterbank.

Note

This function uses scipy.signal.lfilter().

order = None

Filter order of Butterworth filter.

plot_power(signal)[source]

Plot power in signal.

plot_response()[source]

Plot frequency response.

Note

The follow phase response is obtained in case lfilter() is used. The method filtfilt() results in a zero-phase response.

power(signal)[source]

Power per band in signal.

sample_frequency

Sample frequency.

acoustics.signal.bandpass_filter(lowcut, highcut, fs, order=8, output='sos')[source]

Band-pass filter.

Parameters:
  • lowcut – Lower cut-off frequency
  • highcut – Upper cut-off frequency
  • fs – Sample frequency
  • order – Filter order
  • output – Output type. {‘ba’, ‘zpk’, ‘sos’}. Default is ‘sos’. See also scipy.signal.butter().
Returns:

Returned value depends on output.

A Butterworth filter is used.

acoustics.signal.octave_filter(center, fs, fraction, order=8, output='sos')[source]

Fractional-octave band-pass filter.

Parameters:
  • center – Centerfrequency of fractional-octave band.
  • fs – Sample frequency
  • fraction – Fraction of fractional-octave band.
  • order – Filter order
  • output – Output type. {‘ba’, ‘zpk’, ‘sos’}. Default is ‘sos’. See also scipy.signal.butter().

A Butterworth filter is used.

acoustics.signal.bandpass(signal, lowcut, highcut, fs, order=8, zero_phase=False)[source]

Filter signal with band-pass filter.

Parameters:
  • signal – Signal
  • lowcut – Lower cut-off frequency
  • highcut – Upper cut-off frequency
  • fs – Sample frequency
  • order – Filter order
  • zero_phase – Prevent phase error by filtering in both directions (filtfilt)

A Butterworth filter is used. Filtering is done with second-order sections.

See also

bandpass_filter() for the filter that is used.

acoustics.signal.lowpass(signal, cutoff, fs, order=4, zero_phase=False)[source]

Filter signal with low-pass filter.

Parameters:
  • signal – Signal
  • fs – Sample frequency
  • cutoff – Cut-off frequency
  • order – Filter order
  • zero_phase – Prevent phase error by filtering in both directions (filtfilt)

A Butterworth filter is used. Filtering is done with second-order sections.

acoustics.signal.highpass(signal, cutoff, fs, order=4, zero_phase=False)[source]

Filter signal with low-pass filter.

Parameters:
  • signal – Signal
  • fs – Sample frequency
  • cutoff – Cut-off frequency
  • order – Filter order
  • zero_phase – Prevent phase error by filtering in both directions (filtfilt)

A Butterworth filter is used. Filtering is done with second-order sections.

acoustics.signal.octavepass(signal, center, fs, fraction, order=8, zero_phase=True)[source]

Filter signal with fractional-octave bandpass filter.

Parameters:
  • signal – Signal
  • center – Centerfrequency of fractional-octave band.
  • fs – Sample frequency
  • fraction – Fraction of fractional-octave band.
  • order – Filter order
  • zero_phase – Prevent phase error by filtering in both directions (filtfilt)

A Butterworth filter is used. Filtering is done with second-order sections.

See also

octave_filter()

acoustics.signal.convolve(signal, ltv, mode='full')[source]

Perform convolution of signal with linear time-variant system ltv.

Parameters:
  • signal – Vector representing input signal \(u\).
  • ltv – 2D array where each column represents an impulse response
  • mode – ‘full’, ‘valid’, or ‘same’. See np.convolve() for an explanation of the options.

The convolution of two sequences is given by

\[\mathbf{y} = \mathbf{t} \star \mathbf{u}\]

This can be written as a matrix-vector multiplication

\[\mathbf{y} = \mathbf{T} \cdot \mathbf{u}\]

where \(T\) is a Toeplitz matrix in which each column represents an impulse response. In the case of a linear time-invariant (LTI) system, each column represents a time-shifted copy of the first column. In the time-variant case (LTV), every column can contain a unique impulse response, both in values as in size.

This function assumes all impulse responses are of the same size. The input matrix ltv thus represents the non-shifted version of the Toeplitz matrix.

See also

np.convolve(), scipy.signal.convolve() and scipy.signal.fftconvolve() for convolution with LTI system.

Windowing

acoustics.signal.window_scaling_factor(window, axis=-1)[source]

Calculate window scaling factor.

Parameters:window – Window.

When analysing broadband (filtered noise) signals it is common to normalize the windowed signal so that it has the same power as the un-windowed one.

\[S = \sqrt{\frac{\sum_{i=0}^N w_i^2}{N}}\]
acoustics.signal.apply_window(x, window)[source]

Apply window to signal.

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • window – Vector representing window.
Returns:

Signal with window applied to it.

\[x_s(t) = x(t) / S\]

where \(S\) is the window scaling factor.

Spectra

Different types of spectra exist.

acoustics.signal.amplitude_spectrum(x, fs, N=None)[source]

Amplitude spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

The amplitude spectrum gives the amplitudes of the sinusoidal the signal is built up from, and the RMS (root-mean-square) amplitudes can easily be found by dividing these amplitudes with \(\sqrt{2}\).

The amplitude spectrum is double-sided.

acoustics.signal.auto_spectrum(x, fs, N=None)[source]

Auto-spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

The auto-spectrum contains the squared amplitudes of the signal. Squared amplitudes are used when presenting data as it is a measure of the power/energy in the signal.

\[S_{xx} (f_n) = \overline{X (f_n)} \cdot X (f_n)\]

The auto-spectrum is double-sided.

acoustics.signal.power_spectrum(x, fs, N=None)[source]

Power spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

The power spectrum, or single-sided autospectrum, contains the squared RMS amplitudes of the signal.

A power spectrum is a spectrum with squared RMS values. The power spectrum is calculated from the autospectrum of the signal.

Warning

Does not include scaling to reference value!

See also

auto_spectrum()

acoustics.signal.density_spectrum(x, fs, N=None)[source]

Density spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

A density spectrum considers the amplitudes per unit frequency. Density spectra are used to compare spectra with different frequency resolution as the magnitudes are not influenced by the resolution because it is per Hertz. The amplitude spectra on the other hand depend on the chosen frequency resolution.

acoustics.signal.angle_spectrum(x, fs, N=None)[source]

Phase angle spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

This function returns a single-sided wrapped phase angle spectrum.

See also

phase_spectrum() for unwrapped phase spectrum.

acoustics.signal.phase_spectrum(x, fs, N=None)[source]

Phase spectrum of instantaneous signal \(x(t)\).

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency \(f_s\).
  • N – Amount of FFT bins.

This function returns a single-sided unwrapped phase spectrum.

See also

angle_spectrum() for wrapped phase angle.

Frequency bands

class acoustics.signal.Frequencies(center, lower, upper, bandwidth=None)[source]

Object describing frequency bands.

__init__(center, lower, upper, bandwidth=None)[source]

Initialize self. See help(type(self)) for accurate signature.

angular()[source]

Angular center frequency in radians per second.

bandwidth = None

Bandwidth.

center = None

Center frequencies.

lower = None

Lower frequencies.

upper = None

Upper frequencies.

class acoustics.signal.EqualBand(center=None, fstart=None, fstop=None, nbands=None, bandwidth=None)[source]

Equal bandwidth spectrum. Generally used for narrowband data.

__init__(center=None, fstart=None, fstop=None, nbands=None, bandwidth=None)[source]
Parameters:
  • center – Vector of center frequencies.
  • fstart – First center frequency.
  • fstop – Last center frequency.
  • nbands – Amount of frequency bands.
  • bandwidth – Bandwidth of bands.
class acoustics.signal.OctaveBand(center=None, fstart=None, fstop=None, nbands=None, fraction=1, reference=1000.0)[source]

Fractional-octave band spectrum.

__init__(center=None, fstart=None, fstop=None, nbands=None, fraction=1, reference=1000.0)[source]

Initialize self. See help(type(self)) for accurate signature.

fraction = None

Fraction of fractional-octave filter.

nominal = None

Nominal center frequencies.

reference = None

Reference center frequency.

acoustics.signal.integrate_bands(data, a, b)[source]

Reduce frequency resolution of power spectrum. Merges frequency bands by integration.

Parameters:

Note

Needs rewriting so that the summation goes over axis=1.

acoustics.signal.octaves(p, fs, density=False, frequencies=array([ 16., 31.5, 63., 125., 250., 500., 1000., 2000., 4000., 8000., 16000. ]), ref=2e-05)[source]

Calculate level per 1/1-octave in frequency domain using the FFT.

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency.
  • density – Power density instead of power.
  • frequencies – Frequencies.
  • ref – Reference value.
Returns:

Tuple. First element is an instance of OctaveBand. The second element an array.

Note

Based on power spectrum (FFT)

See also

acoustics.bands.OCTAVE_CENTER_FREQUENCIES

Note

Exact center frequencies are always calculated.

acoustics.signal.third_octaves(p, fs, density=False, frequencies=array([1.00e+01, 1.25e+01, 1.60e+01, 2.00e+01, 2.50e+01, 3.15e+01, 4.00e+01, 5.00e+01, 6.30e+01, 8.00e+01, 1.00e+02, 1.25e+02, 1.60e+02, 2.00e+02, 2.50e+02, 3.15e+02, 4.00e+02, 5.00e+02, 6.30e+02, 8.00e+02, 1.00e+03, 1.25e+03, 1.60e+03, 2.00e+03, 2.50e+03, 3.15e+03, 4.00e+03, 5.00e+03, 6.30e+03, 8.00e+03, 1.00e+04, 1.25e+04, 1.60e+04, 2.00e+04]), ref=2e-05)[source]

Calculate level per 1/3-octave in frequency domain using the FFT.

Parameters:
  • x – Instantaneous signal \(x(t)\).
  • fs – Sample frequency.
  • density – Power density instead of power.
Returns:

Tuple. First element is an instance of OctaveBand. The second element an array.

Note

Based on power spectrum (FFT)

See also

acoustics.bands.THIRD_OCTAVE_CENTER_FREQUENCIES

Note

Exact center frequencies are always calculated.

Hilbert transform

acoustics.signal.amplitude_envelope(signal, fs, axis=-1)[source]

Instantaneous amplitude of tone.

The instantaneous amplitude is the magnitude of the analytic signal.

Parameters:
  • signal – Signal.
  • fs – Sample frequency.
  • axis – Axis.
Returns:

Amplitude envelope of signal.

acoustics.signal.instantaneous_phase(signal, fs, axis=-1)[source]

Instantaneous phase of tone.

Parameters:
  • signal – Signal.
  • fs – Sample frequency.
  • axis – Axis.
Returns:

Instantaneous phase of signal.

The instantaneous phase is the angle of the analytic signal. This function returns a wrapped angle.

acoustics.signal.instantaneous_frequency(signal, fs, axis=-1)[source]

Determine instantaneous frequency of tone.

Parameters:
  • signal – Signal.
  • fs – Sample frequency.
  • axis – Axis.
Returns:

Instantaneous frequency of signal.

The instantaneous frequency can be obtained by differentiating the unwrapped instantaneous phase.

Conversion

acoustics.signal.decibel_to_neper(decibel)[source]

Convert decibel to neper.

Parameters:decibel – Value in decibel (dB).
Returns:Value in neper (Np).

The conversion is done according to

\[\mathrm{dB} = \frac{\log{10}}{20} \mathrm{Np}\]
acoustics.signal.neper_to_decibel(neper)[source]

Convert neper to decibel.

Parameters:neper – Value in neper (Np).
Returns:Value in decibel (dB).

The conversion is done according to

\[\mathrm{Np} = \frac{20}{\log{10}} \mathrm{dB}\]

Other

acoustics.signal.isolate(signals)[source]

Isolate signals.

Parameters:signals – Array of shape N x M where N is the amount of samples and M the amount of signals. Thus, each column is a signal.
Returns:Array of isolated signals. Each column is a signal.

Isolate signals using Singular Value Decomposition.

acoustics.signal.zero_crossings(data)[source]

Determine the positions of zero crossings in data.

Parameters:data – Vector
Returns:Vector with indices of samples before the zero crossing.
acoustics.signal.rms(x)[source]

Root mean squared of signal x.

Parameters:x – Dynamic quantity.
\[\begin{split}x_{rms} = lim_{T \\to \\infty} \\sqrt{\\frac{1}{T} \int_0^T |f(x)|^2 \\mathrm{d} t }\end{split}\]
Seealso:ms().
acoustics.signal.ms(x)[source]

Mean value of signal x squared.

Parameters:x – Dynamic quantity.
Returns:Mean squared of x.
acoustics.signal.normalize(y, x=None)[source]

normalize power in y to a (standard normal) white noise signal.

Optionally normalize to power in signal x.

#The mean power of a Gaussian with \(\mu=0\) and \(\sigma=1\) is 1.

acoustics.signal.ir2fr(ir, fs, N=None)[source]

Convert impulse response into frequency response. Returns single-sided RMS spectrum.

Parameters:
  • ir – Impulser response
  • fs – Sample frequency
  • N – Blocks

Calculates the positive frequencies using np.fft.rfft(). Corrections are then applied to obtain the single-sided spectrum.

Note

Single-sided spectrum. Therefore, the amount of bins returned is either N/2 or N/2+1.

acoustics.signal.wvd(signal, fs, analytic=True)[source]

Wigner-Ville Distribution

Parameters:
  • signal – Signal
  • fs – Sample frequency
  • analytic – Use the analytic signal, calculated using Hilbert transform.
\[W_z(n, \omega) = 2 \sum_k z^*[n-k]z[n+k] e^{-j\omega 2kT}\]

Includes positive and negative frequencies.