v2.0.0
Loading...
Searching...
No Matches
rtsourcedatacontroller.h
Go to the documentation of this file.
1//=============================================================================================================
34
35#ifndef BRAINVIEW_RTSOURCEDATACONTROLLER_H
36#define BRAINVIEW_RTSOURCEDATACONTROLLER_H
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
43
44#include <QObject>
45#include <QSharedPointer>
46#include <QVector>
47#include <QList>
48#include <QString>
49#include <vector>
50#include <Eigen/Core>
51#include <Eigen/SparseCore>
52
53//=============================================================================================================
54// FORWARD DECLARATIONS
55//=============================================================================================================
56
57class QThread;
58class QTimer;
59
60namespace FSLIB {
61 class Label;
62}
63
64namespace DISP3DRHILIB {
67}
68
69//=============================================================================================================
100{
101 Q_OBJECT
102
103public:
104 //=========================================================================================================
110 explicit RtSourceDataController(QObject *parent = nullptr);
111
112 //=========================================================================================================
116 ~RtSourceDataController() override;
117
118 //=========================================================================================================
125 void addData(const Eigen::VectorXd &data);
126
127 //=========================================================================================================
133 void setInterpolationMatrixLeft(QSharedPointer<Eigen::SparseMatrix<float>> mat);
134
135 //=========================================================================================================
141 void setInterpolationMatrixRight(QSharedPointer<Eigen::SparseMatrix<float>> mat);
142
143 //=========================================================================================================
149 void setStreamingState(bool state);
150
151 //=========================================================================================================
157 bool isStreaming() const;
158
159 //=========================================================================================================
165 void setTimeInterval(int msec);
166
167 //=========================================================================================================
173 void setNumberAverages(int numAvr);
174
175 //=========================================================================================================
181 void setColormapType(const QString &name);
182
183 //=========================================================================================================
191 void setThresholds(double min, double mid, double max);
192
193 //=========================================================================================================
199 void setLoopState(bool enabled);
200
201 //=========================================================================================================
207 void setSFreq(double sFreq);
208
209 //=========================================================================================================
213 void clearData();
214
215 //=========================================================================================================
224 void setSurfaceColor(const QVector<uint32_t> &baseColorsLh,
225 const QVector<uint32_t> &baseColorsRh);
226
227 //=========================================================================================================
234 void setStreamSmoothedData(bool bStreamSmoothedData);
235
236 //=========================================================================================================
237 // ── On-the-fly interpolation matrix computation ─────────────────────
238 //=========================================================================================================
239
246 void setInterpolationFunction(const QString &sInterpolationFunction);
247
248 //=========================================================================================================
255 void setCancelDistance(double dCancelDist);
256
257 //=========================================================================================================
266 void setInterpolationInfoLeft(const Eigen::MatrixX3f &matVertices,
267 const std::vector<Eigen::VectorXi> &vecNeighborVertices,
268 const Eigen::VectorXi &vecSourceVertices);
269
270 //=========================================================================================================
279 void setInterpolationInfoRight(const Eigen::MatrixX3f &matVertices,
280 const std::vector<Eigen::VectorXi> &vecNeighborVertices,
281 const Eigen::VectorXi &vecSourceVertices);
282
283 //=========================================================================================================
289
290 //=========================================================================================================
296 void setVisualizationType(int iVisType);
297
298 //=========================================================================================================
306 void setAnnotationInfoLeft(const Eigen::VectorXi &vecLabelIds,
307 const QList<FSLIB::Label> &lLabels,
308 const Eigen::VectorXi &vecVertNo);
309
310 //=========================================================================================================
318 void setAnnotationInfoRight(const Eigen::VectorXi &vecLabelIds,
319 const QList<FSLIB::Label> &lLabels,
320 const Eigen::VectorXi &vecVertNo);
321
322signals:
323 //=========================================================================================================
330 void newSmoothedDataAvailable(const QVector<uint32_t> &colorsLh,
331 const QVector<uint32_t> &colorsRh);
332
333 //=========================================================================================================
340 void newRawDataAvailable(const Eigen::VectorXd &dataLh,
341 const Eigen::VectorXd &dataRh);
342
343 //=========================================================================================================
349 void newInterpolationMatrixLeftAvailable(QSharedPointer<Eigen::SparseMatrix<float>> interpMat);
350
351 //=========================================================================================================
357 void newInterpolationMatrixRightAvailable(QSharedPointer<Eigen::SparseMatrix<float>> interpMat);
358
359private slots:
360 void onNewInterpolationMatrixLeft(QSharedPointer<Eigen::SparseMatrix<float>> interpMat);
361 void onNewInterpolationMatrixRight(QSharedPointer<Eigen::SparseMatrix<float>> interpMat);
362
363private:
364 QThread *m_pWorkerThread = nullptr;
365 DISP3DRHILIB::RtSourceDataWorker *m_pWorker = nullptr;
366 QTimer *m_pTimer = nullptr;
367 bool m_bIsStreaming = false;
368 int m_iTimeInterval = 17;
369
370 QThread *m_pInterpThread = nullptr;
371 DISP3DRHILIB::RtSourceInterpolationMatWorker *m_pInterpWorker = nullptr;
372};
373
374#endif // BRAINVIEW_RTSOURCEDATACONTROLLER_H
disp3D_rhi library export/import macros.
#define DISP3DRHISHARED_EXPORT
FreeSurfer surface and annotation I/O.
3-D brain visualisation using the Qt RHI rendering backend.
void setAnnotationInfoRight(const Eigen::VectorXi &vecLabelIds, const QList< FSLIB::Label > &lLabels, const Eigen::VectorXi &vecVertNo)
void setThresholds(double min, double mid, double max)
void setInterpolationMatrixLeft(QSharedPointer< Eigen::SparseMatrix< float > > mat)
void addData(const Eigen::VectorXd &data)
void setStreamSmoothedData(bool bStreamSmoothedData)
void setAnnotationInfoLeft(const Eigen::VectorXi &vecLabelIds, const QList< FSLIB::Label > &lLabels, const Eigen::VectorXi &vecVertNo)
void newInterpolationMatrixRightAvailable(QSharedPointer< Eigen::SparseMatrix< float > > interpMat)
void setSurfaceColor(const QVector< uint32_t > &baseColorsLh, const QVector< uint32_t > &baseColorsRh)
void setCancelDistance(double dCancelDist)
void setInterpolationInfoRight(const Eigen::MatrixX3f &matVertices, const std::vector< Eigen::VectorXi > &vecNeighborVertices, const Eigen::VectorXi &vecSourceVertices)
void newRawDataAvailable(const Eigen::VectorXd &dataLh, const Eigen::VectorXd &dataRh)
void newSmoothedDataAvailable(const QVector< uint32_t > &colorsLh, const QVector< uint32_t > &colorsRh)
void setVisualizationType(int iVisType)
void setInterpolationFunction(const QString &sInterpolationFunction)
void setColormapType(const QString &name)
void newInterpolationMatrixLeftAvailable(QSharedPointer< Eigen::SparseMatrix< float > > interpMat)
RtSourceDataController(QObject *parent=nullptr)
void setInterpolationMatrixRight(QSharedPointer< Eigen::SparseMatrix< float > > mat)
void setInterpolationInfoLeft(const Eigen::MatrixX3f &matVertices, const std::vector< Eigen::VectorXi > &vecNeighborVertices, const Eigen::VectorXi &vecSourceVertices)
Background worker for real-time source estimate streaming.
Background worker for computing source interpolation matrices.
Freesurfer/MNE label.
Definition label.h:81