LASP 1.0
Library for Acoustic Signal Processing
Loading...
Searching...
No Matches
lasp_slm.cpp
Go to the documentation of this file.
1/* #define DEBUGTRACE_ENABLED */
2#include "debugtrace.hpp"
3#include "lasp_slm.h"
4#include "lasp_thread.h"
5#include <algorithm>
6#include <cmath>
7#include <future>
8#include <memory>
9
10using std::cerr;
11using std::endl;
12using rte = std::runtime_error;
13using std::unique_ptr;
14
15SLM::SLM(const d fs, const d Lref, const us downsampling_fac, const d tau,
16 std::unique_ptr<Filter> pre_filter,
17 std::vector<std::unique_ptr<Filter>> bandpass)
18 : _pre_filter(std::move(pre_filter)), _bandpass(std::move(bandpass)),
19 _alpha(exp(-1 / (fs * tau))),
20 _sp_storage(_bandpass.size(), arma::fill::zeros), // Storage for
21 // components of
22 // single pole low pass
23 // filter
24 Lrefsq(Lref * Lref), // Reference level
25 downsampling_fac(downsampling_fac),
26
27 // Initalize mean square
28 Pm(_bandpass.size(), arma::fill::zeros),
29
30 // Initalize max
31 Pmax(_bandpass.size(), arma::fill::zeros),
32
33 // Initalize peak
34 Ppeak(_bandpass.size(), arma::fill::zeros)
35
36{
37 DEBUGTRACE_ENTER;
38 DEBUGTRACE_PRINT(_alpha);
39
40 if (Lref <= 0) {
41 throw rte("Invalid reference level");
42 }
43 if (tau < 0) {
44 throw rte("Invalid time constant for Single pole lowpass filter");
45 }
46 if (fs <= 0) {
47 throw rte("Invalid sampling frequency");
48 }
49}
51
59std::vector<unique_ptr<Filter>> createBandPass(const dmat &coefs) {
60 DEBUGTRACE_ENTER;
61 std::vector<unique_ptr<Filter>> bf;
62 for (us colno = 0; colno < coefs.n_cols; colno++) {
63 bf.emplace_back(std::make_unique<SeriesBiquad>(coefs.col(colno)));
64 }
65 return bf;
66}
67us SLM::suggestedDownSamplingFac(const d fs, const d tau) {
68 if (fs <= 0)
69 throw rte("Invalid sampling frequency");
70 // A reasonable 'framerate' for the sound level meter, based on the
71 // filtering time constant.
72 if (tau > 0) {
73 d fs_slm = 10 / tau;
74 if (fs_slm < 30) {
75 fs_slm = 30;
76 }
77 return std::max((us)1, static_cast<us>(fs / fs_slm));
78 } else {
79 return 1;
80 }
81}
82
83SLM SLM::fromBiquads(const d fs, const d Lref, const us downsampling_fac,
84 const d tau, const vd &pre_filter_coefs,
85 const dmat &bandpass_coefs) {
86 DEBUGTRACE_ENTER;
87
88 return SLM(fs, Lref, downsampling_fac, tau,
89 std::make_unique<SeriesBiquad>(pre_filter_coefs),
90 createBandPass(bandpass_coefs));
91}
92SLM SLM::fromBiquads(const d fs, const d Lref, const us downsampling_fac,
93 const d tau, const dmat &bandpass_coefs) {
94
95 DEBUGTRACE_ENTER;
96
97 return SLM(fs, Lref, downsampling_fac, tau,
98 nullptr, // Pre-filter
99 createBandPass(bandpass_coefs) // Bandpass coefficients
100 );
101}
102
103vd SLM::run_single(vd work, const us i) {
104
105 // Filter input in-place
106 _bandpass[i]->filter(work);
107
108 /* cerr << "Filter done" << endl; */
109
110 // Square input --> Signal powers
111 /* work.transform([](d j) { return j * j; }); */
112 work %= work;
113
114 // Compute peak level, that is before single-pole low pass filter
115 Ppeak(i) = std::max(Ppeak(i), arma::max(work));
116
117 // Create copy of N, as we run this in multiple threads.
118 us N_local = N;
119 DEBUGTRACE_PRINT(N);
120
121 // Obtain storage of single_pole low pass filter
122 d cur_storage = _sp_storage(i);
123
124 for (us j = 0; j < work.n_rows; j++) {
125 // Update mean square of signal, work is here still signal power
126 Pm(i) = (Pm(i) * static_cast<d>(N_local) + work(j)) /
127 (static_cast<d>(N_local) + 1);
128
129 N_local++;
130
131 cur_storage = _alpha * cur_storage + (1 - _alpha) * work(j);
132
133 // Now work is single-pole lowpassed signal power
134 work(j) = cur_storage;
135 }
136
137 // And update storage of low-pass filter
138 _sp_storage(i) = cur_storage;
139
140 Pmax(i) = std::max(Pmax(i), arma::max(work));
141
142 // Convert to levels in dB
143 work = 10 * arma::log10((work + arma::datum::eps) / Lrefsq);
144
145 return work;
146}
147
148dmat SLM::run(const vd &input_orig) {
149
150 DEBUGTRACE_ENTER;
151 vd input = input_orig;
152
153 // _pre_filter filters in-place
154 if (_pre_filter) {
155 _pre_filter->filter(input);
156 }
157
158 // Fan out over multiple threads, as it is typically a heavy load
159 dmat res(input.n_rows, _bandpass.size());
160
161 // Perform operations in-place.
162
163#pragma omp parallel for
164 for (us i = 0; i < _bandpass.size(); i++) {
165 res.col(i) = run_single(input, i);
166 /* DEBUGTRACE_PRINT(_bandpass.size()); */
167 /* DEBUGTRACE_PRINT(res.n_cols); */
168 /* DEBUGTRACE_PRINT(res.n_rows); */
169 /* DEBUGTRACE_PRINT(futs.size()); */
170
171 // Update the total number of samples harvested so far. NOTE: *This should
172 // be done AFTER the threads are done!!!*
173 }
174 N += input.n_rows;
175
176 // Downsample, if applicable
177 if (downsampling_fac > 1) {
178 dmat res_ds;
179 us rowno = 0;
180 while (cur_offset < res.n_rows) {
181 res_ds.insert_rows(rowno, res.row(cur_offset));
182 rowno++;
183
184 cur_offset += downsampling_fac;
185 }
186 cur_offset -= res.n_rows;
187 // Instead, return a downsampled version
188 return res_ds;
189 }
190
191 return res;
192}
194 Pm.zeros();
195 Pmax.zeros();
196 Ppeak.zeros();
197 for (auto &f : _bandpass) {
198 f.reset();
199 }
200 if (_pre_filter) {
201 _pre_filter->reset();
202 }
203 _sp_storage.zeros();
204 cur_offset = 0;
205
206 N = 0;
207}
Sound Level Meter implementation that gives a result for each channel. A channel is the result of a f...
Definition lasp_slm.h:17
vd Pmax
Public storage for the maximum signal power, after single pole low-pass filter.
Definition lasp_slm.h:48
us N
Definition lasp_slm.h:55
vd Pm
Public storage for the mean of the square of the signal.
Definition lasp_slm.h:43
static us suggestedDownSamplingFac(const d fs, const d tw)
Comput a 'suggested' downsampling factor, i.e. a lower frame rate at which sound level meter values a...
Definition lasp_slm.cpp:67
dmat run(const vd &input)
Run the sound level meter on given input data. Return downsampled level data for each filterbank chan...
Definition lasp_slm.cpp:148
void reset()
Reset state related to samples acquired. All filters reset to zero. Start again from no history.
Definition lasp_slm.cpp:193
vd Ppeak
Storage for maximum computed signal power so far.
Definition lasp_slm.h:53
static SLM fromBiquads(const d fs, const d Lref, const us downsampling_fac, const d tau, const vd &pre_filter_coefs, const dmat &bandpass_coefs)
Convenience function to create a Sound Level meter from Biquad filters only.
Definition lasp_slm.cpp:83
SLM(const d fs, const d Lref, const us downsampling_fac, const d tau, std::unique_ptr< Filter > pre_filter, std::vector< std::unique_ptr< Filter > > bandpass)
Initialize a Sound Level Meter.
Definition lasp_slm.cpp:15
~SLM()
Definition lasp_slm.cpp:50
std::runtime_error rte
Definition lasp_daq.cpp:16
arma::Col< d > vd
arma::Mat< d > dmat
std::vector< unique_ptr< Filter > > createBandPass(const dmat &coefs)
Create set bandpass filters from filter coefficients.
Definition lasp_slm.cpp:59
size_t us
We often use boolean values.
Definition lasp_types.h:29