MNE-CPP  0.1.9
A Framework for Electrophysiology
fiff_cov.cpp
Go to the documentation of this file.
1 //=============================================================================================================
37 //=============================================================================================================
38 // INCLUDES
39 //=============================================================================================================
40 
41 #include "fiff_cov.h"
42 #include "fiff_stream.h"
43 #include "fiff_info_base.h"
44 #include "fiff_dir_node.h"
45 
46 #include <utils/mnemath.h>
47 
48 //=============================================================================================================
49 // QT INCLUDES
50 //=============================================================================================================
51 
52 #include <QPair>
53 
54 //=============================================================================================================
55 // EIGEN INCLUDES
56 //=============================================================================================================
57 
58 #include <Eigen/SVD>
59 
60 //=============================================================================================================
61 // USED NAMESPACES
62 //=============================================================================================================
63 
64 using namespace FIFFLIB;
65 using namespace UTILSLIB;
66 using namespace Eigen;
67 
68 //=============================================================================================================
69 // DEFINE MEMBER METHODS
70 //=============================================================================================================
71 
73 : kind(-1)
74 , diag(false)
75 , dim(-1)
76 , nfree(-1)
77 {
78  qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
79  qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
80 }
81 
82 //=============================================================================================================
83 
84 FiffCov::FiffCov(QIODevice &p_IODevice)
85 : kind(-1)
86 , diag(false)
87 , dim(-1)
88 , nfree(-1)
89 {
90  FiffStream::SPtr t_pStream(new FiffStream(&p_IODevice));
91 
92  if(!t_pStream->open())
93  {
94  printf("\tNot able to open IODevice.\n");//ToDo Throw here
95  return;
96  }
97 
98  if(!t_pStream->read_cov(t_pStream->dirtree(), FIFFV_MNE_NOISE_COV, *this))
99  printf("\tFiff covariance not found.\n");//ToDo Throw here
100 
101  qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
102  qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
103 }
104 
105 //=============================================================================================================
106 
107 FiffCov::FiffCov(const FiffCov &p_FiffCov)
108 : QSharedData(p_FiffCov)
109 , kind(p_FiffCov.kind)
110 , diag(p_FiffCov.diag)
111 , dim(p_FiffCov.dim)
112 , names(p_FiffCov.names)
113 , data(p_FiffCov.data)
114 , projs(p_FiffCov.projs)
115 , bads(p_FiffCov.bads)
116 , nfree(p_FiffCov.nfree)
117 , eig(p_FiffCov.eig)
118 , eigvec(p_FiffCov.eigvec)
119 {
120  qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
121  qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
122 }
123 
124 //=============================================================================================================
125 
127 {
128 }
129 
130 //=============================================================================================================
131 
133 {
134  kind = -1;
135  diag = false;
136  dim = -1;
137  names.clear();
138  data = MatrixXd();
139  projs.clear();
140  bads.clear();
141  nfree = -1;
142  eig = VectorXd();
143  eigvec = MatrixXd();
144 }
145 
146 //=============================================================================================================
147 
148 FiffCov FiffCov::pick_channels(const QStringList &p_include, const QStringList &p_exclude)
149 {
150  RowVectorXi sel = FiffInfoBase::pick_channels(this->names, p_include, p_exclude);
151  FiffCov res;//No deep copy here - since almost everything else is adapted anyway
152 
153  res.kind = this->kind;
154  res.diag = this->diag;
155  res.dim = sel.size();
156 
157  for(qint32 k = 0; k < sel.size(); ++k)
158  res.names << this->names[sel(k)];
159 
160  res.data.resize(res.dim, res.dim);
161  for(qint32 i = 0; i < res.dim; ++i)
162  for(qint32 j = 0; j < res.dim; ++j)
163  res.data(i, j) = this->data(sel(i), sel(j));
164  res.projs = this->projs;
165 
166  for(qint32 k = 0; k < this->bads.size(); ++k)
167  if(res.names.contains(this->bads[k]))
168  res.bads << this->bads[k];
169  res.nfree = this->nfree;
170 
171  return res;
172 }
173 
174 //=============================================================================================================
175 
176 FiffCov FiffCov::prepare_noise_cov(const FiffInfo &p_Info, const QStringList &p_ChNames) const
177 {
178  FiffCov p_NoiseCov(*this);
179 
180  VectorXi C_ch_idx = VectorXi::Zero(p_NoiseCov.names.size());
181  qint32 count = 0;
182  for(qint32 i = 0; i < p_ChNames.size(); ++i)
183  {
184  qint32 idx = p_NoiseCov.names.indexOf(p_ChNames[i]);
185  if(idx > -1)
186  {
187  C_ch_idx[count] = idx;
188  ++count;
189  }
190  }
191  C_ch_idx.conservativeResize(count);
192 
193  MatrixXd C(count, count);
194 
195  if(!p_NoiseCov.diag)
196  for(qint32 i = 0; i < count; ++i)
197  for(qint32 j = 0; j < count; ++j)
198  C(i,j) = p_NoiseCov.data(C_ch_idx(i), C_ch_idx(j));
199  else
200  {
201  qWarning("Warning in FiffCov::prepare_noise_cov: This has to be debugged - not done before!");
202  C = MatrixXd::Zero(count, count);
203  for(qint32 i = 0; i < count; ++i)
204  C.diagonal()[i] = p_NoiseCov.data(C_ch_idx(i),0);
205  }
206 
207  MatrixXd proj;
208  qint32 ncomp = p_Info.make_projector(proj, p_ChNames);
209 
210  //Create the projection operator
211  if (ncomp > 0 && proj.rows() == count)
212  {
213  printf("Created an SSP operator (subspace dimension = %d)\n", ncomp);
214  C = proj * (C * proj.transpose());
215  } else {
216  qWarning("Warning in FiffCov::prepare_noise_cov: No projections applied since no projectors specified or projector dimensions do not match!");
217  }
218 
219  RowVectorXi pick_meg = p_Info.pick_types(true, false, false, defaultQStringList, p_Info.bads);
220  RowVectorXi pick_eeg = p_Info.pick_types(false, true, false, defaultQStringList, p_Info.bads);
221 
222  QStringList meg_names, eeg_names;
223 
224  for(qint32 i = 0; i < pick_meg.size(); ++i)
225  meg_names << p_Info.chs[pick_meg[i]].ch_name;
226  VectorXi C_meg_idx = VectorXi::Zero(p_NoiseCov.names.size());
227  count = 0;
228  for(qint32 k = 0; k < C.rows(); ++k)
229  {
230  if(meg_names.indexOf(p_ChNames[k]) > -1)
231  {
232  C_meg_idx[count] = k;
233  ++count;
234  }
235  }
236  if(count > 0)
237  C_meg_idx.conservativeResize(count);
238  else
239  C_meg_idx = VectorXi();
240 
241  //
242  for(qint32 i = 0; i < pick_eeg.size(); ++i)
243  eeg_names << p_Info.chs[pick_eeg(0,i)].ch_name;
244  VectorXi C_eeg_idx = VectorXi::Zero(p_NoiseCov.names.size());
245  count = 0;
246  for(qint32 k = 0; k < C.rows(); ++k)
247  {
248  if(eeg_names.indexOf(p_ChNames[k]) > -1)
249  {
250  C_eeg_idx[count] = k;
251  ++count;
252  }
253  }
254 
255  if(count > 0)
256  C_eeg_idx.conservativeResize(count);
257  else
258  C_eeg_idx = VectorXi();
259 
260  bool has_meg = C_meg_idx.size() > 0;
261  bool has_eeg = C_eeg_idx.size() > 0;
262 
263  MatrixXd C_meg, C_eeg;
264  VectorXd C_meg_eig, C_eeg_eig;
265  MatrixXd C_meg_eigvec, C_eeg_eigvec;
266  if (has_meg)
267  {
268  count = C_meg_idx.rows();
269  C_meg = MatrixXd(count,count);
270  for(qint32 i = 0; i < count; ++i)
271  for(qint32 j = 0; j < count; ++j)
272  C_meg(i,j) = C(C_meg_idx(i), C_meg_idx(j));
273  MNEMath::get_whitener(C_meg, false, QString("MEG"), C_meg_eig, C_meg_eigvec);
274  }
275 
276  if (has_eeg)
277  {
278  count = C_eeg_idx.rows();
279  C_eeg = MatrixXd(count,count);
280  for(qint32 i = 0; i < count; ++i)
281  for(qint32 j = 0; j < count; ++j)
282  C_eeg(i,j) = C(C_eeg_idx(i), C_eeg_idx(j));
283  MNEMath::get_whitener(C_eeg, false, QString("EEG"), C_eeg_eig, C_eeg_eigvec);
284  }
285 
286  qint32 n_chan = p_ChNames.size();
287  p_NoiseCov.eigvec = MatrixXd::Zero(n_chan, n_chan);
288  p_NoiseCov.eig = VectorXd::Zero(n_chan);
289 
290  if(has_meg)
291  {
292  for(qint32 i = 0; i < C_meg_idx.rows(); ++i)
293  for(qint32 j = 0; j < C_meg_idx.rows(); ++j)
294  p_NoiseCov.eigvec(C_meg_idx[i], C_meg_idx[j]) = C_meg_eigvec(i, j);
295  for(qint32 i = 0; i < C_meg_idx.rows(); ++i)
296  p_NoiseCov.eig(C_meg_idx[i]) = C_meg_eig[i];
297  }
298  if(has_eeg)
299  {
300  for(qint32 i = 0; i < C_eeg_idx.rows(); ++i)
301  for(qint32 j = 0; j < C_eeg_idx.rows(); ++j)
302  p_NoiseCov.eigvec(C_eeg_idx[i], C_eeg_idx[j]) = C_eeg_eigvec(i, j);
303  for(qint32 i = 0; i < C_eeg_idx.rows(); ++i)
304  p_NoiseCov.eig(C_eeg_idx[i]) = C_eeg_eig[i];
305  }
306 
307  if (C_meg_idx.size() + C_eeg_idx.size() != n_chan)
308  {
309  printf("Error in FiffCov::prepare_noise_cov: channel sizes do no match!\n");//ToDo Throw here
310  return FiffCov();
311  }
312 
313  p_NoiseCov.data = C;
314  p_NoiseCov.dim = p_ChNames.size();
315  p_NoiseCov.diag = false;
316  p_NoiseCov.names = p_ChNames;
317 
318  return p_NoiseCov;
319 }
320 
321 //=============================================================================================================
322 
323 FiffCov FiffCov::regularize(const FiffInfo& p_info, double p_fRegMag, double p_fRegGrad, double p_fRegEeg, bool p_bProj, QStringList p_exclude) const
324 {
325  FiffCov cov(*this);
326 
327  if(p_exclude.size() == 0)
328  {
329  p_exclude = p_info.bads;
330  for(qint32 i = 0; i < cov.bads.size(); ++i)
331  if(!p_exclude.contains(cov.bads[i]))
332  p_exclude << cov.bads[i];
333  }
334 
335  //Allways exclude all STI channels from covariance computation
336  int iNoStimCh = 0;
337 
338  for(int i=0; i<p_info.chs.size(); i++) {
339  if(p_info.chs[i].kind == FIFFV_STIM_CH) {
340  p_exclude << p_info.chs[i].ch_name;
341  iNoStimCh++;
342  }
343  }
344 
345  RowVectorXi sel_eeg = p_info.pick_types(false, true, false, defaultQStringList, p_exclude);
346  RowVectorXi sel_mag = p_info.pick_types(QString("mag"), false, false, defaultQStringList, p_exclude);
347  RowVectorXi sel_grad = p_info.pick_types(QString("grad"), false, false, defaultQStringList, p_exclude);
348 
349  QStringList info_ch_names = p_info.ch_names;
350  QStringList ch_names_eeg, ch_names_mag, ch_names_grad;
351  for(qint32 i = 0; i < sel_eeg.size(); ++i)
352  ch_names_eeg << info_ch_names[sel_eeg(i)];
353  for(qint32 i = 0; i < sel_mag.size(); ++i)
354  ch_names_mag << info_ch_names[sel_mag(i)];
355  for(qint32 i = 0; i < sel_grad.size(); ++i)
356  ch_names_grad << info_ch_names[sel_grad(i)];
357 
358  // This actually removes bad channels from the cov, which is not backward
359  // compatible, so let's leave all channels in
360  FiffCov cov_good = cov.pick_channels(info_ch_names, p_exclude);
361  QStringList ch_names = cov_good.names;
362 
363  std::vector<qint32> idx_eeg, idx_mag, idx_grad;
364  for(qint32 i = 0; i < ch_names.size(); ++i)
365  {
366  if(ch_names_eeg.contains(ch_names[i]))
367  idx_eeg.push_back(i);
368  else if(ch_names_mag.contains(ch_names[i]))
369  idx_mag.push_back(i);
370  else if(ch_names_grad.contains(ch_names[i]))
371  idx_grad.push_back(i);
372  }
373 
374  MatrixXd C(cov_good.data);
375 
376  //Subtract number of found stim channels because they are still in C but not the idx_eeg, idx_mag or idx_grad
377  if((unsigned) C.rows() - iNoStimCh != idx_eeg.size() + idx_mag.size() + idx_grad.size()) {
378  printf("Error in FiffCov::regularize: Channel dimensions do not fit.\n");//ToDo Throw
379  }
380 
381  QList<FiffProj> t_listProjs;
382  if(p_bProj)
383  {
384  t_listProjs = p_info.projs + cov_good.projs;
385  FiffProj::activate_projs(t_listProjs);
386  }
387 
388  //Build regularization MAP
389  QMap<QString, QPair<double, std::vector<qint32> > > regData;
390  regData.insert("EEG", QPair<double, std::vector<qint32> >(p_fRegEeg, idx_eeg));
391  regData.insert("MAG", QPair<double, std::vector<qint32> >(p_fRegMag, idx_mag));
392  regData.insert("GRAD", QPair<double, std::vector<qint32> >(p_fRegGrad, idx_grad));
393 
394  //
395  //Regularize
396  //
397  QMap<QString, QPair<double, std::vector<qint32> > >::Iterator it;
398  for(it = regData.begin(); it != regData.end(); ++it)
399  {
400  QString desc(it.key());
401  double reg = it.value().first;
402  std::vector<qint32> idx = it.value().second;
403 
404  if(idx.size() == 0 || reg == 0.0)
405  printf("\tNothing to regularize within %s data.\n", desc.toUtf8().constData());
406  else
407  {
408  printf("\tRegularize %s: %f\n", desc.toUtf8().constData(), reg);
409  MatrixXd this_C(idx.size(), idx.size());
410  for(quint32 i = 0; i < idx.size(); ++i)
411  for(quint32 j = 0; j < idx.size(); ++j)
412  this_C(i,j) = cov_good.data(idx[i], idx[j]);
413 
414  MatrixXd U;
415  qint32 ncomp;
416  if(p_bProj)
417  {
418  QStringList this_ch_names;
419  for(quint32 k = 0; k < idx.size(); ++k)
420  this_ch_names << ch_names[idx[k]];
421 
422  MatrixXd P;
423  ncomp = FiffProj::make_projector(t_listProjs, this_ch_names, P); //ToDo: Synchronize with mne-python and debug
424 
425  JacobiSVD<MatrixXd> svd(P, ComputeFullU);
426  //Sort singular values and singular vectors
427  VectorXd t_s = svd.singularValues();
428  MatrixXd t_U = svd.matrixU();
429  MNEMath::sort<double>(t_s, t_U);
430 
431  U = t_U.block(0,0, t_U.rows(), t_U.cols()-ncomp);
432 
433  if (ncomp > 0)
434  {
435  printf("\tCreated an SSP operator for %s (dimension = %d).\n", desc.toUtf8().constData(), ncomp);
436  this_C = U.transpose() * (this_C * U);
437  }
438  }
439 
440  double sigma = this_C.diagonal().mean();
441  this_C.diagonal() = this_C.diagonal().array() + reg * sigma; // modify diag inplace
442  if(p_bProj && ncomp > 0)
443  this_C = U * (this_C * U.transpose());
444 
445  for(qint32 i = 0; i < this_C.rows(); ++i)
446  for(qint32 j = 0; j < this_C.cols(); ++j)
447  C(idx[i],idx[j]) = this_C(i,j);
448  }
449  }
450 
451  // Put data back in correct locations
452  RowVectorXi idx = FiffInfo::pick_channels(cov.names, info_ch_names, p_exclude);
453  for(qint32 i = 0; i < idx.size(); ++i)
454  for(qint32 j = 0; j < idx.size(); ++j)
455  cov.data(idx[i], idx[j]) = C(i, j);
456 
457  return cov;
458 }
459 
460 //=============================================================================================================
461 
463 {
464  if (this != &rhs) // protect against invalid self-assignment
465  {
466  kind = rhs.kind;
467  diag = rhs.diag;
468  dim = rhs.dim;
469  names = rhs.names;
470  data = rhs.data;
471  projs = rhs.projs;
472  bads = rhs.bads;
473  nfree = rhs.nfree;
474  eig = rhs.eig;
475  eigvec = rhs.eigvec;
476  }
477  // to support chained assignment operators (a=b=c), always return *this
478  return *this;
479 }
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
fiff_info_base.h
FiffInfoBase class declaration.
FIFFLIB::FiffProj::activate_projs
static void activate_projs(QList< FiffProj > &p_qListFiffProj)
Definition: fiff_proj.cpp:103
FIFFLIB::FiffStream::SPtr
QSharedPointer< FiffStream > SPtr
Definition: fiff_stream.h:107
FIFFLIB::FiffInfo
FIFF measurement file information.
Definition: fiff_info.h:84
FIFFLIB::FiffCov::names
QStringList names
Definition: fiff_cov.h:194
FIFFLIB::FiffCov::~FiffCov
~FiffCov()
Definition: fiff_cov.cpp:126
fiff_cov.h
FiffCov class declaration.
FIFFLIB::FiffInfo::make_projector
qint32 make_projector(Eigen::MatrixXd &proj) const
Definition: fiff_info.h:278
FIFFLIB::FiffInfoBase::pick_types
Eigen::RowVectorXi pick_types(const QString meg, bool eeg=false, bool stim=false, const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList) const
Definition: fiff_info_base.cpp:132
fiff_stream.h
FiffStream class declaration.
FIFFLIB::FiffCov::FiffCov
FiffCov()
Definition: fiff_cov.cpp:72
FIFFLIB::FiffCov::kind
fiff_int_t kind
Definition: fiff_cov.h:190
FIFFLIB::FiffCov::prepare_noise_cov
FiffCov prepare_noise_cov(const FiffInfo &p_info, const QStringList &p_chNames) const
Definition: fiff_cov.cpp:176
FIFFLIB::FiffInfoBase::pick_channels
static Eigen::RowVectorXi pick_channels(const QStringList &ch_names, const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList)
Definition: fiff_info_base.cpp:198
k
int k
Definition: fiff_tag.cpp:322
FIFFLIB::FiffInfoBase::bads
QStringList bads
Definition: fiff_info_base.h:220
FIFFLIB::FiffCov::operator=
FiffCov & operator=(const FiffCov &rhs)
Definition: fiff_cov.cpp:462
FIFFLIB::FiffInfo::projs
QList< FiffProj > projs
Definition: fiff_info.h:268
FIFFLIB::FiffCov::diag
bool diag
Definition: fiff_cov.h:192
FIFFLIB::FiffCov::eig
Eigen::VectorXd eig
Definition: fiff_cov.h:199
FIFFLIB::FiffCov::clear
void clear()
Definition: fiff_cov.cpp:132
fiff_dir_node.h
FiffDirNode class declaration, which provides fiff dir tree processing methods.
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
FIFFLIB::FiffInfoBase::ch_names
QStringList ch_names
Definition: fiff_info_base.h:224
FIFFLIB::FiffStream
FIFF File I/O routines.
Definition: fiff_stream.h:104
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::pick_channels
FiffCov pick_channels(const QStringList &p_include=defaultQStringList, const QStringList &p_exclude=defaultQStringList)
Definition: fiff_cov.cpp:148
FIFFLIB::FiffCov::dim
fiff_int_t dim
Definition: fiff_cov.h:193
FIFFLIB::FiffProj::make_projector
static fiff_int_t make_projector(const QList< FiffProj > &projs, const QStringList &ch_names, Eigen::MatrixXd &proj, const QStringList &bads=defaultQStringList, Eigen::MatrixXd &U=defaultMatrixXd)
Definition: fiff_proj.cpp:115
mnemath.h
MNEMath class declaration.
FIFFLIB::FiffCov::eigvec
Eigen::MatrixXd eigvec
Definition: fiff_cov.h:200