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 56946290183e636376b32d514974cbe6226d0c12
parent 5c37b7b76602b59c7c5b4891d3689aeb3a952ac7
Author: [email protected] <[email protected]>
Date:   Fri,  9 Aug 2024 16:20:35 +0100

More detailed comments for examples

Diffstat:
Mexamples/biquads.cpp | 134++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++-------------
Mexamples/dft.cpp | 23++++++++++++++++-------
Mexamples/fir.cpp | 75++++++++++++++++++++++++++++++++++++++++++---------------------------------
Mexamples/iir.cpp | 159++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++---------------
Mexamples/sample_rate_conversion.cpp | 49++++++++++++++++++++++++++++++++++++++++++++++---
Mexamples/window.cpp | 19+++++++++++++++++++
6 files changed, 365 insertions(+), 94 deletions(-)

diff --git a/examples/biquads.cpp b/examples/biquads.cpp @@ -14,101 +14,193 @@ using namespace kfr; int main() { + // Print the version of the KFR library being used println(library_version()); + // Define options for plotting DSP data const std::string options = "phaseresp=True"; + // Define a buffer for storing the filter output univector<fbase, 128> output; + + // -------------------------------------------------------------------------------------- + // ------------------------- Biquad Notch Filters Example ------------------------------- + // -------------------------------------------------------------------------------------- { - biquad_section<fbase> bq[] = { biquad_notch(0.1, 0.5), biquad_notch(0.2, 0.5), biquad_notch(0.3, 0.5), - biquad_notch(0.4, 0.5) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Initialize an array of biquad notch filters with different center frequencies + // biquad_notch(frequency, Q) where the frequency is relative to the sample rate + biquad_section<fbase> bq[] = { + biquad_notch(0.1, 0.5), + biquad_notch(0.2, 0.5), + biquad_notch(0.3, 0.5), + biquad_notch(0.4, 0.5), + }; + // Apply the biquad filters to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter responses plot_save("biquad_notch", output, options + ", title='Four Biquad Notch filters'"); + // -------------------------------------------------------------------------------------- + // ------------------------- Biquad Lowpass Filter Example ------------------------------ + // -------------------------------------------------------------------------------------- { + // Initialize a biquad lowpass filter with specific parameters + // biquad_lowpass(frequency, Q) where the frequency is relative to the sample rate biquad_section<fbase> bq[] = { biquad_lowpass(0.2, 0.9) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad lowpass filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_lowpass", output, options + ", title='Biquad Low pass filter (0.2, 0.9)'"); + // -------------------------------------------------------------------------------------- + // ------------------------- Biquad Highpass Filter Example ----------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad highpass filter with specific parameters + // biquad_highpass(frequency, Q) where the frequency is relative to the sample rate biquad_section<fbase> bq[] = { biquad_highpass(0.3, 0.1) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad highpass filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_highpass", output, options + ", title='Biquad High pass filter (0.3, 0.1)'"); + // -------------------------------------------------------------------------------------- + // -------------------------- Biquad Peak Filter Example -------------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad peak filter with specific parameters + // biquad_peak(frequency, Q, gain) where the frequency is relative to the sample rate and the gain is + // in decibels biquad_section<fbase> bq[] = { biquad_peak(0.3, 0.5, +9.0) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad peak filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_peak", output, options + ", title='Biquad Peak filter (0.2, 0.5, +9)'"); + // -------------------------------------------------------------------------------------- + // -------------------------- Biquad Peak Filter Example (2) ---------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize another biquad peak filter with different parameters + // biquad_peak(frequency, Q, gain) where the frequency is relative to the sample rate and the gain is + // in decibels biquad_section<fbase> bq[] = { biquad_peak(0.3, 3.0, -2.0) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad peak filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_peak2", output, options + ", title='Biquad Peak filter (0.3, 3, -2)'"); + // -------------------------------------------------------------------------------------- + // -------------------------- Biquad Low Shelf Filter Example --------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad low shelf filter with specific parameters + // biquad_lowshelf(frequency, gain) where the frequency is relative to the sample rate and the gain is + // in decibels biquad_section<fbase> bq[] = { biquad_lowshelf(0.3, -1.0) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad low shelf filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_lowshelf", output, options + ", title='Biquad low shelf filter (0.3, -1)'"); + // -------------------------------------------------------------------------------------- + // -------------------------- Biquad High Shelf Filter Example -------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad high shelf filter with specific parameters + // biquad_highshelf(frequency, gain) where the frequency is relative to the sample rate and the gain + // is in decibels biquad_section<fbase> bq[] = { biquad_highshelf(0.3, +9.0) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad high shelf filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_highshelf", output, options + ", title='Biquad high shelf filter (0.3, +9)'"); + // -------------------------------------------------------------------------------------- + // ------------------------- Biquad Bandpass Filter Example ----------------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad bandpass filter with specific parameters + // biquad_bandpass(frequency, Q) where the frequency is relative to the sample rate biquad_section<fbase> bq[] = { biquad_bandpass(0.25, 0.2) }; - output = iir(unitimpulse(), iir_params{ bq }); + // Apply the biquad bandpass filter to a unit impulse signal and store the result in 'output' + output = iir(unitimpulse(), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_bandpass", output, options + ", title='Biquad band pass (0.25, 0.2)'"); + // -------------------------------------------------------------------------------------- + // ----------------- Biquad Bandpass Filter Example with std::vector -------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad bandpass filter with specific parameters biquad_section<fbase> bq[] = { biquad_bandpass(0.25, 0.2) }; + // Create a std::vector for the input data and set the first element to 1 (unit impulse) std::vector<fbase> data(output.size(), 0.f); data[0] = 1.f; - output = iir(make_univector(data), iir_params{ bq }); + // Apply the biquad bandpass filter to the input data and store the result in 'output' + output = iir(make_univector(data), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_bandpass_stdvector", output, options + ", title='Biquad band pass (0.25, 0.2)'"); + // -------------------------------------------------------------------------------------- + // ------------------- Biquad Bandpass Filter Example with C array ---------------------- + // -------------------------------------------------------------------------------------- { + // Initialize a biquad bandpass filter with specific parameters biquad_section<fbase> bq[] = { biquad_bandpass(0.25, 0.2) }; - fbase data[output.size()] = { 0 }; // .size() is constexpr - data[0] = 1.f; - output = iir(make_univector(data), iir_params{ bq }); + // Create a C array for the input data and set the first element to 1 (unit impulse) + fbase data[output.size()] = { 0 }; // .size() is constexpr + data[0] = 1.f; + // Apply the biquad bandpass filter to the input data and store the result in 'output' + output = iir(make_univector(data), iir_params{ bq }); } + // Save the plot of the filter response plot_save("biquad_bandpass_carray", output, options + ", title='Biquad band pass (0.25, 0.2)'"); + // -------------------------------------------------------------------------------------- + // ----------------- Custom Biquad Lowpass Filter using expression_filter --------------- + // -------------------------------------------------------------------------------------- { - // filter initialization - biquad_section<fbase> bq[] = { biquad_lowpass(0.2, 0.9) }; + // Initialize a biquad lowpass filter with specific parameters + biquad_section<fbase> bq[] = { biquad_lowpass(0.2, 0.9) }; + // Create a type-erased expression filter for the biquad lowpass filter expression_filter<fbase> filter = to_filter(iir(placeholder<fbase>(), iir_params{ bq })); - // prepare data + // Prepare a unit impulse signal output = unitimpulse(); - // apply filter + // Apply the expression filter to the unit impulse signal filter.apply(output); } + // Save the plot of the filter response plot_save("biquad_custom_filter_lowpass", output, options + ", title='Biquad Low pass filter (0.2, 0.9) (using expression_filter)'"); + // -------------------------------------------------------------------------------------- + // ---------------------- Custom Biquad Lowpass Filter using iir_filter ----------------- + // -------------------------------------------------------------------------------------- { - // filter initialization + // Initialize a biquad lowpass filter with specific parameters biquad_section<fbase> bq[] = { biquad_lowpass(0.2, 0.9) }; + // Create an IIR filter for the biquad lowpass filter iir_filter<fbase> filter(bq); - // prepare data + // Prepare a unit impulse signal output = unitimpulse(); - // apply filter + // Apply the IIR filter to the unit impulse signal filter.apply(output); } + // Save the plot of the filter response plot_save("biquad_filter_lowpass", output, options + ", title='Biquad Low pass filter (0.2, 0.9) (using iir_filter)'"); diff --git a/examples/dft.cpp b/examples/dft.cpp @@ -13,39 +13,48 @@ using namespace kfr; int main() { + // Print the version of the KFR library being used println(library_version()); - // fft size + // Define the size of the Fast Fourier Transform (FFT) const size_t size = 128; - // initialize input & output buffers + // Initialize input and output buffers + // 'in' buffer is filled with a sine wave spanning 4 cycles over the given range + // 'out' buffer is initialized with 'qnan' (quiet NaN) to represent uninitialized state univector<complex<fbase>, size> in = sin(linspace(0.0, c_pi<fbase, 2> * 4.0, size)); univector<complex<fbase>, size> out = scalar(qnan); - // initialize fft + // Create an FFT plan for the defined size const dft_plan<fbase> dft(size); + // Dump the details of the FFT plan (for debugging or information purposes) dft.dump(); - // allocate work buffer for fft (if needed) + // Allocate a temporary buffer for the FFT computation if needed univector<u8> temp(dft.temp_size); - // perform forward fft + // Perform the forward FFT on the input buffer 'in', store the result in the 'out' buffer dft.execute(out, in, temp); - // scale output + // Scale the output of the FFT by dividing by the size out = out / size; - // get magnitude and convert to decibels + // Convert the amplitude of the FFT output to decibels (dB) + // 'cabs' computes the magnitude of the complex numbers in 'out' + // 'amp_to_dB' converts the amplitude values to decibel scale univector<fbase, size> dB = amp_to_dB(cabs(out)); + // Print the maximum, minimum, mean, and root mean square (RMS) of the dB values println("max = ", maxof(dB)); println("min = ", minof(dB)); println("mean = ", mean(dB)); println("rms = ", rms(dB)); + // Print the input buffer 'in' println(in); println(); + // Print the dB values println(dB); return 0; } diff --git a/examples/fir.cpp b/examples/fir.cpp @@ -13,84 +13,91 @@ using namespace kfr; +// Define a macro to check if Python is installed #ifndef PYTHON_IS_INSTALLED #define PYTHON_IS_INSTALLED 1 #endif int main() { - // Show KFR version + // Print the version of the KFR library being used println(library_version()); // -------------------------------------------------------------------------------------- // --------------------- FIR filter design (using window functions) --------------------- // -------------------------------------------------------------------------------------- - // Options for dspplot + // Define options for plotting DSP data const std::string options = "phaseresp=True"; + // Define filter coefficient buffers of different sizes univector<fbase, 15> taps15; univector<fbase, 127> taps127; univector<fbase, 8191> taps8191; - // Prepare window functions (only expression saved here, not data) + // Prepare window function expressions, data is not processed at this point + // Hann window for 15-point filter expression_handle<fbase> hann = to_handle(window_hann(taps15.size())); - + // Kaiser window with alpha = 3.0 for 127-point filter expression_handle<fbase> kaiser = to_handle(window_kaiser(taps127.size(), 3.0)); - + // Blackman-Harris window for 8191-point filter expression_handle<fbase> blackman_harris = to_handle(window_blackman_harris(taps8191.size())); - // Fill taps15 with the low pass FIR filter coefficients using hann window and cutoff=0.15 + // Design a 15-point lowpass FIR filter with a cutoff frequency of 0.15 using Hann window fir_lowpass(taps15, 0.15, hann, true); #if PYTHON_IS_INSTALLED - // Plot filter, frequency and impulse response - // plot_save calls python (matplotlib and numpy must be installed) and saves SVG file + // Plot the filter, frequency, and impulse response + // internally plot_save calls Python (matplotlib and numpy packages used) to save SVG files plot_save("fir_lowpass_hann", taps15, options + ", phasearg='auto', title='15-point lowpass FIR, Hann window'"); #endif - // Fill taps127 with the low pass FIR filter coefficients using kaiser window and cutoff=0.2 + // Design a 127-point lowpass FIR filter with a cutoff frequency of 0.2 using the Kaiser window + // The resulting coefficients are in taps127 fir_lowpass(taps127, 0.2, kaiser, true); #if PYTHON_IS_INSTALLED - // Plot filter, frequency and impulse response + // Plot the filter, frequency, and impulse response plot_save("fir_lowpass_kaiser", taps127, - options + ", phasearg='auto', title=r'127-point lowpass FIR, Kaiser window ($\\alpha=3.0$)'"); + options + ", phasearg='auto', title='127-point lowpass FIR, Kaiser window (\\alpha=3.0)'"); #endif - // Fill taps127 with the high pass FIR filter coefficients using kaiser window and cutoff=0.2 + // Design a 127-point highpass FIR filter with a cutoff frequency of 0.2 using the Kaiser window + // The resulting coefficients are in taps127 fir_highpass(taps127, 0.2, kaiser, true); #if PYTHON_IS_INSTALLED - // Plot filter, frequency and impulse response + // Plot the filter, frequency, and impulse response plot_save("fir_highpass_kaiser", taps127, - options + ", phasearg='auto', title=r'127-point highpass FIR, Kaiser window ($\\alpha=3.0$)'"); + options + ", phasearg='auto', title='127-point highpass FIR, Kaiser window (\\alpha=3.0)'"); #endif - // Fill taps127 with the band pass FIR filter coefficients using kaiser window and cutoff=0.2 and 0.4 + // Design a 127-point bandpass FIR filter with cutoff frequencies of 0.2 and 0.4 using the Kaiser window fir_bandpass(taps127, 0.2, 0.4, kaiser, true); #if PYTHON_IS_INSTALLED - // Plot filter, frequency and impulse response + // Plot the filter, frequency, and impulse response plot_save("fir_bandpass_kaiser", taps127, - options + ", phasearg='auto', title=r'127-point bandpass FIR, Kaiser window ($\\alpha=3.0$)'"); + options + ", phasearg='auto', title='127-point bandpass FIR, Kaiser window (\\alpha=3.0)'"); #endif - // Fill taps127 with the band stop FIR filter coefficients using kaiser window and cutoff=0.2 and 0.4 + // Design a 127-point bandstop FIR filter with cutoff frequencies of 0.2 and 0.4 using the Kaiser window fir_bandstop(taps127, 0.2, 0.4, kaiser, true); #if PYTHON_IS_INSTALLED - // Show filter, frequency and impulse response + // Plot the filter, frequency, and impulse response plot_save("fir_bandstop_kaiser", taps127, - options + ", phasearg='auto', title=r'127-point bandstop FIR, Kaiser window ($\\alpha=3.0$)'"); + options + ", phasearg='auto', title='127-point bandstop FIR, Kaiser window (\\alpha=3.0)'"); #endif - // Fill taps8191 with the low pass FIR filter coefficients using blackman harris window and cutoff=0.15 + // Design an 8191-point lowpass FIR filter with a cutoff frequency of 0.15 using the Blackman-Harris + // window fir_lowpass(taps8191, 0.15, blackman_harris, true); + // Initialize a buffer that is 8191+150 samples long and set it to zero univector<fbase, 8191 + 150> taps8191_150 = scalar(0); - // Shift by 150 samples + // Shift the filter coefficients by 150 samples taps8191_150.slice(150) = taps8191; #if PYTHON_IS_INSTALLED - // Plot filter, frequency and impulse response, pass phasearg to get correct phase shift (phasearg=offset - // to unit impulse in samples) + // Plot the filter, frequency, and impulse response + // phasearg is used to get the correct phase shift (offset by 150 samples) plot_save( "fir_lowpass_blackman", taps8191_150, options + @@ -101,14 +108,14 @@ int main() // -------------------------- Using FIR filter as an expression ------------------------- // -------------------------------------------------------------------------------------- - // Prepare 10000 samples of white noise + // Generate 10,000 samples of white noise univector<float> noise = truncate(gen_random_range(random_init(1, 2, 3, 4), -1.f, +1.f), 10000); - // Apply band stop filter + // Apply the bandstop filter designed earlier to the noise univector<float> filtered_noise = fir(noise, fir_params{ taps127 }); #if PYTHON_IS_INSTALLED - // Plot results + // Plot the original noise and the filtered noise plot_save("noise", noise, "title='Original noise', div_by_N=True"); plot_save("filtered_noise", filtered_noise, "title='Filtered noise', div_by_N=True"); #endif @@ -117,16 +124,18 @@ int main() // ------------------------------- FIR filter as a class -------------------------------- // -------------------------------------------------------------------------------------- + // Redesign the 127-point bandpass FIR filter with cutoff frequencies of 0.2 and 0.4 using the Kaiser + // window fir_bandpass(taps127, 0.2, 0.4, kaiser, true); - // Initialize FIR filter with float input/output and fbase taps + // Initialize an FIR filter class with float input/output and fbase-typed taps filter_fir<fbase, float> fir_filter(taps127); - // Apply to univector, static array, data by pointer or anything + // Apply the FIR filter to the noise data univector<float> filtered_noise2; fir_filter.apply(filtered_noise2, noise); #if PYTHON_IS_INSTALLED - // Plot results + // Plot the results of the filtered noise plot_save("filtered_noise2", filtered_noise2, "title='Filtered noise 2', div_by_N=True"); #endif @@ -135,15 +144,15 @@ int main() // ---------------------- Convolution filter (optimized using DFT) ---------------------- // -------------------------------------------------------------------------------------- - // Initialize FIR filter with float input/output and fbase taps + // Initialize a convolution filter with float input/output and fbase-typed taps convolve_filter<fbase> conv_filter(taps127); - // Apply to univector, static array, data by pointer or anything + // Apply the convolution filter to the noise data univector<fbase> filtered_noise3; conv_filter.apply(filtered_noise3, noise * fbase(1.0)); #if PYTHON_IS_INSTALLED - // Plot results, same as filtered_noise2 + // Plot the results, same as filtered_noise2 plot_save("filtered_noise3", filtered_noise3, "title='Filtered noise 3', div_by_N=True"); #endif #endif diff --git a/examples/iir.cpp b/examples/iir.cpp @@ -12,86 +12,185 @@ using namespace kfr; int main() { + // Print the version of the KFR library being used println(library_version()); + // Define options for plotting DSP data const std::string options = "phaseresp=True, log_freq=True, freq_dB_lim=(-160, 10), padwidth=8192"; + // Define an output univector with 1024 elements univector<fbase, 1024> output; + + // -------------------------------------------------------------------------------------- + // ---------------------------- 24th-Order Bessel Lowpass Filter ------------------------ + // -------------------------------------------------------------------------------------- { + // Create a 24th-order Bessel lowpass filter with a cutoff frequency of 1 kHz and a sample rate of + // 48 kHz zpk<fbase> filt = iir_lowpass(bessel<fbase>(24), 1000, 48000); - output = iir(unitimpulse(), filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } - plot_save("bessel_lowpass24", output, options + ", title='24th-order Bessel filter, lowpass 1khz'"); + plot_save("bessel_lowpass24", output, options + ", title='24th-order Bessel filter, lowpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ---------------------------- 12th-Order Bessel Lowpass Filter ------------------------ + // -------------------------------------------------------------------------------------- { + // Create a 12th-order Bessel lowpass filter with a cutoff frequency of 1 kHz and a sample rate of + // 48 kHz zpk<fbase> filt = iir_lowpass(bessel<fbase>(12), 1000, 48000); - output = iir(unitimpulse(), filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } - plot_save("bessel_lowpass12", output, options + ", title='12th-order Bessel filter, lowpass 1khz'"); + plot_save("bessel_lowpass12", output, options + ", title='12th-order Bessel filter, lowpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ----------------------------- 6th-Order Bessel Lowpass Filter ------------------------ + // -------------------------------------------------------------------------------------- { + // Create a 6th-order Bessel lowpass filter with a cutoff frequency of 1 kHz and a sample rate of + // 48 kHz zpk<fbase> filt = iir_lowpass(bessel<fbase>(6), 1000, 48000); - output = iir(unitimpulse(), filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } - plot_save("bessel_lowpass6", output, options + ", title='6th-order Bessel filter, lowpass 1khz'"); + plot_save("bessel_lowpass6", output, options + ", title='6th-order Bessel filter, lowpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ------------------------ 24th-Order Butterworth Lowpass Filter ----------------------- + // -------------------------------------------------------------------------------------- { + // Create a 24th-order Butterworth lowpass filter with a cutoff frequency of 1 kHz and a sample rate + // of 48 kHz zpk<fbase> filt = iir_lowpass(butterworth<fbase>(24), 1000, 48000); - output = iir(unitimpulse(), filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } plot_save("butterworth_lowpass24", output, - options + ", title='24th-order Butterworth filter, lowpass 1khz'"); + options + ", title='24th-order Butterworth filter, lowpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ------------------------ 12th-Order Butterworth Lowpass Filter ----------------------- + // -------------------------------------------------------------------------------------- { + // Create a 12th-order Butterworth lowpass filter with a cutoff frequency of 1 kHz and a sample rate + // of 48 kHz zpk<fbase> filt = iir_lowpass(butterworth<fbase>(12), 1000, 48000); - output = iir(unitimpulse(), filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } plot_save("butterworth_lowpass12", output, - options + ", title='12th-order Butterworth filter, lowpass 1khz'"); + options + ", title='12th-order Butterworth filter, lowpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ------------------------ 12th-Order Butterworth Highpass Filter ---------------------- + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_highpass(butterworth<fbase>(12), 1000, 48000); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), bqs); + // Create a 12th-order Butterworth highpass filter with a cutoff frequency of 1 kHz and a sample rate + // of 48 kHz + zpk<fbase> filt = iir_highpass(butterworth<fbase>(12), 1000, 48000); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), bqs); } plot_save("butterworth_highpass12", output, - options + ", title='12th-order Butterworth filter, highpass 1khz'"); + options + ", title='12th-order Butterworth filter, highpass 1 kHz'"); + // -------------------------------------------------------------------------------------- + // ---------------------- 12th-Order Butterworth Bandpass Filter ------------------------ + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_bandpass(butterworth<fbase>(12), 0.1, 0.2); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), bqs); + // Create a 12th-order Butterworth bandpass filter with a passband from 0.1 to 0.2 (normalized + // frequency) + zpk<fbase> filt = iir_bandpass(butterworth<fbase>(12), 0.1, 0.2); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), bqs); } plot_save("butterworth_bandpass12", output, options + ", title='12th-order Butterworth filter, bandpass'"); + // -------------------------------------------------------------------------------------- + // ---------------------- 12th-Order Butterworth Bandstop Filter ------------------------ + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_bandstop(butterworth<fbase>(12), 0.1, 0.2); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), bqs); + // Create a 12th-order Butterworth bandstop filter with a stopband from 0.1 to 0.2 (normalized + // frequency) + zpk<fbase> filt = iir_bandstop(butterworth<fbase>(12), 0.1, 0.2); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), bqs); } plot_save("butterworth_bandstop12", output, options + ", title='12th-order Butterworth filter, bandstop'"); + // -------------------------------------------------------------------------------------- + // ------------------------ 4th-Order Butterworth Bandpass Filter ----------------------- + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_bandpass(butterworth<fbase>(4), 0.005, 0.9); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), bqs); + // Create a 4th-order Butterworth bandpass filter with a passband from 0.005 to 0.9 (normalized + // frequency) + zpk<fbase> filt = iir_bandpass(butterworth<fbase>(4), 0.005, 0.9); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), bqs); } plot_save("butterworth_bandpass4", output, options + ", title='4th-order Butterworth filter, bandpass'"); + // -------------------------------------------------------------------------------------- + // ------------------- 8th-Order Chebyshev Type I Lowpass Filter ------------------------ + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_lowpass(chebyshev1<fbase>(8, 2), 0.09); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), bqs); + // Create an 8th-order Chebyshev Type I lowpass filter with a cutoff frequency of 0.09 (normalized + // frequency) and 2 dB ripple in the passband + zpk<fbase> filt = iir_lowpass(chebyshev1<fbase>(8, 2), 0.09); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), bqs); } plot_save("chebyshev1_lowpass8", output, - options + ", title='8th-order Chebyshev type I filter, lowpass'"); + options + ", title='8th-order Chebyshev Type I filter, lowpass'"); + // -------------------------------------------------------------------------------------- + // ------------------- 8th-Order Chebyshev Type II Lowpass Filter ----------------------- + // -------------------------------------------------------------------------------------- { - zpk<fbase> filt = iir_lowpass(chebyshev2<fbase>(8, 80), 0.09); - iir_params<fbase> bqs = to_sos(filt); // to_sos is expensive, keep iir_params if reused - output = iir(unitimpulse(), filt); + // Create an 8th-order Chebyshev Type II lowpass filter with a cutoff frequency of 0.09 (normalized + // frequency) and 80 dB attenuation in the stopband + zpk<fbase> filt = iir_lowpass(chebyshev2<fbase>(8, 80), 0.09); + + // Convert the filter to second-order sections (SOS). + // This is an expensive operation, so keep 'iir_params' if it is reused later + iir_params<fbase> bqs = to_sos(filt); + + // Apply the filter to a unit impulse signal to get the filter's impulse response + output = iir(unitimpulse(), filt); } plot_save("chebyshev2_lowpass8", output, options + ", title='8th-order Chebyshev type II filter, lowpass'"); diff --git a/examples/sample_rate_conversion.cpp b/examples/sample_rate_conversion.cpp @@ -10,65 +10,108 @@ using namespace kfr; -constexpr size_t input_sr = 96000; -constexpr size_t output_sr = 44100; -constexpr size_t len = 96000 * 6; +// Define constants for input and output sample rates and the length of the signal +constexpr size_t input_sr = 96000; // Input sample rate (96 kHz) +constexpr size_t output_sr = 44100; // Output sample rate (44.1 kHz) +constexpr size_t len = 96000 * 6; // Length of the signal (6 seconds at 96 kHz) int main() { + // Print the version of the KFR library being used println(library_version()); + // Generate a swept sine wave signal with a duration of 'len' samples univector<fbase> swept_sine = swept(0.5, len); + // -------------------------------------------------------------------------------------- + // ----------------------------- High Quality Resampling -------------------------------- + // -------------------------------------------------------------------------------------- { + // Create a high-quality resampler from input_sr to output_sr auto r = resampler<fbase>(resample_quality::high, output_sr, input_sr); + + // Create a buffer for the resampled signal, taking the resampler delay into account univector<fbase> resampled(len * output_sr / input_sr + r.get_delay()); + + // Perform the resampling process r.process(resampled, swept_sine); + // Write the resampled signal to a WAV file audio_writer_wav<fbase> writer(open_file_for_writing(KFR_FILEPATH("audio_high_quality.wav")), audio_format{ 1, audio_sample_type::i32, output_sr }); writer.write(resampled.data(), resampled.size()); writer.close(); + // Save a plot of the high-quality resampled audio plot_save("audio_high_quality", "audio_high_quality.wav", ""); } + // -------------------------------------------------------------------------------------- + // ----------------------------- Normal Quality Resampling ------------------------------ + // -------------------------------------------------------------------------------------- { + // Create a normal-quality resampler from input_sr to output_sr auto r = resampler<fbase>(resample_quality::normal, output_sr, input_sr); + + // Create a buffer for the resampled signal, taking the resampler delay into account univector<fbase> resampled(len * output_sr / input_sr + r.get_delay()); + + // Perform the resampling process r.process(resampled, swept_sine); + // Write the resampled signal to a WAV file audio_writer_wav<fbase> writer(open_file_for_writing(KFR_FILEPATH("audio_normal_quality.wav")), audio_format{ 1, audio_sample_type::i32, output_sr }); writer.write(resampled.data(), resampled.size()); writer.close(); + // Save a plot of the normal-quality resampled audio plot_save("audio_normal_quality", "audio_normal_quality.wav", ""); } + // -------------------------------------------------------------------------------------- + // ----------------------------- Low Quality Resampling --------------------------------- + // -------------------------------------------------------------------------------------- { + // Create a low-quality resampler from input_sr to output_sr auto r = resampler<fbase>(resample_quality::low, output_sr, input_sr); + + // Create a buffer for the resampled signal, taking the resampler delay into account univector<fbase> resampled(len * output_sr / input_sr + r.get_delay()); + + // Perform the resampling process r.process(resampled, swept_sine); + // Write the resampled signal to a WAV file audio_writer_wav<fbase> writer(open_file_for_writing(KFR_FILEPATH("audio_low_quality.wav")), audio_format{ 1, audio_sample_type::i32, output_sr }); writer.write(resampled.data(), resampled.size()); writer.close(); + // Save a plot of the low-quality resampled audio plot_save("audio_low_quality", "audio_low_quality.wav", ""); } + // -------------------------------------------------------------------------------------- + // ----------------------------- Draft Quality Resampling ------------------------------- + // -------------------------------------------------------------------------------------- { + // Create a draft-quality resampler from input_sr to output_sr auto r = resampler<fbase>(resample_quality::draft, output_sr, input_sr); + + // Create a buffer for the resampled signal, taking the resampler delay into account univector<fbase> resampled(len * output_sr / input_sr + r.get_delay()); + + // Perform the resampling process r.process(resampled, swept_sine); + // Write the resampled signal to a WAV file audio_writer_wav<fbase> writer(open_file_for_writing(KFR_FILEPATH("audio_draft_quality.wav")), audio_format{ 1, audio_sample_type::i32, output_sr }); writer.write(resampled.data(), resampled.size()); writer.close(); + // Save a plot of the draft-quality resampled audio plot_save("audio_draft_quality", "audio_draft_quality.wav", ""); } diff --git a/examples/window.cpp b/examples/window.cpp @@ -12,54 +12,73 @@ using namespace kfr; int main() { + // Print the version of the KFR library being used println(library_version()); + // Define options for plotting DSP data const std::string options = "freqresp=True, dots=True, padwidth=1024, " "log_freq=False, horizontal=False, normalized_freq=True"; + // Declare a univector of type fbase with a size of 64 to hold the window function output univector<fbase, 64> output; + + // Generate a Hann window function and save the plot output = window_hann(output.size()); plot_save("window_hann", output, options + ", title='Hann window'"); + // Generate a Hamming window function and save the plot output = window_hamming(output.size()); plot_save("window_hamming", output, options + ", title='Hamming window'"); + // Generate a Blackman window function and save the plot output = window_blackman(output.size()); plot_save("window_blackman", output, options + ", title='Blackman window'"); + // Generate a Blackman-Harris window function and save the plot output = window_blackman_harris(output.size()); plot_save("window_blackman_harris", output, options + ", title='Blackman-Harris window'"); + // Generate a Gaussian window function and save the plot output = window_gaussian(output.size()); plot_save("window_gaussian", output, options + ", title='Gaussian window'"); + // Generate a Triangular window function and save the plot output = window_triangular(output.size()); plot_save("window_triangular", output, options + ", title='Triangular window'"); + // Generate a Bartlett window function and save the plot output = window_bartlett(output.size()); plot_save("window_bartlett", output, options + ", title='Bartlett window'"); + // Generate a Cosine window function and save the plot output = window_cosine(output.size()); plot_save("window_cosine", output, options + ", title='Cosine window'"); + // Generate a Cosine window function compatible with numpy and save the plot output = window_cosine_np(output.size()); plot_save("window_cosine_np", output, options + ", title='Cosine window (numpy compatible)'"); + // Generate a Bartlett-Hann window function and save the plot output = window_bartlett_hann(output.size()); plot_save("window_bartlett_hann", output, options + ", title='Bartlett-Hann window'"); + // Generate a Bohman window function and save the plot output = window_bohman(output.size()); plot_save("window_bohman", output, options + ", title='Bohman window'"); + // Generate a Lanczos window function and save the plot output = window_lanczos(output.size()); plot_save("window_lanczos", output, options + ", title='Lanczos window'"); + // Generate a Flat top window function and save the plot output = window_flattop(output.size()); plot_save("window_flattop", output, options + ", title='Flat top window'"); + // Generate a Kaiser window function with a beta value of 2.5 and save the plot output = window_kaiser(output.size(), 2.5); plot_save("window_kaiser", output, options + ", title='Kaiser window'"); + // Generate a Planck-taper window function with an epsilon value of 0.1 and save the plot output = window_planck_taper(output.size(), 0.1); plot_save("window_planck_taper", output, options + ", title='Planck-taper window'");