v2.0.0
Loading...
Searching...
No Matches
filter.h
Go to the documentation of this file.
1//=============================================================================================================
35
36#ifndef FILTER_RTPROCESSING_H
37#define FILTER_RTPROCESSING_H
38
39//=============================================================================================================
40// INCLUDES
41//=============================================================================================================
42
43#include "rtprocessing_global.h"
44
46
47#include <fiff/fiff_info.h>
48
49//=============================================================================================================
50// QT INCLUDES
51//=============================================================================================================
52
53#include <QSharedPointer>
54#include <QtConcurrent/QtConcurrent>
55
56//=============================================================================================================
57// EIGEN INCLUDES
58//=============================================================================================================
59
60#include <Eigen/Core>
61#include <unsupported/Eigen/FFT>
62
63//=============================================================================================================
64// FORWARD DECLARATIONS
65//=============================================================================================================
66
67namespace FIFFLIB {
68 class FiffRawData;
69}
70
71//=============================================================================================================
72// DEFINE NAMESPACE RTPROCESSINGLIB
73//=============================================================================================================
74
75namespace RTPROCESSINGLIB
76{
77
81typedef struct {
83 int iRow;
84 Eigen::RowVectorXd vecData;
86
87//=========================================================================================================
105RTPROCESINGSHARED_EXPORT bool filterFile(QIODevice& pIODevice,
106 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
107 int type,
108 double dCenterfreq,
109 double dBandwidth,
110 double dTransition,
111 double dSFreq,
112 int iOrder = 4096,
113 int designMethod = FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")),
114 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
115 bool bUseThreads = true);
116
117//=========================================================================================================
130RTPROCESINGSHARED_EXPORT bool filterFile(QIODevice& pIODevice,
131 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
132 const RTPROCESSINGLIB::FilterKernel& filterKernel,
133 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
134 bool bUseThreads = false);
135
136//=========================================================================================================
155RTPROCESINGSHARED_EXPORT Eigen::MatrixXd filterData(const Eigen::MatrixXd& matData,
156 int type,
157 double dCenterfreq,
158 double dBandwidth,
159 double dTransition,
160 double dSFreq,
161 int iOrder = 1024,
162 int designMethod = FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")),
163 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
164 bool bUseThreads = true,
165 bool bKeepOverhead = false);
166
167//=========================================================================================================
180RTPROCESINGSHARED_EXPORT Eigen::MatrixXd filterData(const Eigen::MatrixXd& mataData,
181 const RTPROCESSINGLIB::FilterKernel& filterKernel,
182 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
183 bool bUseThreads = true,
184 bool bKeepOverhead = false);
185
186//=========================================================================================================
198RTPROCESINGSHARED_EXPORT Eigen::MatrixXd filterDataBlock(const Eigen::MatrixXd& mataData,
199 const Eigen::RowVectorXi& vecPicks,
200 const RTPROCESSINGLIB::FilterKernel& filterKernel,
201 bool bUseThreads = true);
202
203//=========================================================================================================
210
211//=============================================================================================================
219{
220public:
221 typedef QSharedPointer<FilterOverlapAdd> SPtr;
222 typedef QSharedPointer<const FilterOverlapAdd> ConstSPtr;
223
224 //=========================================================================================================
243 Eigen::MatrixXd calculate(const Eigen::MatrixXd& matData,
244 int type,
245 double dCenterfreq,
246 double dBandwidth,
247 double dTransition,
248 double dSFreq,
249 int iOrder = 1024,
250 int designMethod = FilterKernel::m_designMethods.indexOf(FilterParameter("Cosine")),
251 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
252 bool bFilterEnd = true,
253 bool bUseThreads = true,
254 bool bKeepOverhead = false);
255
256 //=========================================================================================================
269 Eigen::MatrixXd calculate(const Eigen::MatrixXd& mataData,
270 const RTPROCESSINGLIB::FilterKernel& filterKernel,
271 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
272 bool bFilterEnd = true,
273 bool bUseThreads = true,
274 bool bKeepOverhead = false);
275
276 //=========================================================================================================
280 void reset();
281
282private:
283 Eigen::MatrixXd m_matOverlapBack;
284 Eigen::MatrixXd m_matOverlapFront;
285};
286
287//=============================================================================================================
288// INLINE DEFINITIONS
289//=============================================================================================================
290
291} // NAMESPACE
292
293#endif // FILTER_RTPROCESSING_H
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)
Definition filter.cpp:376
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.
Definition filter.h:81
FilterKernel filterKernel
Definition filter.h:82
Eigen::RowVectorXd vecData
Definition filter.h:84
Applies FIR filtering via FFT-based overlap-add convolution for continuous data streams.
Definition filter.h:219
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
Definition filter.h:222
QSharedPointer< FilterOverlapAdd > SPtr
Definition filter.h:221
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