kfr

Fast, modern C++ DSP framework, FFT, Sample Rate Conversion, FIR/IIR/Biquad Filters (SSE, AVX, AVX-512, ARM NEON)
Log | Files | Refs | README

commit 5b2a10aac1cbc7c2c693e91acbd6e985e287440c
parent ebd5d10dd7d7954cf69403344da10c7747ce1f7d
Author: d.levin256@gmail.com <d.levin256@gmail.com>
Date:   Tue, 23 Jan 2024 08:37:43 +0000

iir_design: explicit casts

Diffstat:
Minclude/kfr/dsp/iir_design.hpp | 44++++++++++++++++++++++----------------------
1 file changed, 22 insertions(+), 22 deletions(-)

diff --git a/include/kfr/dsp/iir_design.hpp b/include/kfr/dsp/iir_design.hpp @@ -57,8 +57,8 @@ KFR_FUNCTION zpk<T> chebyshev1(int N, identity<T> rp) return { {}, {}, 1 }; } - T eps = sqrt(exp10(0.1 * rp) - 1.0); - T mu = 1.0 / N * std::asinh(1 / eps); + T eps = sqrt(exp10(T(0.1) * rp) - T(1.0)); + T mu = T(1.0) / N * std::asinh(1 / eps); univector<T> m = linspace(-N + 1, N + 1, N, false, ctrue); @@ -66,7 +66,7 @@ KFR_FUNCTION zpk<T> chebyshev1(int N, identity<T> rp) univector<complex<T>> p = -csinh(make_complex(mu, theta)); T k = product(-p).real(); if (N % 2 == 0) - k = k / sqrt(1.0 + eps * eps); + k = k / sqrt(T(1.0) + eps * eps); return { {}, std::move(p), k }; } @@ -78,8 +78,8 @@ KFR_FUNCTION zpk<T> chebyshev2(int N, identity<T> rs) return { {}, {}, 1 }; } - T de = 1.0 / sqrt(exp10(0.1 * rs) - 1); - T mu = std::asinh(1.0 / de) / N; + T de = T(1.0) / sqrt(exp10(T(0.1) * rs) - 1); + T mu = std::asinh(T(1.0) / de) / N; univector<T> m; @@ -92,12 +92,12 @@ KFR_FUNCTION zpk<T> chebyshev2(int N, identity<T> rs) m = linspace(-N + 1, N + 1, N, false, ctrue); } - univector<complex<T>> z = -cconj(complex<T>(0, 1) / sin(m * c_pi<T> / (2.0 * N))); + univector<complex<T>> z = -cconj(complex<T>(0, 1) / sin(m * c_pi<T> / (T(2.0) * N))); univector<complex<T>> p = -cexp(complex<T>(0, 1) * c_pi<T> * linspace(-N + 1, N + 1, N, false, ctrue) / (2 * N)); p = make_complex(sinh(mu) * real(p), cosh(mu) * imag(p)); - p = 1.0 / p; + p = T(1.0) / p; T k = (product(-p) / product(-z)).real(); @@ -864,7 +864,7 @@ namespace internal template <typename T> KFR_FUNCTION zpk<T> bilinear(const zpk<T>& filter, identity<T> fs) { - const T fs2 = 2.0 * fs; + const T fs2 = T(2.0) * fs; zpk<T> result; result.z = (fs2 + filter.z) / (fs2 - filter.z); result.p = (fs2 + filter.p) / (fs2 - filter.p); @@ -984,8 +984,8 @@ template <typename T> KFR_FUNCTION zpk<T> lp2bp_zpk(const zpk<T>& filter, identity<T> wo, identity<T> bw) { zpk<T> lowpass; - lowpass.z = bw * 0.5 * filter.z; - lowpass.p = bw * 0.5 * filter.p; + lowpass.z = bw * T(0.5) * filter.z; + lowpass.p = bw * T(0.5) * filter.p; zpk<T> result; result.z = concatenate(lowpass.z + csqrt(csqr(lowpass.z) - wo * wo), @@ -1003,8 +1003,8 @@ template <typename T> KFR_FUNCTION zpk<T> lp2bs_zpk(const zpk<T>& filter, identity<T> wo, identity<T> bw) { zpk<T> highpass; - highpass.z = (bw * 0.5) / filter.z; - highpass.p = (bw * 0.5) / filter.p; + highpass.z = (bw * T(0.5)) / filter.z; + highpass.p = (bw * T(0.5)) / filter.p; zpk<T> result; result.z = concatenate(highpass.z + csqrt(csqr(highpass.z) - wo * wo), @@ -1023,7 +1023,7 @@ template <typename T> KFR_FUNCTION T warp_freq(T frequency, T fs) { frequency = 2 * frequency / fs; - fs = 2.0; + fs = T(2.0); T warped = 2 * fs * tan(c_pi<T> * frequency / fs); return warped; } @@ -1031,50 +1031,50 @@ KFR_FUNCTION T warp_freq(T frequency, T fs) } // namespace internal template <typename T> -KFR_FUNCTION zpk<T> iir_lowpass(const zpk<T>& filter, identity<T> frequency, identity<T> fs = 2.0) +KFR_FUNCTION zpk<T> iir_lowpass(const zpk<T>& filter, identity<T> frequency, identity<T> fs = T(2.0)) { T warped = internal::warp_freq(frequency, fs); zpk<T> result = filter; result = internal::lp2lp_zpk(result, warped); - result = internal::bilinear(result, 2.0); // fs = 2.0 + result = internal::bilinear(result, T(2.0)); // fs = 2.0 return result; } template <typename T> -KFR_FUNCTION zpk<T> iir_highpass(const zpk<T>& filter, identity<T> frequency, identity<T> fs = 2.0) +KFR_FUNCTION zpk<T> iir_highpass(const zpk<T>& filter, identity<T> frequency, identity<T> fs = T(2.0)) { T warped = internal::warp_freq(frequency, fs); zpk<T> result = filter; result = internal::lp2hp_zpk(result, warped); - result = internal::bilinear(result, 2.0); // fs = 2.0 + result = internal::bilinear(result, T(2.0)); // fs = 2.0 return result; } template <typename T> KFR_FUNCTION zpk<T> iir_bandpass(const zpk<T>& filter, identity<T> lowfreq, identity<T> highfreq, - identity<T> fs = 2.0) + identity<T> fs = T(2.0)) { T warpedlow = internal::warp_freq(lowfreq, fs); T warpedhigh = internal::warp_freq(highfreq, fs); zpk<T> result = filter; result = internal::lp2bp_zpk(result, std::sqrt(warpedlow * warpedhigh), warpedhigh - warpedlow); - result = internal::bilinear(result, 2.0); // fs = 2.0 + result = internal::bilinear(result, T(2.0)); // fs = 2.0 return result; } template <typename T> KFR_FUNCTION zpk<T> iir_bandstop(const zpk<T>& filter, identity<T> lowfreq, identity<T> highfreq, - identity<T> fs = 2.0) + identity<T> fs = T(2.0)) { T warpedlow = internal::warp_freq(lowfreq, fs); T warpedhigh = internal::warp_freq(highfreq, fs); zpk<T> result = filter; result = internal::lp2bs_zpk(result, std::sqrt(warpedlow * warpedhigh), warpedhigh - warpedlow); - result = internal::bilinear(result, 2.0); // fs = 2.0 + result = internal::bilinear(result, T(2.0)); // fs = 2.0 return result; } @@ -1082,7 +1082,7 @@ template <typename T> KFR_FUNCTION std::vector<biquad_params<T>> to_sos(const zpk<T>& filter) { if (filter.p.empty() && filter.z.empty()) - return { biquad_params<T>(filter.k, 0., 0., 1., 0., 0) }; + return { biquad_params<T>(filter.k, T(0.), T(0.), T(1.), T(0.), 0) }; zpk<T> filt = filter; size_t length = std::max(filter.p.size(), filter.z.size());