36#ifndef FILTER_RTPROCESSING_H
37#define FILTER_RTPROCESSING_H
53#include <QSharedPointer>
54#include <QtConcurrent/QtConcurrent>
61#include <unsupported/Eigen/FFT>
106 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
114 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
115 bool bUseThreads =
true);
131 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
133 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
134 bool bUseThreads =
false);
163 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
164 bool bUseThreads =
true,
165 bool bKeepOverhead =
false);
182 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
183 bool bUseThreads =
true,
184 bool bKeepOverhead =
false);
199 const Eigen::RowVectorXi& vecPicks,
201 bool bUseThreads =
true);
221 typedef QSharedPointer<FilterOverlapAdd>
SPtr;
222 typedef QSharedPointer<const FilterOverlapAdd>
ConstSPtr;
243 Eigen::MatrixXd
calculate(
const Eigen::MatrixXd& matData,
251 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
252 bool bFilterEnd =
true,
253 bool bUseThreads =
true,
254 bool bKeepOverhead =
false);
269 Eigen::MatrixXd
calculate(
const Eigen::MatrixXd& mataData,
271 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
272 bool bFilterEnd =
true,
273 bool bUseThreads =
true,
274 bool bKeepOverhead =
false);
283 Eigen::MatrixXd m_matOverlapBack;
284 Eigen::MatrixXd m_matOverlapFront;
realtime library export/import macros.
#define RTPROCESINGSHARED_EXPORT
The FilterKernel class represents a filter object that generates the FIR filter coefficients using Pa...
FiffInfo class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Real-time signal processing (filtering, averaging, HPI fitting, noise reduction).
Eigen::MatrixXd filterData(const Eigen::MatrixXd &matData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=1024, int designMethod=FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bUseThreads=true, bool bKeepOverhead=false)
void filterChannel(FilterObject &channelDataTime)
Eigen::MatrixXd filterDataBlock(const Eigen::MatrixXd &mataData, const Eigen::RowVectorXi &vecPicks, const RTPROCESSINGLIB::FilterKernel &filterKernel, bool bUseThreads=true)
bool filterFile(QIODevice &pIODevice, QSharedPointer< FIFFLIB::FiffRawData > pFiffRawData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=4096, int designMethod=FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bUseThreads=true)
FIFF raw measurement data.
Lightweight filter configuration holding kernel coefficients and overlap-add state for one channel.
FilterKernel filterKernel
Eigen::RowVectorXd vecData
Applies FIR filtering via FFT-based overlap-add convolution for continuous data streams.
Eigen::MatrixXd calculate(const Eigen::MatrixXd &mataData, const RTPROCESSINGLIB::FilterKernel &filterKernel, const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bFilterEnd=true, bool bUseThreads=true, bool bKeepOverhead=false)
QSharedPointer< const FilterOverlapAdd > ConstSPtr
QSharedPointer< FilterOverlapAdd > SPtr
Eigen::MatrixXd calculate(const Eigen::MatrixXd &matData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=1024, int designMethod=FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bFilterEnd=true, bool bUseThreads=true, bool bKeepOverhead=false)
Named filter-design parameter descriptor holding a human-readable name and description (e....
The FilterKernel class provides methods to create/design a FIR filter kernel.
static QVector< FilterParameter > m_designMethods