Removed carma usage. This bugs as it changes the default allocator for armadillo, which does not work properly when native C++ threads do not register themselves as a thread for Python. We implemented the required convertors between Python and Armadillo ourselves. This also fixes the problem of Carma not being officially Python3.10 proof.

This commit is contained in:
Anne de Jong 2022-10-11 09:43:36 +02:00
parent ce402acd15
commit e4ab05d36a
6 changed files with 116 additions and 98 deletions

View File

@ -11,6 +11,7 @@ using std::cerr;
using std::endl;
using rte = std::runtime_error;
InDataHandler::InDataHandler(StreamMgr &mgr) : _mgr(mgr) { DEBUGTRACE_ENTER; }
void InDataHandler::start() {
DEBUGTRACE_ENTER;
@ -59,6 +60,7 @@ void StreamMgr::checkRightThread() const {
void StreamMgr::rescanDAQDevices(bool background,
std::function<void()> callback) {
DEBUGTRACE_ENTER;
auto& pool = getPool();
checkRightThread();
if (!_devices_mtx.try_lock()) {
@ -74,7 +76,7 @@ void StreamMgr::rescanDAQDevices(bool background,
if (!background) {
rescanDAQDevices_impl(callback);
} else {
/* pool.push_task(&StreamMgr::rescanDAQDevices_impl, this, callback); */
pool.push_task(&StreamMgr::rescanDAQDevices_impl, this, callback);
}
}
void StreamMgr::rescanDAQDevices_impl(std::function<void()> callback) {

View File

@ -13,56 +13,56 @@ class StreamMgr;
class InDataHandler {
protected:
StreamMgr &_mgr;
protected:
StreamMgr &_mgr;
#if LASP_DEBUG == 1
std::atomic<bool> stopCalled{false};
std::atomic<bool> stopCalled{false};
#endif
public:
virtual ~InDataHandler();
public:
virtual ~InDataHandler();
/**
* @brief When constructed, the handler is added to the stream manager, which
* will call the handlers's inCallback() until stop() is called.
*
* @param mgr Stream manager.
*/
InDataHandler(StreamMgr &mgr);
/**
* @brief When constructed, the handler is added to the stream manager, which
* will call the handlers's inCallback() until stop() is called.
*
* @param mgr Stream manager.
*/
InDataHandler(StreamMgr &mgr);
/**
* @brief This function is called when input data from a DAQ is available.
*
* @param daqdata Input data from DAQ
*
* @return true if no error. False to stop the stream from running.
*/
virtual bool inCallback(const DaqData &daqdata) = 0;
/**
* @brief This function is called when input data from a DAQ is available.
*
* @param daqdata Input data from DAQ
*
* @return true if no error. False to stop the stream from running.
*/
virtual bool inCallback(const DaqData &daqdata) = 0;
/**
* @brief Reset in-data handler.
*
* @param daq New DAQ configuration of inCallback(). If nullptr is given,
* it means that the stream is stopped.
*/
virtual void reset(const Daq* daq = nullptr) = 0;
/**
* @brief Reset in-data handler.
*
* @param daq New DAQ configuration of inCallback(). If nullptr is given,
* it means that the stream is stopped.
*/
virtual void reset(const Daq *daq = nullptr) = 0;
/**
* @brief This function should be called from the constructor of the
* implementation of InDataHandler. It will start the stream's calling of
* inCallback().
*/
void start();
/**
* @brief This function should be called from the constructor of the
* implementation of InDataHandler. It will start the stream's calling of
* inCallback().
*/
void start();
/**
* @brief This function should be called from the destructor of derived
* classes, to disable the calls to inCallback(), such that proper
* destruction of the object is allowed and no other threads call methods
* from the object. It removes the inCallback() from the callback list of the
* StreamMgr(). **Failing to call this function results in deadlocks, errors
* like "pure virtual function called", or other**.
*/
void stop();
/**
* @brief This function should be called from the destructor of derived
* classes, to disable the calls to inCallback(), such that proper
* destruction of the object is allowed and no other threads call methods
* from the object. It removes the inCallback() from the callback list of the
* StreamMgr(). **Failing to call this function results in deadlocks, errors
* like "pure virtual function called", or other**.
*/
void stop();
};
/**
@ -152,8 +152,8 @@ class StreamMgr {
* @param callback Function to call when complete.
*/
void
rescanDAQDevices(bool background = false,
std::function<void()> callback = std::function<void()>());
rescanDAQDevices(bool background = false,
std::function<void()> callback = std::function<void()>());
/**
* @brief Start a stream based on given configuration.
@ -173,13 +173,13 @@ class StreamMgr {
return getStreamStatus(type).runningOK();
}
bool isStreamRunning(const StreamType type) const {
switch(type) {
case(StreamType::input):
return bool(_inputStream);
break;
case(StreamType::output):
return bool(_outputStream);
break;
switch (type) {
case (StreamType::input):
return bool(_inputStream);
break;
case (StreamType::output):
return bool(_outputStream);
break;
}
return false;
}
@ -225,7 +225,7 @@ class StreamMgr {
*/
void setSiggen(std::shared_ptr<Siggen> s);
private:
private:
bool inCallback(const DaqData &data);
bool outCallback(DaqData &data);

View File

@ -125,7 +125,7 @@ class SLM:
if prefilter is not None:
self.slm = cppSLM.fromBiquads(fs, level_ref_value, dsfac,
tw[0],
prefilter, sos)
prefilter.flatten(), sos)
else:
self.slm = cppSLM.fromBiquads(fs, level_ref_value, dsfac,
tw[0],
@ -145,11 +145,7 @@ class SLM:
data: one-dimensional input data
"""
assert data.ndim == 2
assert data.shape[1] == 1, "invalid number of channels, should be 1"
if data.shape[0] == 0:
return {}
assert data.ndim == 1
levels = self.slm.run(data)
@ -197,10 +193,10 @@ class SLM:
logging.debug(output['mid'])
if self.include_overall and self.fbdesigner is not None:
output['overall'] = dat[-1,0]
output['y'] = np.asarray(dat[:-1,0])
output['overall'] = dat[-1]
output['y'] = np.asarray(dat[:-1])
else:
output['y'] = np.asarray(dat[:,0])
output['y'] = np.asarray(dat[:])
return output
def Leq(self):

View File

@ -1,15 +1,16 @@
#include <carma>
#include "arma_npy.h"
#include "lasp_avpowerspectra.h"
#include "lasp_biquadbank.h"
#include "lasp_fft.h"
#include "lasp_streammgr.h"
#include "lasp_filter.h"
#include "lasp_slm.h"
#include "lasp_streammgr.h"
#include "lasp_window.h"
#include <iostream>
#include <pybind11/pybind11.h>
using std::cerr;
using std::endl;
namespace py = pybind11;
/**
@ -23,6 +24,7 @@ namespace py = pybind11;
*
* @param m The Python module to add classes and methods to
*/
void init_dsp(py::module &m) {
py::class_<Fft> fft(m, "Fft");
@ -53,22 +55,30 @@ void init_dsp(py::module &m) {
py::class_<SeriesBiquad, std::shared_ptr<SeriesBiquad>> sbq(m, "SeriesBiquad",
filter);
sbq.def(py::init<const vd &>());
sbq.def("filter", [](SeriesBiquad &s, const vd &input) {
vd res = input;
sbq.def("filter", [](SeriesBiquad &s, dpyarray input) {
vd res = NpyToCol<d, true>(input);
s.filter(res);
return res;
return ColToNpy<d>(res);
});
/// BiquadBank
py::class_<BiquadBank, std::shared_ptr<BiquadBank>> bqb(m, "BiquadBank");
bqb.def(py::init<const dmat &, const vd *>());
bqb.def("setGains", &BiquadBank::setGains);
bqb.def("filter", &BiquadBank::filter);
bqb.def("setGains", [](BiquadBank &b, dpyarray gains) {
b.setGains(NpyToCol(gains));
});
bqb.def("filter", [](BiquadBank &b, dpyarray input) {
vd inout = NpyToCol<d, true>(input);
b.filter(inout);
return ColToNpy(inout);
});
/// PowerSpectra
py::class_<PowerSpectra> ps(m, "PowerSpectra");
ps.def(py::init<const us, const Window::WindowType>());
ps.def("compute", &PowerSpectra::compute);
ps.def("compute", [](PowerSpectra &ps, dpyarray input) {
return CubeToNpy<c>(ps.compute(NpyToMat<d, false>(input)));
});
/// AvPowerSpectra
py::class_<AvPowerSpectra> aps(m, "AvPowerSpectra");
@ -77,29 +87,37 @@ void init_dsp(py::module &m) {
py::arg("windowType") = Window::WindowType::Hann,
py::arg("overlap_percentage") = 50.0, py::arg("time_constant") = -1);
aps.def("compute", [](AvPowerSpectra &aps, const dmat &timedata) {
std::optional<arma::cx_cube> res = aps.compute(timedata);
return res.value_or(arma::cx_cube(0, 0, 0));
aps.def("compute", [](AvPowerSpectra &aps, dpyarray timedata) {
std::optional<arma::cx_cube> res =
aps.compute(NpyToMat<d, false>(timedata));
return CubeToNpy<c>(res.value_or(cube(0, 0, 0)));
});
py::class_<SLM> slm(m, "cppSLM");
slm.def_static(
"fromBiquads",
py::overload_cast<const d, const d, const us, const d, const dmat &>(
&SLM::fromBiquads));
slm.def_static(
"fromBiquads",
py::overload_cast<const d, const d, const us, const d, const vd &,
const dmat &>(&SLM::fromBiquads));
slm.def("run", &SLM::run);
slm.def_readonly("Pm", &SLM::Pm);
slm.def_readonly("Pmax", &SLM::Pmax);
slm.def_readonly("Ppeak", &SLM::Ppeak);
slm.def("Lpeak", &SLM::Lpeak);
slm.def("Leq", &SLM::Leq);
slm.def("Lmax", &SLM::Lmax);
slm.def_static("suggestedDownSamplingFac", &SLM::suggestedDownSamplingFac);
slm.def_static("fromBiquads", [](const d fs, const d Lref, const us ds,
const d tau, dpyarray bandpass) {
return SLM::fromBiquads(fs, Lref, ds, tau, NpyToMat<d, false>(bandpass));
});
slm.def_static("fromBiquads", [](const d fs, const d Lref, const us ds,
const d tau, dpyarray prefilter,
py::array_t<d> bandpass) {
return SLM::fromBiquads(fs, Lref, ds, tau, NpyToCol<d, false>(prefilter),
NpyToMat<d, false>(bandpass));
});
slm.def("run", [](SLM &slm, dpyarray in) {
return MatToNpy<d>(slm.run(NpyToCol<d, false>(in)));
});
slm.def("Pm", [](const SLM &slm) { return ColToNpy<d>(slm.Pm); });
slm.def("Pmax", [](const SLM &slm) { return ColToNpy<d>(slm.Pmax); });
slm.def("Ppeak", [](const SLM &slm) { return ColToNpy<d>(slm.Ppeak); });
slm.def("Leq", [](const SLM &slm) { return ColToNpy<d>(slm.Leq()); });
slm.def("Lmax", [](const SLM &slm) { return ColToNpy<d>(slm.Lmax()); });
slm.def("Lpeak", [](const SLM &slm) { return ColToNpy<d>(slm.Lpeak()); });
slm.def_static("suggestedDownSamplingFac", &SLM::suggestedDownSamplingFac);
}
/** @} */

View File

@ -1,6 +1,6 @@
#include <carma>
#include <armadillo>
#define DEBUGTRACE_ENABLED
#include "arma_npy.h"
#include "debugtrace.hpp"
#include "lasp_ppm.h"
#include "lasp_rtaps.h"
@ -25,7 +25,8 @@ namespace py = pybind11;
*
* @return Numpy array
*/
template <typename T,bool copy=false> py::array_t<T> getPyArrayNoCpy(const DaqData &d) {
template <typename T, bool copy = false>
py::array_t<T> getPyArrayNoCpy(const DaqData &d) {
// https://github.com/pybind/pybind11/issues/323
//
// When a valid object is passed as 'base', it tells pybind not to take
@ -57,7 +58,8 @@ template <typename T,bool copy=false> py::array_t<T> getPyArrayNoCpy(const DaqDa
);
}
template <typename T,bool copy=false> py::array_t<d> dmat_to_ndarray(const DaqData &d) {
template <typename T, bool copy = false>
py::array_t<d> dmat_to_ndarray(const DaqData &d) {
// https://github.com/pybind/pybind11/issues/323
//
// When a valid object is passed as 'base', it tells pybind not to take
@ -183,10 +185,11 @@ void init_datahandler(py::module &m) {
ppm.def(py::init<StreamMgr &>());
ppm.def("getCurrentValue", [](const PPMHandler &ppm) {
std::tuple<vd, arma::uvec> tp = ppm.getCurrentValue();
return py::make_tuple(carma::col_to_arr(std::get<0>(tp)),
carma::col_to_arr(std::get<1>(tp)));
return py::make_tuple(ColToNpy<d>(std::get<0>(tp)),
ColToNpy<arma::uword>(std::get<1>(tp)));
});
/// Real time Aps
@ -214,8 +217,8 @@ void init_datahandler(py::module &m) {
rtaps.def("getCurrentValue", [](RtAps &rt) {
std::unique_ptr<cube> val = rt.getCurrentValue();
if (val) {
return carma::cube_to_arr<c>(std::move(*val));
return CubeToNpy<c> (*val);
}
return carma::cube_to_arr<c>(cube(1, 0, 0));
return CubeToNpy<c>(cube(1,0,0));
});
}

View File

@ -1,4 +1,3 @@
#include <carma>
#include "lasp_avpowerspectra.h"
#include "lasp_biquadbank.h"
#include "lasp_fft.h"
@ -33,7 +32,7 @@ void init_siggen(py::module &m) {
siggen.def("setMute", &Siggen::setMute);
siggen.def("setLevel", &Siggen::setLevel);
siggen.def("setLevel", &Siggen::setLevel, py::arg("newLevel"), py::arg("dB") = true);
siggen.def("genSignal", &Siggen::genSignal);
/* siggen.def("genSignal", &Siggen::genSignal); */
siggen.def("setFilter", &Siggen::setFilter);
py::class_<Noise,std::shared_ptr<Noise>> noise(m, "Noise", siggen);