MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
filterkernel.h
Go to the documentation of this file.
1//=============================================================================================================
47#ifndef FILTERKERNEL_H
48#define FILTERKERNEL_H
49
50//=============================================================================================================
51// INCLUDES
52//=============================================================================================================
53
54#include "../rtprocessing_global.h"
55//#include "filter.h"
56
57//=============================================================================================================
58// EIGEN INCLUDES
59//=============================================================================================================
60
61#include <Eigen/Core>
62
63//=============================================================================================================
64// QT INCLUDES
65//=============================================================================================================
66
67#include <QString>
68#include <QMetaType>
69#include <QVector>
70#include <QDebug>
71
72//=============================================================================================================
73// DEFINE NAMESPACE RTPROCESSINGLIB
74//=============================================================================================================
75
76namespace RTPROCESSINGLIB
77{
78
79//=============================================================================================================
84
85public:
86 //=========================================================================================================
90 explicit FilterParameter();
91
92 //=========================================================================================================
98 explicit FilterParameter(QString sName);
99
100 //=========================================================================================================
107 explicit FilterParameter(QString sName, QString sDescription);
108
109 //=========================================================================================================
115 QString getName() const;
116
117 friend bool operator == (const FilterParameter& in1, const FilterParameter& in2){
118 //qDebug() << in1.getName() << in2.getName();
119 return (in1.getName() == in2.getName());
120 }
121protected:
122 QString m_sName;
124};
125
126//=============================================================================================================
133{
134
135public:
136 //=========================================================================================================
140 FilterKernel();
141
142 //=========================================================================================================
155 FilterKernel(const QString &sFilterName,
156 int iFilterType,
157 int iOrder,
158 double dCenterfreq,
159 double dBandwidth,
160 double dParkswidth,
161 double dSFreq,
162 int iDesignMethod);
163
164 //=========================================================================================================
174 void prepareFilter(int iDataSize);
175
176 //=========================================================================================================
186 Eigen::RowVectorXd applyConvFilter(const Eigen::RowVectorXd& vecData,
187 bool bKeepOverhead = false) const;
188
189 //=========================================================================================================
199 void applyFftFilter(Eigen::RowVectorXd& vecData,
200 bool bKeepOverhead = false);
201
202 QString getName() const;
203 void setName(const QString& sFilterName);
204
205 double getSamplingFrequency() const;
206 void setSamplingFrequency(double dSFreq);
207
208 int getFilterOrder() const;
209 void setFilterOrder(int iOrder);
210
211 double getCenterFrequency() const;
212 void setCenterFrequency(double dCenterFreq);
213
214 double getBandwidth() const;
215 void setBandwidth(double dBandwidth);
216
217 double getParksWidth() const;
218 void setParksWidth(double dParksWidth);
219
220 double getHighpassFreq() const;
221 void setHighpassFreq(double dHighpassFreq);
222
223 double getLowpassFreq() const;
224 void setLowpassFreq(double dLowpassFreq);
225
226 Eigen::RowVectorXd getCoefficients() const;
227 void setCoefficients(const Eigen::RowVectorXd& vecCoeff);
228
229 Eigen::RowVectorXcd getFftCoefficients() const;
230 void setFftCoefficients(const Eigen::RowVectorXcd& vecFftCoeff);
231
232 FilterParameter getDesignMethod() const;
233 void setDesignMethod(int iDesignMethod);
234
235 FilterParameter getFilterType() const;
236 void setFilterType(int iFilterType);
237
238 QString getShortDescription() const;
239
240 static QVector<FilterParameter> m_designMethods;
241 static QVector<FilterParameter> m_filterTypes;
243private:
244 //=========================================================================================================
250 bool fftTransformCoeffs(int iFftLength);
251
252 //=========================================================================================================
256 void designFilter();
257
258 double m_sFreq;
259 double m_dCenterFreq;
260 double m_dBandwidth;
261 double m_dParksWidth;
262 double m_dLowpassFreq;
263 double m_dHighpassFreq;
265 int m_iFilterOrder;
266 int m_iDesignMethod;
267 int m_iFilterType;
269 QString m_sFilterName;
270 QString m_sFilterShortDescription;
272 Eigen::RowVectorXd m_vecCoeff;
273 Eigen::RowVectorXcd m_vecFftCoeff;
274};
275
276} // NAMESPACE RTPROCESSINGLIB
277
278#ifndef metatype_filterkernel
279#define metatype_filterkernel
281#endif
282
283#ifndef metatype_filterparameter
284#define metatype_filterkernel
286#endif
287
288#endif // FILTERKERNEL_H
#define RTPROCESINGSHARED_EXPORT
Q_DECLARE_METATYPE(Eigen::MatrixXf)
The FilterParameter class.
The FilterKernel class provides methods to create/design a FIR filter kernel.
static QVector< FilterParameter > m_designMethods
static QVector< FilterParameter > m_filterTypes