v2.0.0
Loading...
Searching...
No Matches
rtfiffrawviewmodel.h
Go to the documentation of this file.
1//=============================================================================================================
37
38#ifndef RTFIFFRAWVIEWMODEL_H
39#define RTFIFFRAWVIEWMODEL_H
40
41//=============================================================================================================
42// INCLUDES
43//=============================================================================================================
44
45#include "../../disp_global.h"
46
47#include <fiff/fiff_types.h>
48#include <fiff/fiff_proj.h>
49
50#include <dsp/filterkernel.h>
51
52//=============================================================================================================
53// QT INCLUDES
54//=============================================================================================================
55
56#include <QAbstractTableModel>
57#include <QSharedPointer>
58#include <QColor>
59
60#include <functional>
61#include <vector>
62
63//=============================================================================================================
64// EIGEN INCLUDES
65//=============================================================================================================
66
67#include <Eigen/Core>
68#include <Eigen/SparseCore>
69
70//=============================================================================================================
71// FORWARD DECLARATIONS
72//=============================================================================================================
73
74namespace FIFFLIB {
75 class FiffInfo;
76}
77
78//=============================================================================================================
79// DEFINE NAMESPACE DISPLIB
80//=============================================================================================================
81
82namespace DISPLIB
83{
84
85//=============================================================================================================
86// DISPLIB FORWARD DECLARATIONS
87//=============================================================================================================
88
89//=============================================================================================================
90// DEFINE TYPEDEFS
91//=============================================================================================================
92
93typedef QPair<const double*,qint32> RowVectorPair;
94typedef Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic,Eigen::RowMajor> MatrixXdR;
95
96//=============================================================================================================
102class DISPSHARED_EXPORT RtFiffRawViewModel : public QAbstractTableModel
103{
104 Q_OBJECT
105
106public:
107 typedef QSharedPointer<RtFiffRawViewModel> SPtr;
108 typedef QSharedPointer<const RtFiffRawViewModel> ConstSPtr;
109
110 //=========================================================================================================
116 RtFiffRawViewModel(QObject *parent = 0);
117
118 //=========================================================================================================
123
125 //=========================================================================================================
133 virtual int rowCount(const QModelIndex &parent = QModelIndex()) const ;
134
135 //=========================================================================================================
143 virtual int columnCount(const QModelIndex &parent = QModelIndex()) const;
144
145 //=========================================================================================================
154 virtual QVariant data(const QModelIndex &index, int role = Qt::DisplayRole) const;
155
156 //=========================================================================================================
166 virtual QVariant headerData(int section, Qt::Orientation orientation, int role = Qt::DisplayRole) const;
167
168 //=========================================================================================================
174 void setFiffInfo(QSharedPointer<FIFFLIB::FiffInfo>& p_pFiffInfo);
175
176 //=========================================================================================================
184 void setSamplingInfo(float sps, int T, bool bSetZero = false);
185
186 //=========================================================================================================
192 Eigen::MatrixXd getLastBlock();
193
194 //=========================================================================================================
200 void addData(const QList<Eigen::MatrixXd> &data);
201
202 //=========================================================================================================
210 FIFFLIB::fiff_int_t getKind(qint32 row) const;
211
212 //=========================================================================================================
220 FIFFLIB::fiff_int_t getUnit(qint32 row) const;
221
222 //=========================================================================================================
230 FIFFLIB::fiff_int_t getCoil(qint32 row) const;
231
232 //=========================================================================================================
238 inline qint32 getMaxSamples() const;
239
240 //=========================================================================================================
246 inline qint32 getCurrentSampleIndex() const;
247
248 //=========================================================================================================
256 inline double getLastBlockFirstValue(int row) const;
257
258 //=========================================================================================================
264 inline const QMap<qint32,qint32>& getIdxSelMap() const;
265
266 //=========================================================================================================
272 void selectRows(const QList<qint32> &selection);
273
274 //=========================================================================================================
280 void hideRows(const QList<qint32> &selection);
281
282 //=========================================================================================================
286 void resetSelection();
287
288 //=========================================================================================================
294 void toggleFreeze(const QModelIndex &index);
295
296 //=========================================================================================================
302 void setScaling(const QMap< qint32,float >& p_qMapChScaling);
303
304 //=========================================================================================================
310 void updateProjection(const QList<FIFFLIB::FiffProj>& projs);
311
312 //=========================================================================================================
318 void updateCompensator(int to);
319
320 //=========================================================================================================
326 void updateSpharaActivation(bool state);
327
328 //=========================================================================================================
336 void updateSpharaOptions(const QString& sSytemType, int nBaseFctsFirst, int nBaseFctsSecond);
337
338 //=========================================================================================================
344 void setFilter(QList<UTILSLIB::FilterKernel> filterData);
345
346 //=========================================================================================================
352 void setFilterActive(bool state);
353
354 //=========================================================================================================
360 void setBackgroundColor(const QColor& color);
361
362 //=========================================================================================================
368 void setFilterChannelType(const QString &channelType);
369
370 //=========================================================================================================
376 void createFilterChannelList(QStringList channelNames);
377
378 //=========================================================================================================
385 void markChBad(QModelIndex ch, bool status);
386
387 //=========================================================================================================
394 void markChBad(QModelIndexList chlist, bool status);
395
396 //=========================================================================================================
405 void triggerInfoChanged(const QMap<double, QColor>& colorMap, bool active, QString triggerCh, double threshold);
406
407 //=========================================================================================================
413 void distanceTimeSpacerChanged(int value);
414
415 //=========================================================================================================
419 void resetTriggerCounter();
420
421 //=========================================================================================================
427 inline qint32 numVLines() const;
428
429 //=========================================================================================================
435 inline bool isFreezed() const;
436
437 //=========================================================================================================
443 inline const QMap< qint32,float >& getScaling() const;
444
445 //=========================================================================================================
451 inline QList<QPair<int,double> > getDetectedTriggers() const;
452
453 //=========================================================================================================
459 inline QList<QPair<int, double> > getDetectedTriggersOld() const;
460
461 //=========================================================================================================
467 inline QMap<double, QColor> getTriggerColor() const;
468
469 //=========================================================================================================
475 inline int getNumberOfTimeSpacers() const;
476
477 //=========================================================================================================
483 inline double getTriggerThreshold() const;
484
485 //=========================================================================================================
491 inline QString getTriggerName() const;
492
493 //=========================================================================================================
499 inline int getCurrentTriggerIndex() const;
500
501 //=========================================================================================================
507 inline bool triggerDetectionActive() const;
508
509 //=========================================================================================================
515 inline int getCurrentOverlapAddDelay() const;
516
517 //=========================================================================================================
523 inline int getFirstSampleOffset() const;
524
525 //=========================================================================================================
532 double getMaxValueFromRawViewModel(int row) const;
533
534 //=========================================================================================================
540 void addEvent(int iSample);
541
542 //=========================================================================================================
551 std::vector<int> getEventsToDisplay(int iBegin, int iEnd) const;
552
553 //=========================================================================================================
560 void setEventCallbacks(std::function<void(int)> addFn,
561 std::function<std::vector<int>(int, int)> getFn);
562
563private:
564 //=========================================================================================================
568 void initSphara();
569
570 static void doFilterPerChannelRTMSA(QPair<QList<UTILSLIB::FilterKernel>,QPair<int,Eigen::RowVectorXd> > &channelDataTime);
571
572 //=========================================================================================================
576 void filterDataBlock();
577
578 //=========================================================================================================
585 void filterDataBlock(const Eigen::MatrixXd &data, int iDataIndex);
586
587 //=========================================================================================================
591 void clearModel();
592
593 bool m_bProjActivated;
594 bool m_bCompActivated;
595 bool m_bSpharaActivated;
596 bool m_bIsFreezed;
597 bool m_bDrawFilterFront;
598 bool m_bPerformFiltering;
599 bool m_bTriggerDetectionActive;
600 float m_fSps;
601 double m_dTriggerThreshold;
602 qint32 m_iT;
603 qint32 m_iDownsampling;
604 qint32 m_iMaxSamples;
605 qint32 m_iCurrentSample;
606 qint32 m_iCurrentStartingSample;
607 qint32 m_iCurrentSampleFreeze;
608 qint32 m_iMaxFilterLength;
609 qint32 m_iCurrentBlockSize;
610 qint32 m_iResidual;
611 int m_iCurrentTriggerChIndex;
612 int m_iDistanceTimerSpacer;
613 int m_iDetectedTriggers;
614
615 QString m_sCurrentTriggerCh;
616 QString m_sFilterChannelType;
617
618 QSharedPointer<FIFFLIB::FiffInfo> m_pFiffInfo;
619
620 Eigen::RowVectorXi m_vecBadIdcs;
621 Eigen::VectorXd m_vecLastBlockFirstValuesFiltered;
622 Eigen::VectorXd m_vecLastBlockFirstValuesRaw;
623
624 MatrixXdR m_matDataRaw;
625 MatrixXdR m_matDataFiltered;
626 MatrixXdR m_matDataRawFreeze;
627 MatrixXdR m_matDataFilteredFreeze;
628 Eigen::MatrixXd m_matOverlap;
629
630 Eigen::VectorXi m_vecIndicesFirstVV;
631 Eigen::VectorXi m_vecIndicesSecondVV;
632 Eigen::VectorXi m_vecIndicesFirstBabyMEG;
633 Eigen::VectorXi m_vecIndicesSecondBabyMEG;
634 Eigen::VectorXi m_vecIndicesFirstEEG;
635
636 Eigen::SparseMatrix<double> m_matSparseSpharaMult;
637 Eigen::SparseMatrix<double> m_matSparseProjCompMult;
638 Eigen::SparseMatrix<double> m_matSparseProjMult;
639 Eigen::SparseMatrix<double> m_matSparseCompMult;
640
641 Eigen::MatrixXd m_matProj;
642 Eigen::MatrixXd m_matComp;
643
644 Eigen::MatrixXd m_matSpharaVVGradLoaded;
645 Eigen::MatrixXd m_matSpharaVVMagLoaded;
646 Eigen::MatrixXd m_matSpharaBabyMEGInnerLoaded;
647 Eigen::MatrixXd m_matSpharaBabyMEGOuterLoaded;
648 Eigen::MatrixXd m_matSpharaEEGLoaded;
649
650 QMap<double, QColor> m_qMapTriggerColor;
651 QMap<int,QList<QPair<int,double> > >m_qMapDetectedTrigger;
652 QList<int> m_lTriggerChannelIndices;
653 QMap<int,QList<QPair<int,double> > >m_qMapDetectedTriggerFreeze;
654 QMap<int,QList<QPair<int,double> > >m_qMapDetectedTriggerOld;
655 QMap<int,QList<QPair<int,double> > >m_qMapDetectedTriggerOldFreeze;
656 QMap<qint32,float> m_qMapChScaling;
657 QList<UTILSLIB::FilterKernel>m_filterKernel;
658 QStringList m_filterChannelList;
659 QStringList m_visibleChannelList;
660 QMap<qint32,qint32> m_qMapIdxRowSelection;
661
662 QColor m_colBackground;
663
664 std::function<void(int)> m_fnAddEvent;
665 std::function<std::vector<int>(int, int)> m_fnGetEventSamples;
666
667signals:
668 //=========================================================================================================
674 void newSelection(const QList<qint32>& selection);
675
676 //=========================================================================================================
682 void windowSizeChanged(int windowSize);
683
684 //=========================================================================================================
688 void triggerDetected(int numberDetectedTriggers, const QMap<int,QList<QPair<int,double> > >& mapDetectedTriggers);
689};
690
691//=============================================================================================================
692// INLINE DEFINITIONS
693//=============================================================================================================
694
696{
697 return m_iMaxSamples;
698}
699
700//=============================================================================================================
701
703{
704 if(m_bIsFreezed) {
705 return m_iCurrentSampleFreeze;
706 }
707
708 if(!m_filterKernel.isEmpty() && m_bPerformFiltering) {
709 return m_iCurrentSample-m_iMaxFilterLength/2;
710 }
711
712 return m_iCurrentSample;
713}
714
715//=============================================================================================================
716
718{
719 if(row>m_vecLastBlockFirstValuesFiltered.rows() || row>m_vecLastBlockFirstValuesRaw.rows())
720 return 0;
721
722 if(!m_filterKernel.isEmpty())
723 return m_vecLastBlockFirstValuesFiltered[row];
724
725 return m_vecLastBlockFirstValuesRaw[row];
726}
727
728//=============================================================================================================
729
730inline const QMap<qint32,qint32>& RtFiffRawViewModel::getIdxSelMap() const
731{
732 return m_qMapIdxRowSelection;
733}
734
735//=============================================================================================================
736
737inline qint32 RtFiffRawViewModel::numVLines() const
738{
739 return (m_iT - 1);
740}
741
742//=============================================================================================================
743
745{
746 return m_bIsFreezed;
747}
748
749//=============================================================================================================
750
751inline const QMap<qint32,float>& RtFiffRawViewModel::getScaling() const
752{
753 return m_qMapChScaling;
754}
755
756//=============================================================================================================
757
758inline QMap<double, QColor> RtFiffRawViewModel::getTriggerColor() const
759{
760 if(m_bTriggerDetectionActive) {
761 return m_qMapTriggerColor;
762 }
763
764 QMap<double, QColor> map;
765 return map;
766}
767
768//=============================================================================================================
769
770inline QList<QPair<int,double> > RtFiffRawViewModel::getDetectedTriggers() const
771{
772 QList<QPair<int,double> > triggerIndices;
773
774 if(m_bIsFreezed)
775 return m_qMapDetectedTriggerFreeze[m_iCurrentTriggerChIndex];
776
777 if(m_bTriggerDetectionActive) {
778 return m_qMapDetectedTrigger[m_iCurrentTriggerChIndex];
779 }
780 else
781 return triggerIndices;
782}
783
784//=============================================================================================================
785
786inline QList<QPair<int,double> > RtFiffRawViewModel::getDetectedTriggersOld() const
787{
788 QList<QPair<int,double> > triggerIndices;
789
790 if(m_bIsFreezed)
791 return m_qMapDetectedTriggerOldFreeze[m_iCurrentTriggerChIndex];
792
793 if(m_bTriggerDetectionActive) {
794 return m_qMapDetectedTriggerOld[m_iCurrentTriggerChIndex];
795 }
796 else
797 return triggerIndices;
798}
799
800//=============================================================================================================
801
803{
804 //qDebug()<<((m_iT*1000)/m_iDistanceTimerSpacer)-1;
805 return ((1000)/m_iDistanceTimerSpacer)-1;
806}
807
808//=============================================================================================================
809
811{
812 return m_dTriggerThreshold;
813}
814
815//=============================================================================================================
816
818{
819 return m_sCurrentTriggerCh;
820}
821
822//=============================================================================================================
823
825{
826 return m_iCurrentTriggerChIndex;
827}
828
829//=============================================================================================================
830
832{
833 return m_bTriggerDetectionActive;
834}
835
836//=============================================================================================================
837
839{
840 if(!m_filterKernel.isEmpty())
841 return m_iMaxFilterLength/2;
842 else
843 return 0;
844}
845
846//=============================================================================================================
847
849{
850 return m_iCurrentStartingSample;
851}
852} // NAMESPACE
853
854#if QT_VERSION < QT_VERSION_CHECK(6, 0, 0)
855#ifndef metatype_rowvectorpair
856#define metatype_rowvectorpair
858#endif
859#endif
860
861#endif // RTFIFFRAWVIEWMODEL_H
disp library export/import macros.
#define DISPSHARED_EXPORT
Definition disp_global.h:51
FiffProj class declaration.
Old fiff_type declarations - replace them.
The FilterKernel class represents a filter object that generates the FIR filter coefficients using Pa...
Q_DECLARE_METATYPE(Eigen::MatrixXf)
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
2-D display widgets and visualisation helpers (charts, topography, colour maps).
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > MatrixXdR
QPair< const double *, qint32 > RowVectorPair
double getLastBlockFirstValue(int row) const
const QMap< qint32, float > & getScaling() const
void setBackgroundColor(const QColor &color)
void setFilter(QList< UTILSLIB::FilterKernel > filterData)
void selectRows(const QList< qint32 > &selection)
void newSelection(const QList< qint32 > &selection)
FIFFLIB::fiff_int_t getKind(qint32 row) const
QSharedPointer< RtFiffRawViewModel > SPtr
QSharedPointer< const RtFiffRawViewModel > ConstSPtr
void setSamplingInfo(float sps, int T, bool bSetZero=false)
void markChBad(QModelIndex ch, bool status)
RtFiffRawViewModel(QObject *parent=0)
void windowSizeChanged(int windowSize)
void createFilterChannelList(QStringList channelNames)
void setFilterChannelType(const QString &channelType)
void setFiffInfo(QSharedPointer< FIFFLIB::FiffInfo > &p_pFiffInfo)
virtual int rowCount(const QModelIndex &parent=QModelIndex()) const
QList< QPair< int, double > > getDetectedTriggers() const
void toggleFreeze(const QModelIndex &index)
void hideRows(const QList< qint32 > &selection)
void setEventCallbacks(std::function< void(int)> addFn, std::function< std::vector< int >(int, int)> getFn)
void updateSpharaOptions(const QString &sSytemType, int nBaseFctsFirst, int nBaseFctsSecond)
void setScaling(const QMap< qint32, float > &p_qMapChScaling)
virtual int columnCount(const QModelIndex &parent=QModelIndex()) const
FIFFLIB::fiff_int_t getCoil(qint32 row) const
virtual QVariant headerData(int section, Qt::Orientation orientation, int role=Qt::DisplayRole) const
const QMap< qint32, qint32 > & getIdxSelMap() const
void triggerDetected(int numberDetectedTriggers, const QMap< int, QList< QPair< int, double > > > &mapDetectedTriggers)
QMap< double, QColor > getTriggerColor() const
FIFFLIB::fiff_int_t getUnit(qint32 row) const
std::vector< int > getEventsToDisplay(int iBegin, int iEnd) const
void triggerInfoChanged(const QMap< double, QColor > &colorMap, bool active, QString triggerCh, double threshold)
void updateProjection(const QList< FIFFLIB::FiffProj > &projs)
virtual QVariant data(const QModelIndex &index, int role=Qt::DisplayRole) const
QList< QPair< int, double > > getDetectedTriggersOld() const
void addData(const QList< Eigen::MatrixXd > &data)
double getMaxValueFromRawViewModel(int row) const
FIFF measurement file information.
Definition fiff_info.h:86