v2.0.0
Loading...
Searching...
No Matches
rtcov.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "rtcov.h"
42
43//=============================================================================================================
44// QT INCLUDES
45//=============================================================================================================
46
47#include <QDebug>
48#include <QtConcurrent>
49
50//=============================================================================================================
51// STL INCLUDES
52//=============================================================================================================
53
54#include <iostream>
55
56//=============================================================================================================
57// USED NAMESPACES
58//=============================================================================================================
59
60using namespace RTPROCESSINGLIB;
61using namespace FIFFLIB;
62using namespace Eigen;
63
64//=============================================================================================================
65// DEFINE MEMBER METHODS
66//=============================================================================================================
67
68RtCov::RtCov(QSharedPointer<FIFFLIB::FiffInfo> pFiffInfo)
69: m_fiffInfo(*pFiffInfo)
70, m_iSamples(0)
71{
72}
73
74//=============================================================================================================
75
76FiffCov RtCov::estimateCovariance(const Eigen::MatrixXd& matData,
77 int iNewMaxSamples)
78{
79 if(m_fiffInfo.chs.isEmpty()) {
80 qWarning() << "[RtCov::estimateCovariance] FiffInfo was not set. Returning empty covariance estimation.";
81 return FiffCov();
82 }
83
84 m_lData.append(matData);
85 m_iSamples += matData.cols();
86
87 if(m_iSamples < iNewMaxSamples) {
88 return FiffCov();
89 }
90
91 QFuture<RtCovComputeResult> result = QtConcurrent::mappedReduced(m_lData,
92 compute,
93 reduce);
94
95 result.waitForFinished();
96
97 RtCovComputeResult finalResult = result.result();
98
99 //Final computation
100 FiffCov computedCov;
101 computedCov.data = finalResult.matData;
102
103 QStringList exclude;
104 for(int i = 0; i<m_fiffInfo.chs.size(); i++) {
105 if(m_fiffInfo.chs.at(i).kind != FIFFV_MEG_CH &&
106 m_fiffInfo.chs.at(i).kind != FIFFV_EEG_CH) {
107 exclude << m_fiffInfo.chs.at(i).ch_name;
108 }
109 }
110 bool doProj = true;
111
112 if(m_iSamples > 0) {
113 finalResult.mu /= (float)m_iSamples;
114 computedCov.data.array() -= m_iSamples * (finalResult.mu * finalResult.mu.transpose()).array();
115 computedCov.data.array() /= (m_iSamples - 1);
116
117 computedCov.kind = FIFFV_MNE_NOISE_COV;
118 computedCov.diag = false;
119 computedCov.dim = computedCov.data.rows();
120
121 //ToDo do picks
122 computedCov.names = m_fiffInfo.ch_names;
123 computedCov.projs = m_fiffInfo.projs;
124 computedCov.bads = m_fiffInfo.bads;
125 computedCov.nfree = m_iSamples;
126
127 // regularize noise covariance
128 computedCov = computedCov.regularize(m_fiffInfo, 0.05, 0.05, 0.1, doProj, exclude);
129
130// qint32 samples = rawSegment.cols();
131// VectorXf mu = rawSegment.rowwise().sum().array() / (float)samples;
132
133// MatrixXf noise_covariance = rawSegment * rawSegment.transpose();// noise_covariance == raw_covariance
134// noise_covariance.array() -= samples * (mu * mu.transpose()).array();
135// noise_covariance.array() /= (samples - 1);
136
137// std::cout << "Noise Covariance:\n" << noise_covariance.block(0,0,10,10) << std::endl;
138
139// printf("%d raw buffer (%d x %d) generated\r\n", count, tmp.rows(), tmp.cols());
140
141 m_lData.clear();
142 m_iSamples = 0;
143
144 return computedCov;
145 } else {
146 qWarning() << "[RtCov::estimateCovariance] Number of samples equals zero. Regularization not possible. Returning empty covariance estimation.";
147 return FiffCov();
148 }
149
150}
151
152//=============================================================================================================
153
154RtCovComputeResult RtCov::compute(const MatrixXd &matData)
155{
156 RtCovComputeResult result;
157 result.mu = matData.rowwise().sum();
158 result.matData = matData * matData.transpose();
159 return result;
160}
161
162//=============================================================================================================
163
164void RtCov::reduce(RtCovComputeResult& finalResult, const RtCovComputeResult &tempResult)
165{
166 if(finalResult.matData.size() == 0 || finalResult.mu.size() == 0) {
167 finalResult.mu = tempResult.mu;
168 finalResult.matData = tempResult.matData;
169 } else {
170 finalResult.mu += tempResult.mu;
171 finalResult.matData += tempResult.matData;
172 }
173}
RtCov class declaration.
#define FIFFV_EEG_CH
#define FIFFV_MNE_NOISE_COV
#define FIFFV_MEG_CH
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Real-time signal processing (filtering, averaging, HPI fitting, noise reduction).
covariance data
Definition fiff_cov.h:84
QList< FiffProj > projs
Definition fiff_cov.h:252
fiff_int_t nfree
Definition fiff_cov.h:254
fiff_int_t dim
Definition fiff_cov.h:249
FiffCov regularize(const FiffInfo &p_info, double p_fMag=0.1, double p_fGrad=0.1, double p_fEeg=0.1, bool p_bProj=true, QStringList p_exclude=defaultQStringList) const
Definition fiff_cov.cpp:326
fiff_int_t kind
Definition fiff_cov.h:246
QStringList bads
Definition fiff_cov.h:253
QStringList names
Definition fiff_cov.h:250
Eigen::MatrixXd data
Definition fiff_cov.h:251
Bundled output of a real-time covariance computation step containing the covariance matrix and sample...
Definition rtcov.h:83
QList< Eigen::MatrixXd > m_lData
Definition rtcov.h:132
RtCov(QSharedPointer< FIFFLIB::FiffInfo > pFiffInfo)
Definition rtcov.cpp:68
FIFFLIB::FiffCov estimateCovariance(const Eigen::MatrixXd &matData, int iNewMaxSamples)
Definition rtcov.cpp:76
FIFFLIB::FiffInfo m_fiffInfo
Definition rtcov.h:134
static void reduce(RtCovComputeResult &finalResult, const RtCovComputeResult &tempResult)
Definition rtcov.cpp:164
static RtCovComputeResult compute(const Eigen::MatrixXd &matData)
Definition rtcov.cpp:154