v2.0.0
Loading...
Searching...
No Matches
rt_filter.h
Go to the documentation of this file.
1//=============================================================================================================
35
36#ifndef RT_FILTER_RT_H
37#define RT_FILTER_RT_H
38
39//=============================================================================================================
40// INCLUDES
41//=============================================================================================================
42
43#include "../dsp_global.h"
44
45#include "../filterkernel.h"
46
47#include <fiff/fiff_info.h>
48#include <fiff/fiff_evoked.h>
49
50//=============================================================================================================
51// QT INCLUDES
52//=============================================================================================================
53
54#include <QSharedPointer>
55#include <QtConcurrent/QtConcurrent>
56
57//=============================================================================================================
58// EIGEN INCLUDES
59//=============================================================================================================
60
61#include <Eigen/Core>
62#include <unsupported/Eigen/FFT>
63
64//=============================================================================================================
65// FORWARD DECLARATIONS
66//=============================================================================================================
67
68namespace FIFFLIB {
69 class FiffRawData;
70}
71
72//=============================================================================================================
73// DEFINE NAMESPACE RTPROCESSINGLIB
74//=============================================================================================================
75
76namespace RTPROCESSINGLIB
77{
78
87
88//=========================================================================================================
106DSPSHARED_EXPORT bool filterFile(QIODevice& pIODevice,
107 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
108 int type,
109 double dCenterfreq,
110 double dBandwidth,
111 double dTransition,
112 double dSFreq,
113 int iOrder = 4096,
114 int designMethod = UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")),
115 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
116 bool bUseThreads = true);
117
118//=========================================================================================================
131DSPSHARED_EXPORT bool filterFile(QIODevice& pIODevice,
132 QSharedPointer<FIFFLIB::FiffRawData> pFiffRawData,
133 const UTILSLIB::FilterKernel& filterKernel,
134 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
135 bool bUseThreads = false);
136
137//=========================================================================================================
156DSPSHARED_EXPORT Eigen::MatrixXd filterData(const Eigen::MatrixXd& matData,
157 int type,
158 double dCenterfreq,
159 double dBandwidth,
160 double dTransition,
161 double dSFreq,
162 int iOrder = 1024,
163 int designMethod = UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")),
164 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
165 bool bUseThreads = true,
166 bool bKeepOverhead = false);
167
168//=========================================================================================================
181DSPSHARED_EXPORT Eigen::MatrixXd filterData(const Eigen::MatrixXd& matData,
182 const UTILSLIB::FilterKernel& filterKernel,
183 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
184 bool bUseThreads = true,
185 bool bKeepOverhead = false);
186
187//=========================================================================================================
199DSPSHARED_EXPORT Eigen::MatrixXd filterDataBlock(const Eigen::MatrixXd& matData,
200 const Eigen::RowVectorXi& vecPicks,
201 const UTILSLIB::FilterKernel& filterKernel,
202 bool bUseThreads = true);
203
204//=========================================================================================================
210DSPSHARED_EXPORT void filterChannel(FilterObject &channelDataTime);
211
212//=========================================================================================================
234 const Eigen::MatrixXi& matEvents,
235 float fTMinS,
236 float fTMaxS,
237 qint32 eventType,
238 bool bApplyBaseline,
239 float fTBaselineFromS,
240 float fTBaselineToS,
241 const QMap<QString,double>& mapReject,
242 const UTILSLIB::FilterKernel& filterKernel,
243 const QStringList &lExcludeChs = QStringList(),
244 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi());
245
246//=============================================================================================================
254{
255public:
256 typedef QSharedPointer<FilterOverlapAdd> SPtr;
257 typedef QSharedPointer<const FilterOverlapAdd> ConstSPtr;
258
259 //=========================================================================================================
278 Eigen::MatrixXd calculate(const Eigen::MatrixXd& matData,
279 int type,
280 double dCenterfreq,
281 double dBandwidth,
282 double dTransition,
283 double dSFreq,
284 int iOrder = 1024,
285 int designMethod = UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")),
286 const Eigen::RowVectorXi &vecPicks = Eigen::RowVectorXi(),
287 bool bFilterEnd = true,
288 bool bUseThreads = true,
289 bool bKeepOverhead = false);
290
291 //=========================================================================================================
304 Eigen::MatrixXd calculate(const Eigen::MatrixXd& matData,
305 const UTILSLIB::FilterKernel& filterKernel,
306 const Eigen::RowVectorXi& vecPicks = Eigen::RowVectorXi(),
307 bool bFilterEnd = true,
308 bool bUseThreads = true,
309 bool bKeepOverhead = false);
310
311 //=========================================================================================================
315 void reset();
316
317private:
318 Eigen::MatrixXd m_matOverlapBack;
319 Eigen::MatrixXd m_matOverlapFront;
320};
321
322//=============================================================================================================
323// INLINE DEFINITIONS
324//=============================================================================================================
325
326} // NAMESPACE
327
328#endif // RT_FILTER_RT_H
FiffInfo class declaration.
FiffEvoked class declaration.
The FilterKernel class represents a filter object that generates the FIR filter coefficients using Pa...
dsp library export/import macros.
#define DSPSHARED_EXPORT
Definition dsp_global.h:56
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
DSPSHARED_EXPORT FIFFLIB::FiffEvoked computeFilteredAverage(const FIFFLIB::FiffRawData &raw, const Eigen::MatrixXi &matEvents, float fTMinS, float fTMaxS, qint32 eventType, bool bApplyBaseline, float fTBaselineFromS, float fTBaselineToS, const QMap< QString, double > &mapReject, const UTILSLIB::FilterKernel &filterKernel, const QStringList &lExcludeChs=QStringList(), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi())
DSPSHARED_EXPORT Eigen::MatrixXd filterData(const Eigen::MatrixXd &matData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=1024, int designMethod=UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bUseThreads=true, bool bKeepOverhead=false)
DSPSHARED_EXPORT Eigen::MatrixXd filterDataBlock(const Eigen::MatrixXd &matData, const Eigen::RowVectorXi &vecPicks, const UTILSLIB::FilterKernel &filterKernel, bool bUseThreads=true)
DSPSHARED_EXPORT void filterChannel(FilterObject &channelDataTime)
DSPSHARED_EXPORT bool filterFile(QIODevice &pIODevice, QSharedPointer< FIFFLIB::FiffRawData > pFiffRawData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=4096, int designMethod=UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bUseThreads=true)
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
Lightweight filter configuration holding kernel coefficients and overlap-add state for one channel.
Definition rt_filter.h:82
UTILSLIB::FilterKernel filterKernel
Definition rt_filter.h:83
Eigen::RowVectorXd vecData
Definition rt_filter.h:85
Applies FIR filtering via FFT-based overlap-add convolution for continuous data streams.
Definition rt_filter.h:254
QSharedPointer< const FilterOverlapAdd > ConstSPtr
Definition rt_filter.h:257
Eigen::MatrixXd calculate(const Eigen::MatrixXd &matData, const UTILSLIB::FilterKernel &filterKernel, const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bFilterEnd=true, bool bUseThreads=true, bool bKeepOverhead=false)
Eigen::MatrixXd calculate(const Eigen::MatrixXd &matData, int type, double dCenterfreq, double dBandwidth, double dTransition, double dSFreq, int iOrder=1024, int designMethod=UTILSLIB::FilterKernel::m_designMethods.indexOf(UTILSLIB::FilterParameter("Cosine")), const Eigen::RowVectorXi &vecPicks=Eigen::RowVectorXi(), bool bFilterEnd=true, bool bUseThreads=true, bool bKeepOverhead=false)
QSharedPointer< FilterOverlapAdd > SPtr
Definition rt_filter.h:256
FIFF raw measurement data.