| filter {gsignal} | R Documentation | 
Filter a signal
Description
Apply a 1-D digital filter compatible with 'Matlab' and 'Octave'.
Usage
filter(filt, ...)
## Default S3 method:
filter(filt, a, x, zi = NULL, ...)
## S3 method for class 'Arma'
filter(filt, x, ...)
## S3 method for class 'Ma'
filter(filt, x, ...)
## S3 method for class 'Sos'
filter(filt, x, ...)
## S3 method for class 'Zpg'
filter(filt, x, ...)
## S3 method for class 'Zpg'
filter(filt, x, ...)
Arguments
| filt | For the default case, the moving-average coefficients of an ARMA
filter (normally called  | 
| ... | additional arguments (ignored). | 
| a | the autoregressive (recursive) coefficients of an ARMA filter,
specified as a numeric or complex vector. If  | 
| x | the input signal to be filtered, specified as a numeric or complex
vector or matrix. If  | 
| zi | If  | 
Details
The filter is a direct form II transposed implementation of the standard linear time-invariant difference equation:
N M SUM a(k+1)y(n-k) + SUM b(k+1)x(n-k) = 0; 1 <= n <= length(x) k=0 k=0
where N = length(a) - 1 and M = length(b) - 1.
The initial and final conditions for filter delays can be used to filter data in sections, especially if memory limitations are a consideration. See the examples.
Value
The filtered signal, of the same dimensions as the input signal. In
case the zi input argument was specified, a list with two elements
is returned containing the variables y, which represents the output
signal, and zf, which contains the final state vector or matrix.
Author(s)
Geert van Boxtel, G.J.M.vanBoxtel@gmail.com.
See Also
filter_zi, sosfilt (preferred because it
avoids numerical problems).
Examples
bf <- butter(3, 0.1)                                 # 10 Hz low-pass filter
t <- seq(0, 1, len = 100)                            # 1 second sample
x <- sin(2* pi * t * 2.3) + 0.25 * rnorm(length(t))  # 2.3 Hz sinusoid+noise
z <- filter(bf, x)                                   # apply filter
plot(t, x, type = "l")
lines(t, z, col = "red")
## specify initial conditions
## from Python scipy.signal.lfilter() documentation
t <- seq(-1, 1, length.out =  201)
x <- (sin(2 * pi * 0.75 * t * (1 - t) + 2.1)
      + 0.1 * sin(2 * pi * 1.25 * t + 1)
      + 0.18 * cos(2 * pi * 3.85 * t))
h <- butter(3, 0.05)
lab <- max(length(h$b), length(h$a)) - 1
zi <- filtic(h$b, h$a, rep(1, lab), rep(1, lab))
z1 <- filter(h, x)
z2 <- filter(h, x, zi * x[1])
plot(t, x, type = "l")
lines(t, z1, col = "red")
lines(t, z2$y, col = "green")
legend("bottomright", legend = c("Original signal",
        "Filtered without initial conditions",
        "Filtered with initial conditions"),
       lty = 1, col = c("black", "red", "green"))