v2.0.0
Loading...
Searching...
No Matches
rt_filter.cpp
Go to the documentation of this file.
1//=============================================================================================================
35
36//=============================================================================================================
37// INCLUDES
38//=============================================================================================================
39
40#include "rt_filter.h"
41
42#include <fiff/fiff_raw_data.h>
43#include <fiff/fiff_file.h>
44#include <mne/mne_epoch_data.h>
46
47//=============================================================================================================
48// QT INCLUDES
49//=============================================================================================================
50
51#include <QDebug>
52
53//=============================================================================================================
54// EIGEN INCLUDES
55//=============================================================================================================
56
57#include <Eigen/Dense>
58#include <Eigen/Core>
59
60//=============================================================================================================
61// USED NAMESPACES
62//=============================================================================================================
63
64using namespace RTPROCESSINGLIB;
65using namespace Eigen;
66using namespace FIFFLIB;
67using namespace UTILSLIB;
68using namespace MNELIB;
69
70//=============================================================================================================
71// DEFINE GLOBAL RTPROCESSINGLIB METHODS
72//=============================================================================================================
73
74bool RTPROCESSINGLIB::filterFile(QIODevice &pIODevice,
75 QSharedPointer<FiffRawData> pFiffRawData,
76 int type,
77 double dCenterfreq,
78 double bandwidth,
79 double dTransition,
80 double dSFreq,
81 int iOrder,
82 int designMethod,
83 const RowVectorXi& vecPicks,
84 bool bUseThreads)
85{
86 // Normalize cut off frequencies to nyquist
87 dCenterfreq = dCenterfreq/(dSFreq/2.0);
88 bandwidth = bandwidth/(dSFreq/2.0);
89 dTransition = dTransition/(dSFreq/2.0);
90
91 // create filter
92 FilterKernel filter = FilterKernel("filter_kernel",
93 type,
94 iOrder,
95 dCenterfreq,
96 bandwidth,
97 dTransition,
98 dSFreq,
99 designMethod);
100
101 return filterFile(pIODevice,
102 pFiffRawData,
103 filter,
104 vecPicks,
105 bUseThreads);
106}
107
108//=============================================================================================================
109
110bool RTPROCESSINGLIB::filterFile(QIODevice &pIODevice,
111 QSharedPointer<FiffRawData> pFiffRawData,
112 const FilterKernel& filterKernel,
113 const RowVectorXi& vecPicks,
114 bool bUseThreads)
115{
116 int iOrder = filterKernel.getFilterOrder();
117
118 RowVectorXd cals;
119 SparseMatrix<double> mult;
120 RowVectorXi sel;
121 FiffStream::SPtr outfid = FiffStream::start_writing_raw(pIODevice, pFiffRawData->info, cals);
122
123 //Setup reading parameters
124 fiff_int_t from = pFiffRawData->first_samp;
125 fiff_int_t to = pFiffRawData->last_samp;
126
127 // slice input data into data junks with proper length so that the slices are always >= the filter order
128 float fFactor = 2.0f;
129 int iSize = fFactor * iOrder;
130 int residual = (to - from) % iSize;
131 while(residual < iOrder) {
132 fFactor = fFactor - 0.1f;
133 iSize = fFactor * iOrder;
134 residual = (to - from) % iSize;
135
136 if((iSize < iOrder)) {
137 qInfo() << "[Filter::filterData] Sliced data block size is too small. Filtering whole block at once.";
138 iSize = to - from;
139 break;
140 }
141 }
142
143 float quantum_sec = iSize/pFiffRawData->info.sfreq;
144 fiff_int_t quantum = ceil(quantum_sec*pFiffRawData->info.sfreq);
145
146 // Read, filter and write the data
147 bool first_buffer = true;
148
149 fiff_int_t first, last;
150 MatrixXd matData, matDataOverlap;
151 MatrixXd times;
152
153 for(first = from; first < to; first+=quantum) {
154 last = first+quantum-1;
155 if (last > to) {
156 last = to;
157 }
158
159 if (!pFiffRawData->read_raw_segment(matData, times, mult, first, last, sel)) {
160 qWarning("[Filter::filterData] Error during read_raw_segment\n");
161 return false;
162 }
163
164 qInfo() << "Filtering and writing block" << first << "to" << last;
165
166 if (first_buffer) {
167 if (first > 0) {
168 outfid->write_int(FIFF_FIRST_SAMPLE,&first);
169 }
170 first_buffer = false;
171 }
172
173 matData = filterDataBlock(matData,
174 vecPicks,
175 filterKernel,
176 bUseThreads);
177
178 if(first == from) {
179 outfid->write_raw_buffer(matData.block(0,iOrder/2,matData.rows(),matData.cols()-iOrder), cals);
180 } else if(first + quantum >= to) {
181 matData.block(0,0,matData.rows(),iOrder) += matDataOverlap;
182 outfid->write_raw_buffer(matData.block(0,0,matData.rows(),matData.cols()-iOrder), cals);
183 } else {
184 matData.block(0,0,matData.rows(),iOrder) += matDataOverlap;
185 outfid->write_raw_buffer(matData.block(0,0,matData.rows(),matData.cols()-iOrder), cals);
186 }
187
188 matDataOverlap = matData.block(0,matData.cols()-iOrder,matData.rows(),iOrder);
189 }
190
191 outfid->finish_writing_raw();
192
193 return true;
194}
195
196//=============================================================================================================
197
198MatrixXd RTPROCESSINGLIB::filterData(const MatrixXd& matData,
199 int type,
200 double dCenterfreq,
201 double bandwidth,
202 double dTransition,
203 double dSFreq,
204 int iOrder,
205 int designMethod,
206 const RowVectorXi& vecPicks,
207 bool bUseThreads,
208 bool bKeepOverhead)
209{
210 // Check for size of data
211 if(matData.cols() < iOrder){
212 qWarning() << QString("[Filter::filterData] Filter length/order is bigger than data length. Returning.");
213 return matData;
214 }
215
216 // Normalize cut off frequencies to nyquist
217 dCenterfreq = dCenterfreq/(dSFreq/2.0);
218 bandwidth = bandwidth/(dSFreq/2.0);
219 dTransition = dTransition/(dSFreq/2.0);
220
221 // create filter
222 FilterKernel filter = FilterKernel("filter_kernel",
223 type,
224 iOrder,
225 dCenterfreq,
226 bandwidth,
227 dTransition,
228 dSFreq,
229 designMethod);
230
231 return filterData(matData,
232 filter,
233 vecPicks,
234 bUseThreads,
235 bKeepOverhead);
236}
237
238//=============================================================================================================
239
240MatrixXd RTPROCESSINGLIB::filterData(const MatrixXd& matData,
241 const FilterKernel& filterKernel,
242 const RowVectorXi& vecPicks,
243 bool bUseThreads,
244 bool bKeepOverhead)
245{
246 int iOrder = filterKernel.getFilterOrder();
247
248 // Check for size of data
249 if(matData.cols() < iOrder){
250 qWarning() << "[Filter::filterData] Filter length/order is bigger than data length. Returning.";
251 return matData;
252 }
253
254 // Create output matrix with size of input matrix
255 MatrixXd matDataOut(matData.rows(), matData.cols()+iOrder);
256 matDataOut.setZero();
257 MatrixXd sliceFiltered;
258
259 // slice input data into data junks with proper length so that the slices are always >= the filter order
260 float fFactor = 2.0f;
261 int iSize = fFactor * iOrder;
262 int residual = matData.cols() % iSize;
263 while(residual < iOrder) {
264 fFactor = fFactor - 0.1f;
265 iSize = fFactor * iOrder;
266 residual = matData.cols() % iSize;
267
268 if(iSize < iOrder) {
269 iSize = matData.cols();
270 break;
271 }
272 }
273
274 if(matData.cols() > iSize) {
275 int from = 0;
276 int numSlices = ceil(float(matData.cols())/float(iSize)); //calculate number of data slices
277
278 for (int i = 0; i < numSlices; i++) {
279 if(i == numSlices-1) {
280 //catch the last one that might be shorter than the other blocks
281 iSize = matData.cols() - (iSize * (numSlices -1));
282 }
283
284 // Filter the data block. This will return data with a fitler delay of iOrder/2 in front and back
285 sliceFiltered = filterDataBlock(matData.block(0,from,matData.rows(),iSize),
286 vecPicks,
287 filterKernel,
288 bUseThreads);
289
290 // Perform overlap add
291 if(i == 0) {
292 matDataOut.block(0,0,matData.rows(),sliceFiltered.cols()) += sliceFiltered;
293 } else {
294 matDataOut.block(0,from,matData.rows(),sliceFiltered.cols()) += sliceFiltered;
295 }
296
297 from += iSize;
298 }
299 } else {
300 matDataOut = filterDataBlock(matData,
301 vecPicks,
302 filterKernel,
303 bUseThreads);
304 }
305
306 if(bKeepOverhead) {
307 return matDataOut;
308 } else {
309 return matDataOut.block(0,iOrder/2,matDataOut.rows(),matData.cols());
310 }
311}
312
313//=============================================================================================================
314
315MatrixXd RTPROCESSINGLIB::filterDataBlock(const MatrixXd& matData,
316 const RowVectorXi& vecPicks,
317 const FilterKernel& filterKernel,
318 bool bUseThreads)
319{
320 int iOrder = filterKernel.getFilterOrder();
321
322 // Check for size of data
323 if(matData.cols() < iOrder){
324 qWarning() << QString("[Filter::filterDataBlock] Filter length/order is bigger than data length. Returning.");
325 return matData;
326 }
327
328 // Setup filters to the correct length, so we do not have to do this everytime we call the FFT filter function
329 FilterKernel filterKernelSetup = filterKernel;
330 filterKernelSetup.prepareFilter(matData.cols());
331
332 // Do the concurrent filtering
333 RowVectorXi vecPicksNew = vecPicks;
334 if(vecPicksNew.cols() == 0) {
335 vecPicksNew = RowVectorXi::LinSpaced(matData.rows(), 0, matData.rows());
336 }
337
338 // Generate QList structure which can be handled by the QConcurrent framework
339 QList<FilterObject> timeData;
340
341 // Only select channels specified in vecPicksNew
342 FilterObject data;
343 for(qint32 i = 0; i < vecPicksNew.cols(); ++i) {
344 data.filterKernel = filterKernelSetup;
345 data.iRow = vecPicksNew[i];
346 data.vecData = matData.row(vecPicksNew[i]);
347 timeData.append(data);
348 }
349
350 // Copy in data from last data block. This is necessary in order to also delay channels which are not filtered
351 MatrixXd matDataOut(matData.rows(), matData.cols()+iOrder);
352 matDataOut.setZero();
353 matDataOut.block(0, iOrder/2, matData.rows(), matData.cols()) = matData;
354
355 if(bUseThreads) {
356 QFuture<void> future = QtConcurrent::map(timeData,
358 future.waitForFinished();
359 } else {
360 for(int i = 0; i < timeData.size(); ++i) {
361 filterChannel(timeData[i]);
362 }
363 }
364
365 // Do the overlap add method and store in matDataOut
366 RowVectorXd tempData;
367
368 for(int r = 0; r < timeData.size(); r++) {
369 // Write the newly calculated filtered data to the filter data matrix. This data has a delay of iOrder/2 in front and back
370 matDataOut.row(timeData.at(r).iRow) = timeData.at(r).vecData;
371 }
372
373 return matDataOut;
374}
375
376//=============================================================================================================
377
379{
380 //channelDataTime.vecData = channelDataTime.first.at(i).applyConvFilter(channelDataTime.vecData, true);
381 channelDataTime.filterKernel.applyFftFilter(channelDataTime.vecData, true); //FFT Convolution for rt is not suitable. FFT make the signal filtering non causal.
382}
383
384//=============================================================================================================
385// DEFINE MEMBER METHODS
386//=============================================================================================================
387
388MatrixXd FilterOverlapAdd::calculate(const MatrixXd& matData,
389 int type,
390 double dCenterfreq,
391 double bandwidth,
392 double dTransition,
393 double dSFreq,
394 int iOrder,
395 int designMethod,
396 const RowVectorXi& vecPicks,
397 bool bFilterEnd,
398 bool bUseThreads,
399 bool bKeepOverhead)
400{
401 // Check for size of data
402 if(matData.cols() < iOrder){
403 qWarning() << QString("[Filter::filterData] Filter length/order is bigger than data length. Returning.");
404 return matData;
405 }
406
407 // Normalize cut off frequencies to nyquist
408 dCenterfreq = dCenterfreq/(dSFreq/2.0);
409 bandwidth = bandwidth/(dSFreq/2.0);
410 dTransition = dTransition/(dSFreq/2.0);
411
412 // create filter
413 FilterKernel filter = FilterKernel("filter_kernel",
414 type,
415 iOrder,
416 dCenterfreq,
417 bandwidth,
418 dTransition,
419 dSFreq,
420 designMethod);
421
422 return calculate(matData,
423 filter,
424 vecPicks,
425 bFilterEnd,
426 bUseThreads,
427 bKeepOverhead);
428}
429
430//=============================================================================================================
431
432MatrixXd FilterOverlapAdd::calculate(const MatrixXd& matData,
433 const FilterKernel& filterKernel,
434 const RowVectorXi& vecPicks,
435 bool bFilterEnd,
436 bool bUseThreads,
437 bool bKeepOverhead)
438{
439 int iOrder = filterKernel.getFilterOrder();
440
441 // Check for size of data
442 if(matData.cols() < iOrder){
443 qWarning() << "[Filter::filterData] Filter length/order is bigger than data length. Returning.";
444 return matData;
445 }
446
447 // Init overlaps from last block
448 if(m_matOverlapBack.cols() != iOrder || m_matOverlapBack.rows() < matData.rows()) {
449 m_matOverlapBack.resize(matData.rows(), iOrder);
450 m_matOverlapBack.setZero();
451 }
452
453 if(m_matOverlapFront.cols() != iOrder || m_matOverlapFront.rows() < matData.rows()) {
454 m_matOverlapFront.resize(matData.rows(), iOrder);
455 m_matOverlapFront.setZero();
456 }
457
458 // Create output matrix with size of input matrix
459 MatrixXd matDataOut(matData.rows(), matData.cols()+iOrder);
460 matDataOut.setZero();
461 MatrixXd sliceFiltered;
462
463 // slice input data into data junks with proper length so that the slices are always >= the filter order
464 float fFactor = 2.0f;
465 int iSize = fFactor * iOrder;
466 int residual = matData.cols() % iSize;
467 while(residual < iOrder) {
468 fFactor = fFactor - 0.1f;
469 iSize = fFactor * iOrder;
470 residual = matData.cols() % iSize;
471
472 if(iSize < iOrder) {
473 iSize = matData.cols();
474 break;
475 }
476 }
477
478 if(matData.cols() > iSize) {
479 int from = 0;
480 int numSlices = ceil(float(matData.cols())/float(iSize)); //calculate number of data slices
481
482 for (int i = 0; i < numSlices; i++) {
483 if(i == numSlices-1) {
484 //catch the last one that might be shorter than the other blocks
485 iSize = matData.cols() - (iSize * (numSlices -1));
486 }
487
488 // Filter the data block. This will return data with a fitler delay of iOrder/2 in front and back
489 sliceFiltered = filterDataBlock(matData.block(0,from,matData.rows(),iSize),
490 vecPicks,
491 filterKernel,
492 bUseThreads);
493
494 if(i == 0) {
495 matDataOut.block(0,0,matData.rows(),sliceFiltered.cols()) += sliceFiltered;
496 } else {
497 matDataOut.block(0,from,matData.rows(),sliceFiltered.cols()) += sliceFiltered;
498 }
499
500 if(bFilterEnd && (i == 0)) {
501 matDataOut.block(0,0,matDataOut.rows(),iOrder) += m_matOverlapBack;
502 } else if (!bFilterEnd && (i == numSlices-1)) {
503 matDataOut.block(0,matDataOut.cols()-iOrder,matDataOut.rows(),iOrder) += m_matOverlapFront;
504 }
505
506 from += iSize;
507 }
508 } else {
509 matDataOut = filterDataBlock(matData,
510 vecPicks,
511 filterKernel,
512 bUseThreads);
513
514 if(bFilterEnd) {
515 matDataOut.block(0,0,matDataOut.rows(),iOrder) += m_matOverlapBack;
516 } else {
517 matDataOut.block(0,matDataOut.cols()-iOrder,matDataOut.rows(),iOrder) += m_matOverlapFront;
518 }
519 }
520
521 // Refresh the overlap matrix with the new calculated filtered data
522 m_matOverlapBack = matDataOut.block(0,matDataOut.cols()-iOrder,matDataOut.rows(),iOrder);
523 m_matOverlapFront = matDataOut.block(0,0,matDataOut.rows(),iOrder);
524
525 if(bKeepOverhead) {
526 return matDataOut;
527 } else {
528 return matDataOut.block(0,0,matDataOut.rows(),matData.cols());
529 }
530}
531
532//=============================================================================================================
533
535{
536 m_matOverlapBack.resize(0,0);
537 m_matOverlapFront.resize(0,0);
538}
539
540//=============================================================================================================
541
543 const MatrixXi& matEvents,
544 float fTMinS,
545 float fTMaxS,
546 qint32 eventType,
547 bool bApplyBaseline,
548 float fTBaselineFromS,
549 float fTBaselineToS,
550 const QMap<QString,double>& mapReject,
551 const FilterKernel& filterKernel,
552 const QStringList& lExcludeChs,
553 const RowVectorXi& picks)
554{
555 MNEEpochDataList lstEpochDataList;
556
557 // Select the desired events
558 qint32 count = 0;
559 qint32 p;
560 MatrixXi selected = MatrixXi::Zero(1, matEvents.rows());
561 for (p = 0; p < matEvents.rows(); ++p)
562 {
563 if (matEvents(p,1) == 0 && matEvents(p,2) == eventType)
564 {
565 selected(0,count) = p;
566 ++count;
567 }
568 }
569 selected.conservativeResize(1, count);
570 if (count > 0) {
571 qInfo("[RTPROCESSINGLIB::computeFilteredAverage] %d matching events found",count);
572 }
573
574 // If picks are empty, pick all
575 RowVectorXi picksNew = picks;
576 if(picks.cols() <= 0) {
577 picksNew.resize(raw.info.chs.size());
578 for(int i = 0; i < raw.info.chs.size(); ++i) {
579 picksNew(i) = i;
580 }
581 }
582
583 fiff_int_t event_samp, from, to;
584 fiff_int_t dropCount = 0;
585 MatrixXd timesDummy;
586 MatrixXd times;
587
588 std::unique_ptr<MNEEpochData> epoch;
589 int iFilterDelay = filterKernel.getFilterOrder()/2;
590
591 for (p = 0; p < count; ++p) {
592 // Read a data segment
593 event_samp = matEvents(selected(p),0);
594 from = event_samp + fTMinS*raw.info.sfreq;
595 to = event_samp + floor(fTMaxS*raw.info.sfreq + 0.5);
596
597 epoch = std::make_unique<MNEEpochData>();
598
599 if(raw.read_raw_segment(epoch->epoch, timesDummy, from - iFilterDelay, to + iFilterDelay, picksNew)) {
600 // Filter the data
601 epoch->epoch = RTPROCESSINGLIB::filterData(epoch->epoch,filterKernel).block(0, iFilterDelay, epoch->epoch.rows(), to-from);
602
603 if (p == 0) {
604 times.resize(1, to-from+1);
605 for (qint32 i = 0; i < times.cols(); ++i)
606 times(0, i) = ((float)(from-event_samp+i)) / raw.info.sfreq;
607 }
608
609 epoch->event = eventType;
610 epoch->tmin = fTMinS;
611 epoch->tmax = fTMaxS;
612
613 epoch->bReject = MNEEpochDataList::checkForArtifact(epoch->epoch,
614 raw.info,
615 mapReject,
616 lExcludeChs);
617
618 if (epoch->bReject) {
619 dropCount++;
620 }
621
622 //Check if data block has the same size as the previous one
623 if(!lstEpochDataList.isEmpty()) {
624 if(epoch->epoch.size() == lstEpochDataList.last()->epoch.size()) {
625 lstEpochDataList.append(MNEEpochData::SPtr(epoch.release()));
626 }
627 } else {
628 lstEpochDataList.append(MNEEpochData::SPtr(epoch.release()));
629 }
630 } else {
631 qWarning("[MNEEpochDataList::readEpochs] Can't read the event data segments.");
632 }
633 }
634
635 qInfo().noquote() << "[MNEEpochDataList::readEpochs] Read a total of"<< lstEpochDataList.size() <<"epochs of type" << eventType << "and marked"<< dropCount <<"for rejection.";
636
637 if(bApplyBaseline) {
638 QPair<float, float> baselinePair(fTBaselineFromS, fTBaselineToS);
639 lstEpochDataList.applyBaselineCorrection(baselinePair);
640 }
641
642 if(!mapReject.isEmpty()) {
643 lstEpochDataList.dropRejected();
644 }
645
646 return lstEpochDataList.average(raw.info,
647 0,
648 lstEpochDataList.first()->epoch.cols());
649}
FiffRawData class declaration.
Header file describing the numerical values used in fif files.
#define FIFF_FIRST_SAMPLE
Definition fiff_file.h:461
MNEEpochDataList class declaration.
MNEEpochData class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
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)
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
The FilterKernel class provides methods to create/design a FIR filter kernel.
void applyFftFilter(Eigen::RowVectorXd &vecData, bool bKeepOverhead=false)
void prepareFilter(int iDataSize)
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
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)
QList< FiffChInfo > chs
FIFF raw measurement data.
bool read_raw_segment(Eigen::MatrixXd &data, Eigen::MatrixXd &times, fiff_int_t from=-1, fiff_int_t to=-1, const Eigen::RowVectorXi &sel=defaultRowVectorXi, bool do_debug=false) const
QSharedPointer< FiffStream > SPtr
static FiffStream::SPtr start_writing_raw(QIODevice &p_IODevice, const FiffInfo &info, Eigen::RowVectorXd &cals, Eigen::MatrixXi sel=defaultMatrixXi, bool bResetRange=true)
QSharedPointer< MNEEpochData > SPtr
FIFFLIB::FiffEvoked average(const FIFFLIB::FiffInfo &p_info, FIFFLIB::fiff_int_t first, FIFFLIB::fiff_int_t last, Eigen::VectorXi sel=FIFFLIB::defaultVectorXi, bool proj=false)
void applyBaselineCorrection(const QPair< float, float > &baseline)
static bool checkForArtifact(const Eigen::MatrixXd &data, const FIFFLIB::FiffInfo &pFiffInfo, const QMap< QString, double > &mapReject, const QStringList &lExcludeChs=QStringList())