Signal processing

Signal processing toolbox.

arlpy.signal.bb2pb(x, fd, fc, fs=None, axis=-1)

Convert baseband signal to passband.

For communication applications, one may wish to use arlpy.comms.upconvert() instead, as that function supports pulse shaping.

The convention used in that exp(2j*pi*fc*t) is a positive frequency carrier.

Parameters:
  • x – complex baseband signal
  • fd – sampling rate of baseband signal in Hz
  • fc – carrier frequency in passband in Hz
  • fs – sampling rate of passband signal in Hz (None => same as fd)
  • axis – axis of the signal, if multiple signals specified
Returns:

real passband signal, sampled at fs

arlpy.signal.correlate_periodic(a, v=None)

Cross-correlation of two 1-dimensional periodic sequences.

a and v must be sequences with the same length. If v is not specified, it is assumed to be the same as a (i.e. the function computes auto-correlation).

Parameters:
  • a – input sequence #1
  • v – input sequence #2
Returns:

discrete periodic cross-correlation of a and v

arlpy.signal.cw(fc, duration, fs, window=None, complex_output=False)

Generate a sinusoidal pulse.

Parameters:
  • fc – frequency of the pulse in Hz
  • duration – duration of the pulse in s
  • fs – sampling rate in Hz
  • window – window function to use (None means rectangular window)
  • complex_output – True to return complex signal, False for a real signal

For supported window functions, see documentation for scipy.signal.get_window().

>>> import arlpy
>>> x1 = arlpy.signal.cw(fc=27000, duration=0.5, fs=250000)
>>> x2 = arlpy.signal.cw(fc=27000, duration=0.5, fs=250000, window='hamming')
>>> x3 = arlpy.signal.cw(fc=27000, duration=0.5, fs=250000, window=('kaiser', 4.0))
arlpy.signal.detect_impulses(x, fs, k=10, tdist=0.001)

Detect impulses in x

The minimum height of impulses is defined by a+k*b where a is median of the envelope of x and b is median absolute deviation (MAD) of the envelope of x.

Parameters:
  • x – real signal
  • fs – sampling frequency in Hz
  • k – multiple of MAD for the impulse minimum height (default: 10)
  • tdist – minimum time difference between neighbouring impulses in sec (default: 1e-3)
Returns:

indices and heights of detected impulses

>>> nsamp = 1000
>>> ind_impulses = np.array([10, 115, 641, 888])
>>> x = np.zeros((nsamp))
>>> x[ind_impulses] = 1
>>> x += np.random.normal(0, 0.1, nsamp)
>>> ind_pks, h_pks = signal.detect_impulses(x, fs=100000, k=10, tdist=1e-3)
arlpy.signal.envelope(x, axis=-1)

Generate a Hilbert envelope of the real signal x.

Parameters:
  • x – real passband signal
  • axis – axis of the signal, if multiple signals specified
arlpy.signal.gmseq(spec, theta=None)

Generate generalized m-sequence.

Generalized m-sequences are related to m-sequences but have an additional parameter \(\theta\). When \(\theta = \pi/2\), generalized m-sequences become normal m-sequences. When \(\theta < \pi/2\), generalized m-sequences contain a DC-component that leads to an exalted carrier after modulation.

When theta is \(\arctan(\sqrt{n})\) where \(n\) is the length of the m-sequence, the m-sequence is considered to be period matched. Period matched m-sequences are complex sequences with perfect discrete periodic auto-correlation properties, i.e., all non-zero lag periodic auto-correlations are zero. The zero-lag autocorrelation is \(n = 2^m-1\), where m is the shift register length.

This function currently supports shift register lengths between 2 and 30.

Parameters:
  • spec – m-sequence specifier (shift register length or taps)
  • theta – transmission angle (None to use period-matched angle)
>>> import arlpy
>>> x = arlpy.signal.gmseq(7)
>>> len(x)
127
arlpy.signal.goertzel(f, x, fs=2.0, filter=False)

Goertzel algorithm for single tone detection.

The output of the Goertzel algorithm is the same as a single bin DFT if f/(fs/N) is an integer, where N is the number of points in signal x.

The detection metric returned by this function is the magnitude of the output of the Goertzel algorithm at the end of the input block. If filter is set to true, the complex time series at the output of the IIR filter is returned, rather than just the detection metric.

