MNE-CPP  0.1.9
A Framework for Electrophysiology
rtcov.cpp
Go to the documentation of this file.
1 //=============================================================================================================
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 // USED NAMESPACES
52 //=============================================================================================================
53 
54 using namespace RTPROCESSINGLIB;
55 using namespace FIFFLIB;
56 using namespace Eigen;
57 
58 //=============================================================================================================
59 // DEFINE MEMBER METHODS
60 //=============================================================================================================
61 
62 RtCov::RtCov(QSharedPointer<FIFFLIB::FiffInfo> pFiffInfo)
63 : m_fiffInfo(*pFiffInfo)
64 , m_iSamples(0)
65 {
66 }
67 
68 //=============================================================================================================
69 
70 FiffCov RtCov::estimateCovariance(const Eigen::MatrixXd& matData,
71  int iNewMaxSamples)
72 {
73  if(m_fiffInfo.chs.isEmpty()) {
74  qWarning() << "[RtCov::estimateCovariance] FiffInfo was not set. Returning empty covariance estimation.";
75  return FiffCov();
76  }
77 
78  m_lData.append(matData);
79  m_iSamples += matData.cols();
80 
81  if(m_iSamples < iNewMaxSamples) {
82  return FiffCov();
83  }
84 
85  QFuture<RtCovComputeResult> result = QtConcurrent::mappedReduced(m_lData,
86  compute,
87  reduce);
88 
89  result.waitForFinished();
90 
91  RtCovComputeResult finalResult = result.result();
92 
93  //Final computation
94  FiffCov computedCov;
95  computedCov.data = finalResult.matData;
96 
97  QStringList exclude;
98  for(int i = 0; i<m_fiffInfo.chs.size(); i++) {
99  if(m_fiffInfo.chs.at(i).kind != FIFFV_MEG_CH &&
100  m_fiffInfo.chs.at(i).kind != FIFFV_EEG_CH) {
101  exclude << m_fiffInfo.chs.at(i).ch_name;
102  }
103  }
104  bool doProj = true;
105 
106  if(m_iSamples > 0) {
107  finalResult.mu /= (float)m_iSamples;
108  computedCov.data.array() -= m_iSamples * (finalResult.mu * finalResult.mu.transpose()).array();
109  computedCov.data.array() /= (m_iSamples - 1);
110 
111  computedCov.kind = FIFFV_MNE_NOISE_COV;
112  computedCov.diag = false;
113  computedCov.dim = computedCov.data.rows();
114 
115  //ToDo do picks
116  computedCov.names = m_fiffInfo.ch_names;
117  computedCov.projs = m_fiffInfo.projs;
118  computedCov.bads = m_fiffInfo.bads;
119  computedCov.nfree = m_iSamples;
120 
121  // regularize noise covariance
122  computedCov = computedCov.regularize(m_fiffInfo, 0.05, 0.05, 0.1, doProj, exclude);
123 
124 // qint32 samples = rawSegment.cols();
125 // VectorXf mu = rawSegment.rowwise().sum().array() / (float)samples;
126 
127 // MatrixXf noise_covariance = rawSegment * rawSegment.transpose();// noise_covariance == raw_covariance
128 // noise_covariance.array() -= samples * (mu * mu.transpose()).array();
129 // noise_covariance.array() /= (samples - 1);
130 
131 // std::cout << "Noise Covariance:\n" << noise_covariance.block(0,0,10,10) << std::endl;
132 
133 // printf("%d raw buffer (%d x %d) generated\r\n", count, tmp.rows(), tmp.cols());
134 
135  m_lData.clear();
136  m_iSamples = 0;
137 
138  return computedCov;
139  } else {
140  qWarning() << "[RtCov::estimateCovariance] Number of samples equals zero. Regularization not possible. Returning empty covariance estimation.";
141  return FiffCov();
142  }
143 
144 }
145 
146 //=============================================================================================================
147 
148 RtCovComputeResult RtCov::compute(const MatrixXd &matData)
149 {
150  RtCovComputeResult result;
151  result.mu = matData.rowwise().sum();
152  result.matData = matData * matData.transpose();
153  return result;
154 }
155 
156 //=============================================================================================================
157 
158 void RtCov::reduce(RtCovComputeResult& finalResult, const RtCovComputeResult &tempResult)
159 {
160  if(finalResult.matData.size() == 0 || finalResult.mu.size() == 0) {
161  finalResult.mu = tempResult.mu;
162  finalResult.matData = tempResult.matData;
163  } else {
164  finalResult.mu += tempResult.mu;
165  finalResult.matData += tempResult.matData;
166  }
167 }
FIFFLIB::FiffCov::data
Eigen::MatrixXd data
Definition: fiff_cov.h:195
FIFFV_MNE_NOISE_COV
#define FIFFV_MNE_NOISE_COV
Definition: fiff_constants.h:448
FIFFLIB::FiffCov::bads
QStringList bads
Definition: fiff_cov.h:197
rtcov.h
RtCov class declaration.
FIFFLIB::FiffCov::names
QStringList names
Definition: fiff_cov.h:194
RTPROCESSINGLIB::RtCov::reduce
static void reduce(RtCovComputeResult &finalResult, const RtCovComputeResult &tempResult)
Definition: rtcov.cpp:158
FIFFLIB::FiffCov::kind
fiff_int_t kind
Definition: fiff_cov.h:190
FIFFLIB::FiffInfoBase::bads
QStringList bads
Definition: fiff_info_base.h:220
RTPROCESSINGLIB::RtCov::compute
static RtCovComputeResult compute(const Eigen::MatrixXd &matData)
Definition: rtcov.cpp:148
FIFFLIB::FiffInfo::projs
QList< FiffProj > projs
Definition: fiff_info.h:268
RTPROCESSINGLIB::RtCov::m_lData
QList< Eigen::MatrixXd > m_lData
Definition: rtcov.h:129
FIFFLIB::FiffCov::diag
bool diag
Definition: fiff_cov.h:192
RTPROCESSINGLIB::RtCovComputeResult
Definition: rtcov.h:80
FIFFLIB::FiffCov
covariance data
Definition: fiff_cov.h:77
FIFFLIB::FiffInfoBase::chs
QList< FiffChInfo > chs
Definition: fiff_info_base.h:223
FIFFLIB::FiffCov::nfree
fiff_int_t nfree
Definition: fiff_cov.h:198
RTPROCESSINGLIB::RtCov::m_iSamples
int m_iSamples
Definition: rtcov.h:127
FIFFLIB::FiffInfoBase::ch_names
QStringList ch_names
Definition: fiff_info_base.h:224
FIFFLIB::FiffCov::regularize
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:323
FIFFLIB::FiffCov::projs
QList< FiffProj > projs
Definition: fiff_cov.h:196
FIFFLIB::FiffCov::dim
fiff_int_t dim
Definition: fiff_cov.h:193
RTPROCESSINGLIB::RtCov::m_fiffInfo
FIFFLIB::FiffInfo m_fiffInfo
Definition: rtcov.h:131
RTPROCESSINGLIB::RtCov::estimateCovariance
FIFFLIB::FiffCov estimateCovariance(const Eigen::MatrixXd &matData, int iNewMaxSamples)
Definition: rtcov.cpp:70