v2.0.0
Loading...
Searching...
No Matches
fiff_cov.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "fiff_cov.h"
42#include "fiff_stream.h"
43#include "fiff_raw_data.h"
44#include "fiff_info_base.h"
45#include "fiff_dir_node.h"
46#include "fiff_file.h"
47
48#include <utils/mnemath.h>
49
50//=============================================================================================================
51// QT INCLUDES
52//=============================================================================================================
53
54#include <QPair>
55#include <QFile>
56
57//=============================================================================================================
58// EIGEN INCLUDES
59//=============================================================================================================
60
61#include <Eigen/SVD>
62
63//=============================================================================================================
64// USED NAMESPACES
65//=============================================================================================================
66
67using namespace FIFFLIB;
68using namespace UTILSLIB;
69using namespace Eigen;
70
71//=============================================================================================================
72// DEFINE MEMBER METHODS
73//=============================================================================================================
74
76: kind(-1)
77, diag(false)
78, dim(-1)
79, nfree(-1)
80{
81 qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
82 qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
83}
84
85//=============================================================================================================
86
87FiffCov::FiffCov(QIODevice &p_IODevice)
88: kind(-1)
89, diag(false)
90, dim(-1)
91, nfree(-1)
92{
93 FiffStream::SPtr t_pStream(new FiffStream(&p_IODevice));
94
95 if(!t_pStream->open())
96 {
97 printf("\tNot able to open IODevice.\n");//ToDo Throw here
98 return;
99 }
100
101 if(!t_pStream->read_cov(t_pStream->dirtree(), FIFFV_MNE_NOISE_COV, *this))
102 printf("\tFiff covariance not found.\n");//ToDo Throw here
103
104 qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
105 qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
106}
107
108//=============================================================================================================
109
110FiffCov::FiffCov(const FiffCov &p_FiffCov)
111: QSharedData(p_FiffCov)
112, kind(p_FiffCov.kind)
113, diag(p_FiffCov.diag)
114, dim(p_FiffCov.dim)
115, names(p_FiffCov.names)
116, data(p_FiffCov.data)
117, projs(p_FiffCov.projs)
118, bads(p_FiffCov.bads)
119, nfree(p_FiffCov.nfree)
120, eig(p_FiffCov.eig)
121, eigvec(p_FiffCov.eigvec)
122{
123 qRegisterMetaType<QSharedPointer<FIFFLIB::FiffCov> >("QSharedPointer<FIFFLIB::FiffCov>");
124 qRegisterMetaType<FIFFLIB::FiffCov>("FIFFLIB::FiffCov");
125}
126
127//=============================================================================================================
128
132
133//=============================================================================================================
134
136{
137 kind = -1;
138 diag = false;
139 dim = -1;
140 names.clear();
141 data = MatrixXd();
142 projs.clear();
143 bads.clear();
144 nfree = -1;
145 eig = VectorXd();
146 eigvec = MatrixXd();
147}
148
149//=============================================================================================================
150
151FiffCov FiffCov::pick_channels(const QStringList &p_include, const QStringList &p_exclude)
152{
153 RowVectorXi sel = FiffInfoBase::pick_channels(this->names, p_include, p_exclude);
154 FiffCov res;//No deep copy here - since almost everything else is adapted anyway
155
156 res.kind = this->kind;
157 res.diag = this->diag;
158 res.dim = sel.size();
159
160 for(qint32 k = 0; k < sel.size(); ++k)
161 res.names << this->names[sel(k)];
162
163 res.data.resize(res.dim, res.dim);
164 for(qint32 i = 0; i < res.dim; ++i)
165 for(qint32 j = 0; j < res.dim; ++j)
166 res.data(i, j) = this->data(sel(i), sel(j));
167 res.projs = this->projs;
168
169 for(qint32 k = 0; k < this->bads.size(); ++k)
170 if(res.names.contains(this->bads[k]))
171 res.bads << this->bads[k];
172 res.nfree = this->nfree;
173
174 return res;
175}
176
177//=============================================================================================================
178
179FiffCov FiffCov::prepare_noise_cov(const FiffInfo &p_Info, const QStringList &p_ChNames) const
180{
181 FiffCov p_NoiseCov(*this);
182
183 VectorXi C_ch_idx = VectorXi::Zero(p_NoiseCov.names.size());
184 qint32 count = 0;
185 for(qint32 i = 0; i < p_ChNames.size(); ++i)
186 {
187 qint32 idx = p_NoiseCov.names.indexOf(p_ChNames[i]);
188 if(idx > -1)
189 {
190 C_ch_idx[count] = idx;
191 ++count;
192 }
193 }
194 C_ch_idx.conservativeResize(count);
195
196 MatrixXd C(count, count);
197
198 if(!p_NoiseCov.diag)
199 for(qint32 i = 0; i < count; ++i)
200 for(qint32 j = 0; j < count; ++j)
201 C(i,j) = p_NoiseCov.data(C_ch_idx(i), C_ch_idx(j));
202 else
203 {
204 qWarning("Warning in FiffCov::prepare_noise_cov: This has to be debugged - not done before!");
205 C = MatrixXd::Zero(count, count);
206 for(qint32 i = 0; i < count; ++i)
207 C.diagonal()[i] = p_NoiseCov.data(C_ch_idx(i),0);
208 }
209
210 MatrixXd proj;
211 qint32 ncomp = p_Info.make_projector(proj, p_ChNames);
212
213 //Create the projection operator
214 if (ncomp > 0 && proj.rows() == count)
215 {
216 printf("Created an SSP operator (subspace dimension = %d)\n", ncomp);
217 C = proj * (C * proj.transpose());
218 } else {
219 qWarning("Warning in FiffCov::prepare_noise_cov: No projections applied since no projectors specified or projector dimensions do not match!");
220 }
221
222 RowVectorXi pick_meg = p_Info.pick_types(true, false, false, defaultQStringList, p_Info.bads);
223 RowVectorXi pick_eeg = p_Info.pick_types(false, true, false, defaultQStringList, p_Info.bads);
224
225 QStringList meg_names, eeg_names;
226
227 for(qint32 i = 0; i < pick_meg.size(); ++i)
228 meg_names << p_Info.chs[pick_meg[i]].ch_name;
229 VectorXi C_meg_idx = VectorXi::Zero(p_NoiseCov.names.size());
230 count = 0;
231 for(qint32 k = 0; k < C.rows(); ++k)
232 {
233 if(meg_names.indexOf(p_ChNames[k]) > -1)
234 {
235 C_meg_idx[count] = k;
236 ++count;
237 }
238 }
239 if(count > 0)
240 C_meg_idx.conservativeResize(count);
241 else
242 C_meg_idx = VectorXi();
243
244 //
245 for(qint32 i = 0; i < pick_eeg.size(); ++i)
246 eeg_names << p_Info.chs[pick_eeg(0,i)].ch_name;
247 VectorXi C_eeg_idx = VectorXi::Zero(p_NoiseCov.names.size());
248 count = 0;
249 for(qint32 k = 0; k < C.rows(); ++k)
250 {
251 if(eeg_names.indexOf(p_ChNames[k]) > -1)
252 {
253 C_eeg_idx[count] = k;
254 ++count;
255 }
256 }
257
258 if(count > 0)
259 C_eeg_idx.conservativeResize(count);
260 else
261 C_eeg_idx = VectorXi();
262
263 bool has_meg = C_meg_idx.size() > 0;
264 bool has_eeg = C_eeg_idx.size() > 0;
265
266 MatrixXd C_meg, C_eeg;
267 VectorXd C_meg_eig, C_eeg_eig;
268 MatrixXd C_meg_eigvec, C_eeg_eigvec;
269 if (has_meg)
270 {
271 count = C_meg_idx.rows();
272 C_meg = MatrixXd(count,count);
273 for(qint32 i = 0; i < count; ++i)
274 for(qint32 j = 0; j < count; ++j)
275 C_meg(i,j) = C(C_meg_idx(i), C_meg_idx(j));
276 MNEMath::get_whitener(C_meg, false, QString("MEG"), C_meg_eig, C_meg_eigvec);
277 }
278
279 if (has_eeg)
280 {
281 count = C_eeg_idx.rows();
282 C_eeg = MatrixXd(count,count);
283 for(qint32 i = 0; i < count; ++i)
284 for(qint32 j = 0; j < count; ++j)
285 C_eeg(i,j) = C(C_eeg_idx(i), C_eeg_idx(j));
286 MNEMath::get_whitener(C_eeg, false, QString("EEG"), C_eeg_eig, C_eeg_eigvec);
287 }
288
289 qint32 n_chan = p_ChNames.size();
290 p_NoiseCov.eigvec = MatrixXd::Zero(n_chan, n_chan);
291 p_NoiseCov.eig = VectorXd::Zero(n_chan);
292
293 if(has_meg)
294 {
295 for(qint32 i = 0; i < C_meg_idx.rows(); ++i)
296 for(qint32 j = 0; j < C_meg_idx.rows(); ++j)
297 p_NoiseCov.eigvec(C_meg_idx[i], C_meg_idx[j]) = C_meg_eigvec(i, j);
298 for(qint32 i = 0; i < C_meg_idx.rows(); ++i)
299 p_NoiseCov.eig(C_meg_idx[i]) = C_meg_eig[i];
300 }
301 if(has_eeg)
302 {
303 for(qint32 i = 0; i < C_eeg_idx.rows(); ++i)
304 for(qint32 j = 0; j < C_eeg_idx.rows(); ++j)
305 p_NoiseCov.eigvec(C_eeg_idx[i], C_eeg_idx[j]) = C_eeg_eigvec(i, j);
306 for(qint32 i = 0; i < C_eeg_idx.rows(); ++i)
307 p_NoiseCov.eig(C_eeg_idx[i]) = C_eeg_eig[i];
308 }
309
310 if (C_meg_idx.size() + C_eeg_idx.size() != n_chan)
311 {
312 printf("Error in FiffCov::prepare_noise_cov: channel sizes do no match!\n");//ToDo Throw here
313 return FiffCov();
314 }
315
316 p_NoiseCov.data = C;
317 p_NoiseCov.dim = p_ChNames.size();
318 p_NoiseCov.diag = false;
319 p_NoiseCov.names = p_ChNames;
320
321 return p_NoiseCov;
322}
323
324//=============================================================================================================
325
326FiffCov FiffCov::regularize(const FiffInfo& p_info, double p_fRegMag, double p_fRegGrad, double p_fRegEeg, bool p_bProj, QStringList p_exclude) const
327{
328 FiffCov cov(*this);
329
330 if(p_exclude.size() == 0)
331 {
332 p_exclude = p_info.bads;
333 for(qint32 i = 0; i < cov.bads.size(); ++i)
334 if(!p_exclude.contains(cov.bads[i]))
335 p_exclude << cov.bads[i];
336 }
337
338 //Allways exclude all STI channels from covariance computation
339 int iNoStimCh = 0;
340
341 for(int i=0; i<p_info.chs.size(); i++) {
342 if(p_info.chs[i].kind == FIFFV_STIM_CH) {
343 p_exclude << p_info.chs[i].ch_name;
344 iNoStimCh++;
345 }
346 }
347
348 RowVectorXi sel_eeg = p_info.pick_types(false, true, false, defaultQStringList, p_exclude);
349 RowVectorXi sel_mag = p_info.pick_types(QString("mag"), false, false, defaultQStringList, p_exclude);
350 RowVectorXi sel_grad = p_info.pick_types(QString("grad"), false, false, defaultQStringList, p_exclude);
351
352 QStringList info_ch_names = p_info.ch_names;
353 QStringList ch_names_eeg, ch_names_mag, ch_names_grad;
354 for(qint32 i = 0; i < sel_eeg.size(); ++i)
355 ch_names_eeg << info_ch_names[sel_eeg(i)];
356 for(qint32 i = 0; i < sel_mag.size(); ++i)
357 ch_names_mag << info_ch_names[sel_mag(i)];
358 for(qint32 i = 0; i < sel_grad.size(); ++i)
359 ch_names_grad << info_ch_names[sel_grad(i)];
360
361 // This actually removes bad channels from the cov, which is not backward
362 // compatible, so let's leave all channels in
363 FiffCov cov_good = cov.pick_channels(info_ch_names, p_exclude);
364 QStringList ch_names = cov_good.names;
365
366 std::vector<qint32> idx_eeg, idx_mag, idx_grad;
367 for(qint32 i = 0; i < ch_names.size(); ++i)
368 {
369 if(ch_names_eeg.contains(ch_names[i]))
370 idx_eeg.push_back(i);
371 else if(ch_names_mag.contains(ch_names[i]))
372 idx_mag.push_back(i);
373 else if(ch_names_grad.contains(ch_names[i]))
374 idx_grad.push_back(i);
375 }
376
377 MatrixXd C(cov_good.data);
378
379 //Subtract number of found stim channels because they are still in C but not the idx_eeg, idx_mag or idx_grad
380 if((unsigned) C.rows() - iNoStimCh != idx_eeg.size() + idx_mag.size() + idx_grad.size()) {
381 printf("Error in FiffCov::regularize: Channel dimensions do not fit.\n");//ToDo Throw
382 }
383
384 QList<FiffProj> t_listProjs;
385 if(p_bProj)
386 {
387 t_listProjs = p_info.projs + cov_good.projs;
388 FiffProj::activate_projs(t_listProjs);
389 }
390
391 //Build regularization MAP
392 QMap<QString, QPair<double, std::vector<qint32> > > regData;
393 regData.insert("EEG", QPair<double, std::vector<qint32> >(p_fRegEeg, idx_eeg));
394 regData.insert("MAG", QPair<double, std::vector<qint32> >(p_fRegMag, idx_mag));
395 regData.insert("GRAD", QPair<double, std::vector<qint32> >(p_fRegGrad, idx_grad));
396
397 //
398 //Regularize
399 //
400 QMap<QString, QPair<double, std::vector<qint32> > >::Iterator it;
401 for(it = regData.begin(); it != regData.end(); ++it)
402 {
403 QString desc(it.key());
404 double reg = it.value().first;
405 std::vector<qint32> idx = it.value().second;
406
407 if(idx.size() == 0 || reg == 0.0)
408 printf("\tNothing to regularize within %s data.\n", desc.toUtf8().constData());
409 else
410 {
411 printf("\tRegularize %s: %f\n", desc.toUtf8().constData(), reg);
412 MatrixXd this_C(idx.size(), idx.size());
413 for(quint32 i = 0; i < idx.size(); ++i)
414 for(quint32 j = 0; j < idx.size(); ++j)
415 this_C(i,j) = cov_good.data(idx[i], idx[j]);
416
417 MatrixXd U;
418 qint32 ncomp;
419 if(p_bProj)
420 {
421 QStringList this_ch_names;
422 for(quint32 k = 0; k < idx.size(); ++k)
423 this_ch_names << ch_names[idx[k]];
424
425 MatrixXd P;
426 ncomp = FiffProj::make_projector(t_listProjs, this_ch_names, P); //ToDo: Synchronize with mne-python and debug
427
428 JacobiSVD<MatrixXd> svd(P, ComputeFullU);
429 //Sort singular values and singular vectors
430 VectorXd t_s = svd.singularValues();
431 MatrixXd t_U = svd.matrixU();
432 MNEMath::sort<double>(t_s, t_U);
433
434 U = t_U.block(0,0, t_U.rows(), t_U.cols()-ncomp);
435
436 if (ncomp > 0)
437 {
438 printf("\tCreated an SSP operator for %s (dimension = %d).\n", desc.toUtf8().constData(), ncomp);
439 this_C = U.transpose() * (this_C * U);
440 }
441 }
442
443 double sigma = this_C.diagonal().mean();
444 this_C.diagonal() = this_C.diagonal().array() + reg * sigma; // modify diag inplace
445 if(p_bProj && ncomp > 0)
446 this_C = U * (this_C * U.transpose());
447
448 for(qint32 i = 0; i < this_C.rows(); ++i)
449 for(qint32 j = 0; j < this_C.cols(); ++j)
450 C(idx[i],idx[j]) = this_C(i,j);
451 }
452 }
453
454 // Put data back in correct locations
455 RowVectorXi idx = FiffInfo::pick_channels(cov.names, info_ch_names, p_exclude);
456 for(qint32 i = 0; i < idx.size(); ++i)
457 for(qint32 j = 0; j < idx.size(); ++j)
458 cov.data(idx[i], idx[j]) = C(i, j);
459
460 return cov;
461}
462
463//=============================================================================================================
464
466{
467 if (this != &rhs) // protect against invalid self-assignment
468 {
469 kind = rhs.kind;
470 diag = rhs.diag;
471 dim = rhs.dim;
472 names = rhs.names;
473 data = rhs.data;
474 projs = rhs.projs;
475 bads = rhs.bads;
476 nfree = rhs.nfree;
477 eig = rhs.eig;
478 eigvec = rhs.eigvec;
479 }
480 // to support chained assignment operators (a=b=c), always return *this
481 return *this;
482}
483
484//=============================================================================================================
485
487 const MatrixXi &events,
488 const QList<int> &eventCodes,
489 float tmin,
490 float tmax,
491 float bmin,
492 float bmax,
493 bool doBaseline,
494 bool removeMean,
495 unsigned int ignoreMask,
496 float delay)
497{
498 FiffCov cov;
499 float sfreq = raw.info.sfreq;
500 int nchan = raw.info.nchan;
501
502 int minSamp = static_cast<int>(std::round(tmin * sfreq));
503 int maxSamp = static_cast<int>(std::round(tmax * sfreq));
504 int ns = maxSamp - minSamp + 1;
505 int delaySamp = static_cast<int>(std::round(delay * sfreq));
506
507 if (ns <= 0) {
508 qWarning() << "[FiffCov::compute_from_epochs] Invalid time window.";
509 return cov;
510 }
511
512 int bminSamp = 0, bmaxSamp = 0;
513 if (doBaseline) {
514 bminSamp = static_cast<int>(std::round(bmin * sfreq)) - minSamp;
515 bmaxSamp = static_cast<int>(std::round(bmax * sfreq)) - minSamp;
516 }
517
518 MatrixXd covAccum = MatrixXd::Zero(nchan, nchan);
519 VectorXd meanAccum = VectorXd::Zero(nchan);
520 int totalSamples = 0;
521 int nAccepted = 0;
522
523 for (int k = 0; k < events.rows(); ++k) {
524 int evFrom = events(k, 1) & ~static_cast<int>(ignoreMask);
525 int evTo = events(k, 2) & ~static_cast<int>(ignoreMask);
526
527 // Check if event matches any of the desired event codes
528 bool match = false;
529 for (int ec = 0; ec < eventCodes.size(); ++ec) {
530 if (evFrom == 0 && evTo == eventCodes[ec]) {
531 match = true;
532 break;
533 }
534 }
535 if (!match)
536 continue;
537
538 int evSample = events(k, 0);
539 int epochStart = evSample + delaySamp + minSamp;
540 int epochEnd = evSample + delaySamp + maxSamp;
541
542 if (epochStart < raw.first_samp || epochEnd > raw.last_samp)
543 continue;
544
545 MatrixXd epochData;
546 MatrixXd epochTimes;
547 if (!raw.read_raw_segment(epochData, epochTimes, epochStart, epochEnd))
548 continue;
549
550 // Baseline subtraction
551 if (doBaseline) {
552 int bminIdx = qMax(0, bminSamp);
553 int bmaxIdx = qMin(static_cast<int>(epochData.cols()) - 1, bmaxSamp);
554 if (bmaxIdx > bminIdx) {
555 int nBase = bmaxIdx - bminIdx;
556 for (int c = 0; c < nchan; ++c) {
557 double baseVal = epochData.row(c).segment(bminIdx, nBase).mean();
558 epochData.row(c).array() -= baseVal;
559 }
560 }
561 }
562
563 // Accumulate
564 if (removeMean) {
565 VectorXd epochMean = epochData.rowwise().mean();
566 meanAccum += epochMean * static_cast<double>(ns);
567 }
568 covAccum += epochData * epochData.transpose();
569 totalSamples += ns;
570 nAccepted++;
571 }
572
573 if (totalSamples < 2) {
574 qWarning() << "[FiffCov::compute_from_epochs] Not enough data.";
575 return cov;
576 }
577
578 if (removeMean) {
579 VectorXd grandMean = meanAccum / static_cast<double>(totalSamples);
580 cov.data = (covAccum / static_cast<double>(totalSamples - 1))
581 - (grandMean * grandMean.transpose()) * (static_cast<double>(totalSamples) / (totalSamples - 1));
582 } else {
583 cov.data = covAccum / static_cast<double>(totalSamples - 1);
584 }
585
587 cov.dim = nchan;
588 cov.names = raw.info.ch_names;
589 cov.nfree = totalSamples - 1;
590 cov.bads = raw.info.bads;
591 cov.projs = raw.info.projs;
592
593 qInfo() << "[FiffCov::compute_from_epochs] Computed:" << nchan << "channels,"
594 << nAccepted << "epochs," << totalSamples << "total samples.";
595
596 return cov;
597}
598
599//=============================================================================================================
600
601bool FiffCov::save(const QString &fileName) const
602{
603 if (fileName.isEmpty()) {
604 qWarning() << "[FiffCov::save] Output file not specified.";
605 return false;
606 }
607
608 QFile file(fileName);
610 if (!pStream) {
611 qWarning() << "[FiffCov::save] Cannot open" << fileName;
612 return false;
613 }
614
615 pStream->start_block(FIFFB_MEAS);
616 pStream->write_id(FIFF_BLOCK_ID);
617 pStream->write_cov(*this);
618 pStream->end_block(FIFFB_MEAS);
619 pStream->end_file();
620
621 qInfo() << "[FiffCov::save] Saved covariance matrix to" << fileName;
622 return true;
623}
624
625//=============================================================================================================
626
627FiffCov FiffCov::computeGrandAverage(const QList<FiffCov> &covs)
628{
629 FiffCov grandCov;
630
631 if (covs.isEmpty()) {
632 qWarning() << "[FiffCov::computeGrandAverage] No covariance matrices provided.";
633 return grandCov;
634 }
635
636 grandCov = covs[0];
637 MatrixXd sumCov = grandCov.data * static_cast<double>(grandCov.nfree);
638 int totalNfree = grandCov.nfree;
639
640 for (int k = 1; k < covs.size(); ++k) {
641 if (covs[k].dim != grandCov.dim) {
642 qWarning() << "[FiffCov::computeGrandAverage] Dimension mismatch.";
643 return FiffCov();
644 }
645 sumCov += covs[k].data * static_cast<double>(covs[k].nfree);
646 totalNfree += covs[k].nfree;
647 }
648
649 grandCov.data = sumCov / static_cast<double>(totalNfree);
650 grandCov.nfree = totalNfree;
651
652 return grandCov;
653}
#define FIFFV_MNE_NOISE_COV
#define FIFFV_STIM_CH
FiffStream class declaration.
FiffInfoBase class declaration.
FiffRawData class declaration.
FiffDirNode class declaration, which provides fiff dir tree processing methods.
FiffCov class declaration.
Header file describing the numerical values used in fif files.
#define FIFFB_MEAS
Definition fiff_file.h:362
#define FIFF_BLOCK_ID
Definition fiff_file.h:326
MNEMath class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Definition buildinfo.h:45
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
Eigen::MatrixXd eigvec
Definition fiff_cov.h:256
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
static FiffCov compute_from_epochs(const FiffRawData &raw, const Eigen::MatrixXi &events, const QList< int > &eventCodes, float tmin, float tmax, float bmin=0.0f, float bmax=0.0f, bool doBaseline=false, bool removeMean=true, unsigned int ignoreMask=0, float delay=0.0f)
Definition fiff_cov.cpp:486
FiffCov pick_channels(const QStringList &p_include=defaultQStringList, const QStringList &p_exclude=defaultQStringList)
Definition fiff_cov.cpp:151
fiff_int_t kind
Definition fiff_cov.h:246
static FiffCov computeGrandAverage(const QList< FiffCov > &covs)
Definition fiff_cov.cpp:627
FiffCov & operator=(const FiffCov &rhs)
Definition fiff_cov.cpp:465
QStringList bads
Definition fiff_cov.h:253
QStringList names
Definition fiff_cov.h:250
Eigen::VectorXd eig
Definition fiff_cov.h:255
Eigen::MatrixXd data
Definition fiff_cov.h:251
FiffCov prepare_noise_cov(const FiffInfo &p_info, const QStringList &p_chNames) const
Definition fiff_cov.cpp:179
bool save(const QString &fileName) const
Definition fiff_cov.cpp:601
FIFF measurement file information.
Definition fiff_info.h:85
qint32 make_projector(Eigen::MatrixXd &proj) const
Definition fiff_info.h:266
QList< FiffProj > projs
Definition fiff_info.h:256
static Eigen::RowVectorXi pick_channels(const QStringList &ch_names, const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList)
QList< FiffChInfo > chs
Eigen::RowVectorXi pick_types(const QString meg, bool eeg=false, bool stim=false, const QStringList &include=defaultQStringList, const QStringList &exclude=defaultQStringList) const
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)
static void activate_projs(QList< FiffProj > &p_qListFiffProj)
FIFF raw measurement data.
bool read_raw_segment(Eigen::MatrixXd &data, Eigen::MatrixXd &times, fiff_int_t from=-1, fiff_int_t to=-1, const Eigen::RowVectorXi &sel=defaultRowVectorXi, bool do_debug=false) const
FIFF File I/O routines.
static FiffStream::SPtr start_file(QIODevice &p_IODevice)
QSharedPointer< FiffStream > SPtr
static Eigen::VectorXi sort(Eigen::Matrix< T, Eigen::Dynamic, 1 > &v, bool desc=true)
Definition mnemath.h:499
static void get_whitener(Eigen::MatrixXd &A, bool pca, QString ch_type, Eigen::VectorXd &eig, Eigen::MatrixXd &eigvec)