Parameters:
  • f – frequency of tone of interest in Hz
  • x – real or complex input sequence
  • fs – sampling frequency of x in Hz
  • filter – output complex time series if true, detection metric otherwise (default: false)
Returns:

detection metric or complex time series

>>> import arlpy
>>> x1 = arlpy.signal.cw(64, 1, 512)
>>> g1 = arlpy.signal.goertzel(64, x1, 512)
>>> g1
256.0
>>> g2 = arlpy.signal.goertzel(32, x1, 512)
>>> g2
0.0
arlpy.signal.lfilter0(b, a, x, axis=-1)

Filter data with an IIR or FIR filter with zero DC group delay.

scipy.signal.lfilter() provides a way to filter a signal x using a FIR/IIR filter defined by b and a. The resulting output is delayed, as compared to the input by the group delay. This function corrects for the group delay, resulting in an output that is synchronized with the input signal. If the filter as an acausal impulse response, some precursor signal from the output will be lost. To avoid this, pad input signal x with sufficient zeros at the beginning to capture the precursor. Since both, scipy.signal.lfilter() and this function return a signal with the same length as the input, some signal tail is lost at the end. To avoid this, pad input signal x with sufficient zeros at the end.

See documentation for scipy.signal.lfilter() for more details.

>>> import arlpy
>>> import numpy as np
>>> fs = 250000
>>> b = arlpy.uwa.absorption_filter(fs, distance=500)
>>> x = np.pad(arlpy.signal.sweep(20000, 40000, 0.5, fs), (127, 127), 'constant')
>>> y = arlpy.signal.lfilter0(b, 1, x)
arlpy.signal.lfilter_gen(b, a)

Generator form of an FIR/IIR filter.

The filter is a direct form I implementation of the standard difference equation. Data samples can be passed to the filter using the send() method, and the output can be read a sample at a time.

>>> import arlpy
>>> import numpy as np
>>> import scipy.signal as sp
>>> b, a = sp.iirfilter(2, 0.1, btype='lowpass')  # generate a biquad lowpass
>>> f = arlpy.signal.filter_gen(b, a)             # create the filter
>>> x = np.random.normal(0, 1, 1000)              # get some random data
>>> y = [f.send(v) for v in x]                    # filter data by stepping through it
arlpy.signal.mfilter(s, x, complex_output=False, axis=-1)

Matched filter recevied signal using a reference signal.

Parameters:
  • s – reference signal
  • x – recevied signal
  • complex_output – True to return complex signal, False for absolute value of complex signal
  • axis – axis of the signal, if multiple recevied signals specified
arlpy.signal.mseq(spec, n=None)

Generate m-sequence.

m-sequences are sequences of \(\pm 1\) values with near-perfect discrete periodic auto-correlation properties. All non-zero lag periodic auto-correlations are -1. The zero-lag autocorrelation is \(2^m-1\), where m is the shift register length.

This function currently supports shift register lengths between 2 and 30.

Parameters:
  • spec – m-sequence specifier (shift register length or taps)
  • n – length of sequence (None means full length of \(2^m-1\))
>>> import arlpy
>>> x = arlpy.signal.mseq(7)
>>> len(x)
127
arlpy.signal.nco(fc, fs=2.0, phase0=0, wrap=6.283185307179586, func=<function <lambda>>)

Numerically controlled oscillator (NCO).

If fs is specified, fc is given in Hz, otherwise it is specified as normalized frequency (Nyquist = 1).

The default oscillator function is exp(i*phase) to generate a complex sinusoid. Alternate oscillator functions that take in the phase angle and generate other outputs can be specifed. For example, a real sinusoid can be generated by specifying sin as the function. The phase angle can be generated by specifying None as the function.

Parameters:
  • fc – array of instantaneous oscillation frequency
  • fs – sampling frequency in Hz
  • phase0 – initial phase in radians (default: 0)
  • wrap – phase angle to wrap phase around to 0 (default: \(2\pi\))
  • func – oscillator function of phase angle (default: complex sinusoid)
