v2.0.0
Loading...
Searching...
No Matches
iirfilter.cpp
Go to the documentation of this file.
1//=============================================================================================================
33
34//=============================================================================================================
35// INCLUDES
36//=============================================================================================================
37
38#include "iirfilter.h"
39
40//=============================================================================================================
41// QT INCLUDES
42//=============================================================================================================
43
44#include <QDebug>
45
46//=============================================================================================================
47// C++ INCLUDES
48//=============================================================================================================
49
50#include <cmath>
51#include <complex>
52#include <algorithm>
53
54//=============================================================================================================
55// USED NAMESPACES
56//=============================================================================================================
57
58using namespace UTILSLIB;
59using namespace Eigen;
60
61//=============================================================================================================
62// CONSTANTS
63//=============================================================================================================
64
65static constexpr double PI = M_PI;
66
67//=============================================================================================================
68// PRIVATE HELPERS
69//=============================================================================================================
70
71namespace {
72
73//=============================================================================================================
77double evalMagnitude(const QVector<IirBiquad>& sos, double omega)
78{
79 std::complex<double> z(std::cos(omega), std::sin(omega));
80 std::complex<double> H(1.0, 0.0);
81 std::complex<double> z1(1.0, 0.0); // z^{-1} placeholder
82
83 for (const IirBiquad& bq : sos) {
84 std::complex<double> zinv = 1.0 / z;
85 std::complex<double> zinv2 = zinv * zinv;
86 std::complex<double> num = bq.b0 + bq.b1 * zinv + bq.b2 * zinv2;
87 std::complex<double> den = 1.0 + bq.a1 * zinv + bq.a2 * zinv2;
88 H *= num / den;
89 }
90 return std::abs(H);
91}
92
93} // anonymous namespace
94
95//=============================================================================================================
96// MEMBER DEFINITIONS
97//=============================================================================================================
98
99QVector<std::complex<double>> IirFilter::butterworthPrototypePoles(int n)
100{
101 // Butterworth analogue LP prototype poles (normalised cutoff = 1 rad/s):
102 // p_k = exp( j * PI * (2k - 1 + n) / (2n) ) for k = 1 .. n
103 // All poles lie on the unit circle in the left half-plane.
104 QVector<std::complex<double>> poles;
105 poles.reserve(n);
106 for (int k = 1; k <= n; ++k) {
107 double angle = PI * static_cast<double>(2 * k - 1 + n) / static_cast<double>(2 * n);
108 poles.append(std::complex<double>(std::cos(angle), std::sin(angle)));
109 }
110 return poles;
111}
112
113//=============================================================================================================
114
115IirBiquad IirFilter::poleToDigitalBiquad(std::complex<double> pole,
116 double dC,
117 double dGain)
118{
119 // Given complex conjugate pole pair (pole, conj(pole)), the analogue biquad is:
120 // H_a(s) = gain * omega0^2 / (s^2 - 2*Re(pole)*s + |pole|^2)
121 // where omega0 = |pole|
122 //
123 // Bilinear transform s = dC*(z-1)/(z+1):
124 // Let omega0 = |pole|, alpha = -Re(pole) (> 0), c = dC
125 //
126 // H_a(s) numerator coefficients: b2_a=0, b1_a=0, b0_a = gain*omega0^2
127 // H_a(s) denominator: a2_a=1, a1_a = 2*alpha, a0_a = omega0^2
128 //
129 // After bilinear (general formula for b2_a=0, b1_a=0):
130 // d0 = c^2 + 2*alpha*c + omega0^2
131 // b0 = gain * omega0^2 / d0
132 // b1 = 2 * gain * omega0^2 / d0
133 // b2 = gain * omega0^2 / d0
134 // a1 = 2*(omega0^2 - c^2) / d0
135 // a2 = (c^2 - 2*alpha*c + omega0^2) / d0
136
137 double omega0 = std::abs(pole);
138 double alpha = -pole.real(); // positive for stable (left-half-plane) poles
139
140 double omega0sq = omega0 * omega0;
141 double d0 = dC * dC + 2.0 * alpha * dC + omega0sq;
142
143 IirBiquad bq;
144 bq.b0 = dGain * omega0sq / d0;
145 bq.b1 = 2.0 * dGain * omega0sq / d0;
146 bq.b2 = dGain * omega0sq / d0;
147 bq.a1 = 2.0 * (omega0sq - dC * dC) / d0;
148 bq.a2 = (dC * dC - 2.0 * alpha * dC + omega0sq) / d0;
149 return bq;
150}
151
152//=============================================================================================================
153
154IirBiquad IirFilter::realPoleToDigitalSection(double dPoleReal, double dC, double dGain)
155{
156 // H_a(s) = gain * |pole| / (s - dPoleReal) (dPoleReal < 0)
157 // Bilinear: s = dC*(z-1)/(z+1)
158 // H(z) = gain * |pole| * (1 + z^{-1}) / ((dC + |pole|) - (dC - |pole|)*z^{-1})
159 // = [gain*|pole|/(dC+|pole|)] * (1 + z^{-1}) / (1 - [(dC-|pole|)/(dC+|pole|)]*z^{-1})
160
161 double absPole = std::abs(dPoleReal);
162 double sum = dC + absPole;
163
164 IirBiquad bq;
165 bq.b0 = dGain * absPole / sum;
166 bq.b1 = dGain * absPole / sum;
167 bq.b2 = 0.0;
168 bq.a1 = -(dC - absPole) / sum;
169 bq.a2 = 0.0;
170 return bq;
171}
172
173//=============================================================================================================
174
175QVector<IirBiquad> IirFilter::designButterworth(int iOrder,
176 FilterType type,
177 double dCutoffLow,
178 double dCutoffHigh,
179 double dSFreq)
180{
181 if (iOrder < 1) {
182 qWarning() << "IirFilter::designButterworth: order must be >= 1, got" << iOrder;
183 return {};
184 }
185
186 const double dC = 2.0 * dSFreq; // bilinear pre-warp constant
187
188 // Pre-warp analogue cutoff frequencies
189 double dOmegaLow = dC * std::tan(PI * dCutoffLow / dSFreq);
190 double dOmegaHigh = (type == BandPass || type == BandStop)
191 ? dC * std::tan(PI * dCutoffHigh / dSFreq)
192 : 0.0;
193
194 // Butterworth normalised LP prototype poles
195 QVector<std::complex<double>> protoPoles = butterworthPrototypePoles(iOrder);
196
197 QVector<IirBiquad> sos;
198
199 if (type == LowPass || type == HighPass) {
200 // ---- LP / HP -------------------------------------------------------
201 // LP: scale pole by Omega_c → p_k_lp = Omega_c * proto_k
202 // HP: invert-and-scale → p_k_hp = Omega_c / proto_k
203 //
204 // Process complex conjugate pairs (take poles with Im >= 0).
205 // If order is odd there will be one real pole at proto = -1.
206
207 for (int k = 0; k < iOrder; ++k) {
208 std::complex<double> proto = protoPoles[k];
209
210 // Skip the conjugate (imaginary part < 0); it is handled together with its pair.
211 if (proto.imag() < -1e-10) {
212 continue;
213 }
214
215 if (std::abs(proto.imag()) < 1e-10) {
216 // Real pole (only for odd order, proto ≈ -1)
217 double analogPole = (type == LowPass)
218 ? dOmegaLow * proto.real()
219 : dOmegaLow / proto.real();
220
221 // Evaluate passband gain at DC (LP) or Nyquist (HP) for normalisation.
222 // For a single first-order section: |H(0)| = |b0+b1|/|1+a1|
223 // We design with dGain=1 first, then normalise.
224 IirBiquad bq = realPoleToDigitalSection(analogPole, dC, 1.0);
225
226 // Normalise so passband gain ~ 1
227 double w_pb = (type == LowPass) ? 0.0 : PI;
228 double num = std::abs(bq.b0 + bq.b1 * std::cos(-w_pb)); // rough, imaginary part 0
229 double den = std::abs(1.0 + bq.a1 * std::cos(-w_pb));
230 double gain = (den > 1e-12 && num > 1e-12) ? den / num : 1.0;
231 bq.b0 *= gain;
232 bq.b1 *= gain;
233
234 sos.append(bq);
235 } else {
236 // Complex conjugate pair
237 std::complex<double> analogPole;
238 double sectionGain;
239
240 if (type == LowPass) {
241 analogPole = dOmegaLow * proto;
242 sectionGain = dOmegaLow * dOmegaLow; // LP: constant numerator = omega0^2
243 } else {
244 // HP: s -> Omega_c / s transforms pole p to Omega_c / p
245 analogPole = dOmegaLow / proto;
246 sectionGain = 1.0; // HP: numerator has s^2 -> gain handled below
247 }
248
249 if (type == LowPass) {
250 IirBiquad bq = poleToDigitalBiquad(analogPole, dC, sectionGain);
251 // Normalise: LP passband gain at DC (z=1, omega=0) should be 1
252 // At z=1: H = (b0+b1+b2)/(1+a1+a2)
253 double hDC = (bq.b0 + bq.b1 + bq.b2) / (1.0 + bq.a1 + bq.a2);
254 if (std::abs(hDC) > 1e-12) {
255 double scale = 1.0 / std::abs(hDC);
256 bq.b0 *= scale; bq.b1 *= scale; bq.b2 *= scale;
257 }
258 sos.append(bq);
259 } else {
260 // HP biquad: analogue section H(s) = s^2 / (s^2 - 2*Re(p)*s + |p|^2)
261 // b2_a=1, b1_a=0, b0_a=0; a2_a=1, a1_a=-2*Re(p), a0_a=|p|^2
262 // General bilinear for (b2_a=1, b1_a=0, b0_a=0):
263 double omega0 = std::abs(analogPole);
264 double alpha = -analogPole.real();
265 double omega0sq = omega0 * omega0;
266 double d0 = dC * dC + 2.0 * alpha * dC + omega0sq;
267
268 IirBiquad bq;
269 bq.b0 = dC * dC / d0;
270 bq.b1 = -2.0 * dC * dC / d0;
271 bq.b2 = dC * dC / d0;
272 bq.a1 = 2.0 * (omega0sq - dC * dC) / d0;
273 bq.a2 = (dC * dC - 2.0 * alpha * dC + omega0sq) / d0;
274
275 // Normalise: HP passband gain at Nyquist (z=-1, omega=pi) should be 1
276 // At z=-1: H = (b0-b1+b2)/(1-a1+a2)
277 double hNy = (bq.b0 - bq.b1 + bq.b2) / (1.0 - bq.a1 + bq.a2);
278 if (std::abs(hNy) > 1e-12) {
279 double scale = 1.0 / std::abs(hNy);
280 bq.b0 *= scale; bq.b1 *= scale; bq.b2 *= scale;
281 }
282 sos.append(bq);
283 }
284 }
285 }
286
287 } else {
288 // ---- BP / BS -------------------------------------------------------
289 // LP prototype pole p_k → BP/BS analogue poles via:
290 // BP: s^2 - Bw*p_k*s + Omega0^2 = 0
291 // BS: Bw*s / (s^2 + Omega0^2) = p_k (i.e., solve s^2 - (Bw/p_k)*s + Omega0^2 = 0)
292 // This doubles the order, producing 2*iOrder poles → iOrder biquad sections.
293
294 double dOmega0 = std::sqrt(dOmegaLow * dOmegaHigh);
295 double dBw = dOmegaHigh - dOmegaLow;
296
297 for (int k = 0; k < iOrder; ++k) {
298 std::complex<double> proto = protoPoles[k];
299
300 // Solve quadratic for the two BP or BS analogue poles from this prototype pole.
301 // BP quadratic: s^2 - Bw*proto*s + Omega0^2 = 0
302 // BS quadratic: s^2 - (Bw/proto)*s + Omega0^2 = 0
303 std::complex<double> mid;
304 if (type == BandPass) {
305 mid = dBw * proto;
306 } else {
307 mid = dBw / proto;
308 }
309
310 std::complex<double> discriminant = mid * mid - 4.0 * dOmega0 * dOmega0;
311 std::complex<double> sqrtDisc = std::sqrt(discriminant);
312 std::complex<double> pole1 = (mid + sqrtDisc) / 2.0;
313 std::complex<double> pole2 = (mid - sqrtDisc) / 2.0;
314
315 // Each of pole1, pole2 forms a complex conjugate pair with its own conjugate.
316 // Design one biquad section per pole (using poleToDigitalBiquad which handles
317 // the conjugate pair automatically).
318
319 for (auto& pole : {pole1, pole2}) {
320 if (type == BandPass) {
321 // BP analogue biquad numerator: Bw*s (zero at origin, zero at infinity)
322 // i.e., b2_a=0, b1_a=Bw, b0_a=0
323 // Bilinear of (b2_a=0, b1_a=Bw, b0_a=0) / (s^2 - 2Re(p)s + |p|^2):
324 double omega0 = std::abs(pole);
325 double alpha = -pole.real();
326 double omega0sq = omega0 * omega0;
327 double d0 = dC * dC + 2.0 * alpha * dC + omega0sq;
328
329 IirBiquad bq;
330 // num from bilinear(b1_a=Bw): b0=Bw*c, b1=0, b2=-Bw*c
331 bq.b0 = dBw * dC / d0;
332 bq.b1 = 0.0;
333 bq.b2 = -dBw * dC / d0;
334 bq.a1 = 2.0 * (omega0sq - dC * dC) / d0;
335 bq.a2 = (dC * dC - 2.0 * alpha * dC + omega0sq) / d0;
336 sos.append(bq);
337 } else {
338 // BS analogue biquad numerator: s^2 + Omega0^2 (zeros at ±j*Omega0)
339 // b2_a=1, b1_a=0, b0_a=Omega0^2
340 // Bilinear:
341 double omega0 = std::abs(pole);
342 double alpha = -pole.real();
343 double omega0sq = omega0 * omega0;
344 double dOmega0sq = dOmega0 * dOmega0;
345 double d0 = dC * dC + 2.0 * alpha * dC + omega0sq;
346
347 IirBiquad bq;
348 // num from bilinear(b2_a=1, b1_a=0, b0_a=Omega0^2):
349 // b0 = (c^2 + Omega0^2)/d0, b1 = 2*(Omega0^2-c^2)/d0, b2 = (c^2+Omega0^2)/d0
350 bq.b0 = (dC * dC + dOmega0sq) / d0;
351 bq.b1 = 2.0 * (dOmega0sq - dC * dC) / d0;
352 bq.b2 = (dC * dC + dOmega0sq) / d0;
353 bq.a1 = 2.0 * (omega0sq - dC * dC) / d0;
354 bq.a2 = (dC * dC - 2.0 * alpha * dC + omega0sq) / d0;
355 sos.append(bq);
356 }
357 }
358 }
359
360 // Normalise overall passband gain
361 // For BP: evaluate at centre frequency Omega0 (digital: omega_0 = Omega0 / fs * 2pi... use prewarped)
362 // For BS: evaluate at DC (omega = 0)
363 double omegaCheck = (type == BandPass)
364 ? 2.0 * std::atan(dOmega0 / dC) // bilinear inverse: omega_d = 2*atan(Omega_a/c)
365 : 0.0;
366 double totalGain = evalMagnitude(sos, omegaCheck);
367 if (totalGain > 1e-12) {
368 double scale = 1.0 / totalGain;
369 // Distribute scale evenly across sections
370 double perSection = std::pow(scale, 1.0 / sos.size());
371 for (IirBiquad& bq : sos) {
372 bq.b0 *= perSection;
373 bq.b1 *= perSection;
374 bq.b2 *= perSection;
375 }
376 }
377 }
378
379 return sos;
380}
381
382//=============================================================================================================
383
384RowVectorXd IirFilter::applySos(const RowVectorXd& vecData, const QVector<IirBiquad>& sos)
385{
386 if (sos.isEmpty() || vecData.size() == 0) {
387 return vecData;
388 }
389
390 RowVectorXd y = vecData;
391
392 // Cascade all biquad sections using Direct-Form II transposed structure:
393 // w[n] = x[n] - a1*w[n-1] - a2*w[n-2]
394 // y[n] = b0*w[n] + b1*w[n-1] + b2*w[n-2]
395 // Equivalent transposed: uses two state variables (s1, s2)
396
397 for (const IirBiquad& bq : sos) {
398 RowVectorXd x = y;
399 double s1 = 0.0, s2 = 0.0;
400
401 for (int i = 0; i < static_cast<int>(x.size()); ++i) {
402 double xi = x(i);
403 double yi = bq.b0 * xi + s1;
404 s1 = bq.b1 * xi - bq.a1 * yi + s2;
405 s2 = bq.b2 * xi - bq.a2 * yi;
406 y(i) = yi;
407 }
408 }
409
410 return y;
411}
412
413//=============================================================================================================
414
415RowVectorXd IirFilter::applyZeroPhase(const RowVectorXd& vecData, const QVector<IirBiquad>& sos)
416{
417 if (sos.isEmpty() || vecData.size() == 0) {
418 return vecData;
419 }
420
421 // Forward pass
422 RowVectorXd y = applySos(vecData, sos);
423
424 // Reverse
425 RowVectorXd yRev = y.reverse();
426
427 // Backward pass
428 RowVectorXd yRevFilt = applySos(yRev, sos);
429
430 // Reverse again
431 return yRevFilt.reverse();
432}
433
434//=============================================================================================================
435
436MatrixXd IirFilter::applyZeroPhaseMatrix(const MatrixXd& matData, const QVector<IirBiquad>& sos)
437{
438 if (sos.isEmpty() || matData.size() == 0) {
439 return matData;
440 }
441
442 MatrixXd matOut(matData.rows(), matData.cols());
443 for (int i = 0; i < static_cast<int>(matData.rows()); ++i) {
444 matOut.row(i) = applyZeroPhase(matData.row(i), sos);
445 }
446 return matOut;
447}
#define M_PI
Declaration of IirFilter — Butterworth IIR filter design and application using numerically stable sec...
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Coefficients of one second-order IIR section (biquad).
Definition iirfilter.h:78
static Eigen::RowVectorXd applyZeroPhase(const Eigen::RowVectorXd &vecData, const QVector< IirBiquad > &sos)
static QVector< IirBiquad > designButterworth(int iOrder, FilterType type, double dCutoffLow, double dCutoffHigh, double dSFreq)
static Eigen::MatrixXd applyZeroPhaseMatrix(const Eigen::MatrixXd &matData, const QVector< IirBiquad > &sos)
FilterType
Filter type selector.
Definition iirfilter.h:113
static Eigen::RowVectorXd applySos(const Eigen::RowVectorXd &vecData, const QVector< IirBiquad > &sos)