scipp.signal.sosfiltfilt#

scipp.signal.sosfiltfilt(obj, dim, *, sos, **kwargs)#

A forward-backward digital filter using cascaded second-order sections.

This is a wrapper around scipy.signal.sosfiltfilt(). See there for a complete description of parameters. The differences are:

  • Instead of an array x and an optional axis, the input must be a data array (or variable) and dimension label. If it is a variable, the coord used for setting up the second-order sections is used as output coordinate.

  • The array of second-order filter coefficients sos includes the coordinate used for computing the frequency used for creating the coefficients. This is compared to the corresponding coordinate of the input data array to ensure that compatible coefficients are used.

Examples:

>>> from scipp.signal import butter, sosfiltfilt
>>> x = sc.linspace(dim='x', start=1.1, stop=4.0, num=1000, unit='m')
>>> y = sc.sin(x * sc.scalar(1.0, unit='rad/m'))
>>> y += sc.sin(x * sc.scalar(400.0, unit='rad/m'))
>>> da = sc.DataArray(data=y, coords={'x': x})
>>> sos = butter(da.coords['x'], N=4, Wn=20 / x.unit)
>>> out = sosfiltfilt(da, 'x', sos=sos)
>>> sc.plot({'input':da, 'sosfiltfilt':out})
../../_images/scipp-signal-sosfiltfilt-1.png

Instead of calling sosfiltfilt the more convenient filtfilt method of scipp.signal.SOS can be used:

>>> out = butter(da.coords['x'], N=4, Wn=20 / x.unit).filtfilt(da, 'x')
Return type

DataArray