>>> import arlpy
>>> import numpy as np
>>> fc = np.append([27000]*12, [54000]*5)
>>> x = arlpy.signal.nco(fc, 108000, func=np.sin)
>>> x
[0, 1, 0, -1, 0, 1, 0, -1, 0, 1, 0, -1, 1, -1, 1, -1, 1]
arlpy.signal.nco_gen(fc, fs=2.0, phase0=0, wrap=6.283185307179586, func=<function <lambda>>)

Generator form of a numerically controlled oscillator (NCO).

Samples at the output of the oscillator can be generated using the next() method. The oscillator frequency can be modified during operation using the send() method, with fc as the argument.

If fs is specified, fc is given in Hz, otherwise it is specified as normalized frequency (Nyquist = 1).

The default oscillator function is exp(i*phase) to generate a complex sinusoid. Alternate oscillator functions that take in the phase angle and generate other outputs can be specifed. For example, a real sinusoid can be generated by specifying sin as the function. The phase angle can be generated by specifying None as the function.

Parameters:
  • fc – oscillation frequency
  • fs – sampling frequency in Hz
  • phase0 – initial phase in radians (default: 0)
  • wrap – phase angle to wrap phase around to 0 (default: \(2\pi\))
  • func – oscillator function of phase angle (default: complex sinusoid)
>>> import arlpy
>>> import math
>>> nco = arlpy.signal.nco_gen(27000, 108000, func=math.sin)
>>> x = [nco.next() for i in range(12)]
>>> x = np.append(x, nco.send(54000))      # change oscillation frequency
>>> x = np.append(x, [nco.next() for i in range(4)])
>>> x
[0, 1, 0, -1, 0, 1, 0, -1, 0, 1, 0, -1, 1, -1, 1, -1, 1]
arlpy.signal.pb2bb(x, fs, fc, fd=None, flen=127, cutoff=None, axis=-1)

Convert passband signal to baseband.

The baseband conversion uses a low-pass filter after downconversion, with a default cutoff frequency of 0.6*fd, if fd is specified, or 1.1*fc if fd is not specified. Alternatively, the user may specify the cutoff frequency explicitly.

For communication applications, one may wish to use arlpy.comms.downconvert() instead, as that function supports matched filtering with a pulse shape rather than a generic low-pass filter.

The convention used in that exp(2j*pi*fc*t) is a positive frequency carrier.

Parameters:
  • x – passband signal
  • fs – sampling rate of passband signal in Hz
  • fc – carrier frequency in passband in Hz
  • fd – sampling rate of baseband signal in Hz (None => same as fs)
  • flen – number of taps in the low-pass FIR filter
  • cutoff – cutoff frequency in Hz (None means auto-select)
  • axis – axis of the signal, if multiple signals specified
Returns:

complex baseband signal, sampled at fd

arlpy.signal.sweep(f1, f2, duration, fs, method='linear', window=None)

Generate frequency modulated sweep.

Parameters:
  • f1 – starting frequency in Hz
  • f2 – ending frequency in Hz
  • duration – duration of the pulse in s
  • fs – sampling rate in Hz
  • method – type of sweep ('linear', 'quadratic', 'logarithmic', 'hyperbolic')
  • window – window function to use (None means rectangular window)

For supported window functions, see documentation for scipy.signal.get_window().

>>> import arlpy
>>> x1 = arlpy.signal.sweep(20000, 30000, duration=0.5, fs=250000)
>>> x2 = arlpy.signal.sweep(20000, 30000, duration=0.5, fs=250000, window='hamming')
>>> x2 = arlpy.signal.sweep(20000, 30000, duration=0.5, fs=250000, window=('kaiser', 4.0))
arlpy.signal.time(n, fs)

Generate a time vector for time series.

Parameters:
  • n – time series, or number of samples
  • fs – sampling rate in Hz
Returns:

time vector starting at time 0

>>> import arlpy
>>> t = arlpy.signal.time(100000, fs=250000)
>>> t
array([  0.00000000e+00,   4.00000000e-06,   8.00000000e-06, ...,
     3.99988000e-01,   3.99992000e-01,   3.99996000e-01])
>>> x = arlpy.signal.cw(fc=27000, duration=0.5, fs=250000)
>>> t = arlpy.signal.time(x, fs=250000)
>>> t
array([  0.00000000e+00,   4.00000000e-06,   8.00000000e-06, ...,
     4.99988000e-01,   4.99992000e-01,   4.99996000e-01])