v2.0.0
Loading...
Searching...
No Matches
mne_cov_matrix.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "mne_cov_matrix.h"
42#include "mne_sss_data.h"
43#include "mne_proj_item.h"
44#include "mne_proj_op.h"
45
46#include <Eigen/Core>
47#include <Eigen/Eigenvalues>
49#include <fiff/fiff_ch_info.h>
50#include <fiff/fiff_constants.h>
51#include <fiff/fiff_stream.h>
52#include <fiff/fiff_tag.h>
53
54#include <QFile>
55#include <QDebug>
56
57//=============================================================================================================
58// USED NAMESPACES
59//=============================================================================================================
60
61using namespace Eigen;
62using namespace FIFFLIB;
63using namespace MNELIB;
64
65constexpr int FAIL = -1;
66constexpr int OK = 0;
67
68//============================= mne_decompose.c =============================
69
70int mne_decompose_eigen (const VectorXd& mat,
71 VectorXd& lambda,
72 Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>& vectors,
73 int dim)
74/*
75 * Compute the eigenvalue decomposition of
76 * a symmetric matrix using the LAPACK routines
77 *
78 * 'mat' contains the lower triangle of the matrix
79 */
80{
81 int np = dim*(dim+1)/2;
82 int maxi;
83 double scale;
84
85// idamax workaround begin
86 maxi = 0;
87 for(int i = 0; i < np; ++i)
88 if (std::fabs(mat[i]) > std::fabs(mat[maxi]))
89 maxi = i;
90// idamax workaround end
91
92 scale = 1.0/mat[maxi];
93
94// dspev workaround begin
95 MatrixXd dmat_full = MatrixXd::Zero(dim,dim);
96 int idx = 0;
97 for (int i = 0; i < dim; ++i) {
98 for(int j = 0; j <= i; ++j) {
99 double val = mat[idx]*scale;
100 dmat_full(i,j) = val;
101 dmat_full(j,i) = val;
102 ++idx;
103 }
104 }
105 SelfAdjointEigenSolver<MatrixXd> es;
106 es.compute(dmat_full);
107// dspev workaround end
108
109 scale = 1.0/scale;
110 lambda = es.eigenvalues() * scale;
111 vectors = es.eigenvectors().transpose().cast<float>();
112
113 return 0;
114}
115
116//=============================================================================================================
117// DEFINE MEMBER METHODS
118//=============================================================================================================
119
121 int p_ncov,
122 const QStringList& p_names,
123 const VectorXd& p_cov,
124 const VectorXd& p_cov_diag,
125 FiffSparseMatrix* p_cov_sparse)
126:kind(p_kind)
127,ncov(p_ncov)
128,nfree(1)
129,nproj(0)
130,nzero(0)
131,names(p_names)
132,cov(p_cov)
133,cov_diag(p_cov_diag)
134,cov_sparse(p_cov_sparse)
135,proj(nullptr)
136,sss(nullptr)
137,nbad(0)
138{
139}
140
141//=============================================================================================================
142
146
147//=============================================================================================================
148
149std::unique_ptr<MNECovMatrix> MNECovMatrix::read(const QString& name, int kind)
150{
151 QFile file(name);
152 FiffStream::SPtr stream(new FiffStream(&file));
153
154 FiffTag::UPtr t_pTag;
155 QList<FiffDirNode::SPtr> nodes;
156 FiffDirNode::SPtr covnode;
157
158 QStringList names;
159 int nnames = 0;
160 Eigen::VectorXd cov;
161 Eigen::VectorXd cov_diag;
162 std::unique_ptr<FiffSparseMatrix> cov_sparse_owner;
163 Eigen::VectorXd lambda;
164 Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> eigen;
165 MatrixXf tmp_eigen;
166 QStringList bads;
167 int nbad = 0;
168 int ncov = 0;
169 int nfree = 1;
170 std::unique_ptr<MNECovMatrix> res;
171
172 int k,p,nn;
173 const float *f;
174 const double *d;
175 std::unique_ptr<MNEProjOp> op;
176 std::unique_ptr<MNESssData> sss;
177
178 if (!stream->open())
179 return nullptr;
180
181 nodes = stream->dirtree()->dir_tree_find(FIFFB_MNE_COV);
182
183 if (nodes.size() == 0) {
184 qWarning("No covariance matrix available in %s", name.toUtf8().data());
185 stream->close();
186 return nullptr;
187 }
188 /*
189 * Locate the desired matrix
190 */
191 for (k = 0; k < nodes.size(); ++k) {
192 if (!nodes[k]->find_tag(stream, FIFF_MNE_COV_KIND, t_pTag))
193 continue;
194
195 if (*t_pTag->toInt() == kind) {
196 covnode = nodes[k];
197 break;
198 }
199 }
200 if (covnode->isEmpty()) {
201 qWarning("Desired covariance matrix not found from %s", name.toUtf8().data());
202 stream->close();
203 return nullptr;
204 }
205 /*
206 * Read the data
207 */
208 if (!nodes[k]->find_tag(stream, FIFF_MNE_COV_DIM, t_pTag)) {
209 stream->close();
210 return nullptr;
211 }
212 ncov = *t_pTag->toInt();
213
214 if (nodes[k]->find_tag(stream, FIFF_MNE_COV_NFREE, t_pTag)) {
215 nfree = *t_pTag->toInt();
216 }
217 if (covnode->find_tag(stream, FIFF_MNE_ROW_NAMES, t_pTag)) {
218 names = FiffStream::split_name_list(t_pTag->toString());
219 nnames = names.size();
220 if (nnames != ncov) {
221 qCritical("Incorrect number of channel names for a covariance matrix");
222 stream->close();
223 return nullptr;
224 }
225 }
226 if (!nodes[k]->find_tag(stream, FIFF_MNE_COV, t_pTag)) {
227 if (!nodes[k]->find_tag(stream, FIFF_MNE_COV_DIAG, t_pTag)) {
228 stream->close();
229 return nullptr;
230 } else {
231 if (t_pTag->getType() == FIFFT_DOUBLE) {
232 cov_diag.resize(ncov);
233 d = t_pTag->toDouble();
234 for (p = 0; p < ncov; p++)
235 cov_diag[p] = d[p];
236 /*
237 * Check for all-zero data
238 */
239 if (cov_diag.sum() == 0.0) {
240 qCritical("Sum of covariance matrix elements is zero!");
241 stream->close();
242 return nullptr;
243 }
244 } else if (t_pTag->getType() == FIFFT_FLOAT) {
245 cov_diag.resize(ncov);
246 f = t_pTag->toFloat();
247 for (p = 0; p < ncov; p++)
248 cov_diag[p] = f[p];
249 } else {
250 qWarning("Illegal data type for covariance matrix");
251 stream->close();
252 return nullptr;
253 }
254 }
255 } else {
256 nn = ncov * (ncov + 1) / 2;
257 if (t_pTag->getType() == FIFFT_DOUBLE) {
258 cov.resize(nn);
259 d = t_pTag->toDouble();
260 for (p = 0; p < nn; p++)
261 cov[p] = d[p];
262 if (cov.sum() == 0.0) {
263 qCritical("Sum of covariance matrix elements is zero!");
264 stream->close();
265 return nullptr;
266 }
267 } else if (t_pTag->getType() == FIFFT_FLOAT) {
268 cov.resize(nn);
269 f = t_pTag->toFloat();
270 for (p = 0; p < nn; p++)
271 cov[p] = f[p];
272 } else {
273 cov_sparse_owner = FiffSparseMatrix::fiff_get_float_sparse_matrix(t_pTag);
274 if (!cov_sparse_owner) {
275 stream->close();
276 return nullptr;
277 }
278 }
279
280 if (nodes[k]->find_tag(stream, FIFF_MNE_COV_EIGENVALUES, t_pTag)) {
281 const double *lambda_data = static_cast<const double *>(t_pTag->toDouble());
282 lambda = Eigen::Map<const Eigen::VectorXd>(lambda_data, ncov);
283 if (nodes[k]->find_tag(stream, FIFF_MNE_COV_EIGENVECTORS, t_pTag)) {
284 stream->close();
285 return nullptr;
286 }
287
288 tmp_eigen = t_pTag->toFloatMatrix().transpose();
289 eigen.resize(tmp_eigen.rows(), tmp_eigen.cols());
290 for (int r = 0; r < tmp_eigen.rows(); ++r)
291 for (int c = 0; c < tmp_eigen.cols(); ++c)
292 eigen(r, c) = tmp_eigen(r, c);
293 }
294 /*
295 * Read the optional projection operator
296 */
297 op = MNEProjOp::read_from_node(stream, nodes[k]);
298 if (!op) {
299 stream->close();
300 return nullptr;
301 }
302 /*
303 * Read the optional SSS data
304 */
305 sss = MNESssData::read_from_node(stream, nodes[k]);
306 if (!sss) {
307 stream->close();
308 return nullptr;
309 }
310 /*
311 * Read the optional bad channel list
312 */
313 bads = stream->read_bad_channels(nodes[k]);
314 nbad = bads.size();
315 }
316 if (cov_sparse_owner)
317 res = create_sparse(kind, ncov, names, cov_sparse_owner.release());
318 else if (cov.size() > 0)
319 res = create_dense(kind, ncov, names, cov);
320 else if (cov_diag.size() > 0)
322 else {
323 qCritical("MNECovMatrix::read : covariance matrix data is not defined.");
324 stream->close();
325 return nullptr;
326 }
327 res->eigen = std::move(eigen);
328 res->lambda = std::move(lambda);
329 res->nfree = nfree;
330 res->bads = bads;
331 res->nbad = nbad;
332 /*
333 * Count the non-zero eigenvalues
334 */
335 if (res->lambda.size() > 0) {
336 res->nzero = 0;
337 for (k = 0; k < res->ncov; k++, res->nzero++)
338 if (res->lambda[k] > 0)
339 break;
340 }
341
342 if (op && op->nitems > 0) {
343 res->proj = std::move(op);
344 }
345 if (sss && sss->comp_info.size() > 0 && sss->job != FIFFV_SSS_JOB_NOTHING) {
346 res->sss = std::move(sss);
347 }
348
349 stream->close();
350 return res;
351}
352
353//=============================================================================================================
354
355std::unique_ptr<MNECovMatrix> MNECovMatrix::dup() const
356{
357 auto res = cov_diag.size() > 0
358 ? create(kind,ncov,names,VectorXd(),VectorXd(cov_diag))
359 : create(kind,ncov,names,VectorXd(cov),VectorXd());
360 /*
361 * Duplicate additional items
362 */
363 if (ch_class.size() > 0) {
364 res->ch_class = ch_class;
365 }
366 res->bads = bads;
367 res->nbad = nbad;
368 res->proj = proj ? proj->dup() : nullptr;
369 if (sss)
370 res->sss = std::make_unique<MNESssData>(*sss);
371
372 return res;
373}
374
375//=============================================================================================================
376
378{
379 return cov_diag.size() > 0;
380}
381
382//=============================================================================================================
383
385/*
386 * Calculate the inverse square roots for whitening
387 */
388{
389 const VectorXd& src = lambda.size() > 0 ? lambda : cov_diag;
390 int k;
391
392 if (src.size() == 0) {
393 qCritical("Covariance matrix is not diagonal or not decomposed.");
394 return FAIL;
395 }
396 inv_lambda.resize(ncov);
397 for (k = 0; k < ncov; k++) {
398 if (src[k] <= 0.0)
399 inv_lambda[k] = 0.0;
400 else
401 inv_lambda[k] = 1.0/sqrt(src[k]);
402 }
403 return OK;
404}
405
406//=============================================================================================================
407
408int MNECovMatrix::condition(float rank_threshold, int use_rank)
409{
410 VectorXd scale_vec;
411 VectorXd cov_local;
412 VectorXd lambda_local;
413 Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> local_eigen;
414 MatrixXd data1;
415 double magscale,gradscale,eegscale;
416 int nmag,ngrad,neeg,nok;
417 int j,k;
418 int res = FAIL;
419
420 if (cov_diag.size() > 0)
421 return OK;
422 if (ch_class.size() == 0) {
423 qCritical("Channels not classified. Rank cannot be determined.");
424 return FAIL;
425 }
426 magscale = gradscale = eegscale = 0.0;
427 nmag = ngrad = neeg = 0;
428 for (k = 0; k < ncov; k++) {
429 if (ch_class[k] == MNE_COV_CH_MEG_MAG) {
430 magscale += this->cov[lt_packed_index(k,k)]; nmag++;
431 }
432 else if (ch_class[k] == MNE_COV_CH_MEG_GRAD) {
433 gradscale += this->cov[lt_packed_index(k,k)]; ngrad++;
434 }
435 else if (ch_class[k] == MNE_COV_CH_EEG) {
436 eegscale += this->cov[lt_packed_index(k,k)]; neeg++;
437 }
438#ifdef DEBUG
439 fprintf(stdout,"%d ",ch_class[k]);
440#endif
441 }
442#ifdef DEBUG
443 fprintf(stdout,"\n");
444#endif
445 if (nmag > 0)
446 magscale = magscale > 0.0 ? sqrt(nmag/magscale) : 0.0;
447 if (ngrad > 0)
448 gradscale = gradscale > 0.0 ? sqrt(ngrad/gradscale) : 0.0;
449 if (neeg > 0)
450 eegscale = eegscale > 0.0 ? sqrt(neeg/eegscale) : 0.0;
451#ifdef DEBUG
452 fprintf(stdout,"%d %g\n",nmag,magscale);
453 fprintf(stdout,"%d %g\n",ngrad,gradscale);
454 fprintf(stdout,"%d %g\n",neeg,eegscale);
455#endif
456 scale_vec.resize(ncov);
457 for (k = 0; k < ncov; k++) {
459 scale_vec[k] = magscale;
460 else if (ch_class[k] == MNE_COV_CH_MEG_GRAD)
461 scale_vec[k] = gradscale;
462 else if (ch_class[k] == MNE_COV_CH_EEG)
463 scale_vec[k] = eegscale;
464 else
465 scale_vec[k] = 1.0;
466 }
467 cov_local.resize(ncov*(ncov+1)/2);
468 lambda_local.resize(ncov);
469 local_eigen.resize(ncov,ncov);
470 for (j = 0; j < ncov; j++)
471 for (k = 0; k <= j; k++)
472 cov_local[lt_packed_index(j,k)] = this->cov[lt_packed_index(j,k)]*scale_vec[j]*scale_vec[k];
473 if (mne_decompose_eigen(cov_local,lambda_local,local_eigen,ncov) == 0) {
474#ifdef DEBUG
475 for (k = 0; k < ncov; k++)
476 fprintf(stdout,"%g ",lambda_local[k]/lambda_local[ncov-1]);
477 fprintf(stdout,"\n");
478#endif
479 nok = 0;
480 for (k = ncov-1; k >= 0; k--) {
481 if (lambda_local[k] >= rank_threshold*lambda_local[ncov-1])
482 nok++;
483 else
484 break;
485 }
486 qInfo("\n\tEstimated covariance matrix rank = %d (%g)\n",nok,lambda_local[ncov-nok]/lambda_local[ncov-1]);
487 if (use_rank > 0 && use_rank < nok) {
488 nok = use_rank;
489 qInfo("\tUser-selected covariance matrix rank = %d (%g)\n",nok,lambda_local[ncov-nok]/lambda_local[ncov-1]);
490 }
491 /*
492 * Put it back together
493 */
494 for (j = 0; j < ncov-nok; j++)
495 lambda_local[j] = 0.0;
496 data1.resize(ncov,ncov);
497 for (j = 0; j < ncov; j++) {
498#ifdef DEBUG
499 mne_print_vector(stdout,nullptr,local_eigen.row(j).data(),ncov);
500#endif
501 for (k = 0; k < ncov; k++)
502 data1(j,k) = sqrt(lambda_local[j])*local_eigen(j,k);
503 }
504 MatrixXd data2 = data1.transpose() * data1;
505#ifdef DEBUG
506 qInfo(">>>\n");
507 for (j = 0; j < ncov; j++)
508 mne_print_dvector(stdout,nullptr,data2.row(j).data(),ncov);
509 qInfo(">>>\n");
510#endif
511 /*
512 * Scale back
513 */
514 for (k = 0; k < ncov; k++)
515 if (scale_vec[k] > 0.0)
516 scale_vec[k] = 1.0/scale_vec[k];
517 for (j = 0; j < ncov; j++)
518 for (k = 0; k <= j; k++)
519 if (this->cov[lt_packed_index(j,k)] != 0.0)
520 this->cov[lt_packed_index(j,k)] = scale_vec[j]*scale_vec[k]*data2(j,k);
521 res = nok;
522 }
523 return res;
524}
525
526//=============================================================================================================
527
528int MNECovMatrix::decompose_eigen_small(float p_small, int use_rank)
529/*
530 * Do the eigenvalue decomposition
531 */
532{
533 int k,p,rank;
534 float rank_threshold = 1e-6;
535
536 if (p_small < 0)
537 p_small = 1.0;
538
539 if (cov_diag.size() > 0)
540 return add_inv();
541 if (lambda.size() > 0 && eigen.size() > 0) {
542 qInfo("\n\tEigenvalue decomposition had been precomputed.\n");
543 nzero = 0;
544 for (k = 0; k < ncov; k++, nzero++)
545 if (lambda[k] > 0)
546 break;
547 }
548 else {
549 lambda.resize(0);
550 eigen.resize(0,0);
551
552 if ((rank = condition(rank_threshold,use_rank)) < 0)
553 return FAIL;
554
555 lambda.resize(ncov);
556 eigen.resize(ncov,ncov);
558 lambda.resize(0);
559 eigen.resize(0,0);
560 return FAIL;
561 }
562 nzero = ncov - rank;
563 for (k = 0; k < nzero; k++)
564 lambda[k] = 0.0;
565 /*
566 * Find which eigenvectors correspond to EEG/MEG
567 */
568 {
569 float meglike,eeglike;
570 int nmeg,neeg;
571
572 nmeg = neeg = 0;
573 for (k = nzero; k < ncov; k++) {
574 meglike = eeglike = 0.0;
575 for (p = 0; p < ncov; p++) {
576 if (ch_class[p] == MNE_COV_CH_EEG)
577 eeglike += std::fabs(eigen(k,p));
579 meglike += std::fabs(eigen(k,p));
580 }
581 if (meglike > eeglike)
582 nmeg++;
583 else
584 neeg++;
585 }
586 qInfo("\t%d MEG and %d EEG-like channels remain in the whitened data\n",nmeg,neeg);
587 }
588 }
589 return add_inv();
590}
591
592//=============================================================================================================
593
595
596{
597 return decompose_eigen_small(-1.0,-1);
598}
599
600//=============================================================================================================
601
603
604{
605 if (j >= k)
606 return k + j*(j+1)/2;
607 else
608 return j + k*(k+1)/2;
609}
610
611//=============================================================================================================
612
613int MNECovMatrix::classify_channels(const QList<FiffChInfo>& chs, int nchan)
614{
615 int k,p;
616 FiffChInfo ch;
617
618 if (chs.isEmpty()) {
619 qCritical("Channel information not available in classify_channels");
620 return FAIL;
621 }
622 ch_class.resize(ncov);
623 for (k = 0; k < ncov; k++) {
625 for (p = 0; p < nchan; p++) {
626 if (QString::compare(chs[p].ch_name,names[k]) == 0) {
627 ch = chs[p];
628 if (ch.kind == FIFFV_MEG_CH) {
629 if (ch.unit == FIFF_UNIT_T)
631 else
633 }
634 else if (ch.kind == FIFFV_EEG_CH)
636 break;
637 }
638 }
639 }
640 return OK;
641}
642
643//=============================================================================================================
644
645int MNECovMatrix::whiten_vector(Eigen::Ref<Eigen::VectorXf> data, Eigen::Ref<Eigen::VectorXf> whitened_data, int nchan) const
646{
647 if (ncov != nchan) {
648 qWarning("Incompatible covariance matrix. Cannot whiten the data.");
649 return FAIL;
650 }
651 const double *inv = inv_lambda.data();
652 if (is_diag()) {
653 for (int k = 0; k < nchan; k++)
654 whitened_data[k] = data[k]*inv[k];
655 }
656 else {
657 Eigen::VectorXf tmp(nchan);
658 for (int k = nzero; k < nchan; k++)
659 tmp[k] = eigen.row(k).dot(data.cast<float>());
660 for (int k = 0; k < nzero; k++)
661 whitened_data[k] = 0.0;
662 for (int k = nzero; k < nchan; k++)
663 whitened_data[k] = tmp[k]*inv[k];
664 }
665 return OK;
666}
667
668//=============================================================================================================
669
670void MNECovMatrix::regularize(const Eigen::Vector3f& regs)
671{
672 int j;
673 float sums[3],nn[3];
674 int nkind = 3;
675
676 if (cov.size() == 0 || ch_class.size() == 0)
677 return;
678
679 for (j = 0; j < nkind; j++) {
680 sums[j] = 0.0;
681 nn[j] = 0;
682 }
683 for (j = 0; j < ncov; j++) {
684 if (ch_class[j] >= 0) {
685 sums[ch_class[j]] += cov[lt_packed_index(j,j)];
686 nn[ch_class[j]]++;
687 }
688 }
689 qInfo("Average noise-covariance matrix diagonals:");
690 for (j = 0; j < nkind; j++) {
691 if (nn[j] > 0) {
692 sums[j] = sums[j]/nn[j];
693 if (j == MNE_COV_CH_MEG_MAG)
694 qInfo("\tMagnetometers : %-7.2f fT reg = %-6.2f",1e15*sqrt(sums[j]),regs[j]);
695 else if (j == MNE_COV_CH_MEG_GRAD)
696 qInfo("\tPlanar gradiometers : %-7.2f fT/cm reg = %-6.2f",1e13*sqrt(sums[j]),regs[j]);
697 else
698 qInfo("\tEEG : %-7.2f uV reg = %-6.2f",1e6*sqrt(sums[j]),regs[j]);
699 sums[j] = regs[j]*sums[j];
700 }
701 }
702 for (j = 0; j < ncov; j++)
703 if (ch_class[j] >= 0)
704 cov[lt_packed_index(j,j)] += sums[ch_class[j]];
705
706 qInfo("Noise-covariance regularized as requested.");
707}
708
709//=============================================================================================================
710
712{
713 int k,p;
714 if (cov.size() == 0)
715 return;
716
717 cov_diag.resize(ncov);
718
719 for (k = p = 0; k < ncov; k++) {
720 cov_diag[k] = cov[p];
721 p = p + k + 2;
722 }
723 cov.resize(0);
724
725 lambda.resize(0);
726 eigen.resize(0,0);
727}
728
729//=============================================================================================================
730
731std::unique_ptr<MNECovMatrix> MNECovMatrix::pick_chs_omit(const QStringList& new_names,
732 int new_ncov,
733 int omit_meg_eeg,
734 const QList<FiffChInfo>& chs) const
735{
736 int j,k;
737 Eigen::VectorXd cov_local;
738 Eigen::VectorXd cov_diag_local;
739 QStringList picked_names;
740 int from,to;
741 std::unique_ptr<MNECovMatrix> res;
742
743 if (new_ncov == 0) {
744 qCritical("No channels specified for picking in pick_chs_omit");
745 return nullptr;
746 }
747 if (names.isEmpty()) {
748 qCritical("No names in covariance matrix. Cannot do picking.");
749 return nullptr;
750 }
751 Eigen::VectorXi pickVec = Eigen::VectorXi::Constant(new_ncov, -1);
752 for (j = 0; j < new_ncov; j++)
753 for (k = 0; k < ncov; k++)
754 if (QString::compare(names[k],new_names[j]) == 0) {
755 pickVec[j] = k;
756 break;
757 }
758 for (j = 0; j < new_ncov; j++) {
759 if (pickVec[j] < 0) {
760 qWarning("All desired channels not found in the covariance matrix (at least missing %s).", new_names[j].toUtf8().constData());
761 return nullptr;
762 }
763 }
764 Eigen::VectorXi isMegVec;
765 if (omit_meg_eeg) {
766 isMegVec.resize(new_ncov);
767 if (!chs.isEmpty()) {
768 for (j = 0; j < new_ncov; j++)
769 if (chs[j].kind == FIFFV_MEG_CH)
770 isMegVec[j] = true;
771 else
772 isMegVec[j] = false;
773 }
774 else {
775 for (j = 0; j < new_ncov; j++)
776 if (new_names[j].startsWith("MEG"))
777 isMegVec[j] = true;
778 else
779 isMegVec[j] = false;
780 }
781 }
782 if (cov_diag.size() > 0) {
783 cov_diag_local.resize(new_ncov);
784 for (j = 0; j < new_ncov; j++) {
785 cov_diag_local[j] = cov_diag[pickVec[j]];
786 picked_names.append(names[pickVec[j]]);
787 }
788 }
789 else {
790 cov_local.resize(new_ncov*(new_ncov+1)/2);
791 for (j = 0; j < new_ncov; j++) {
792 picked_names.append(names[pickVec[j]]);
793 for (k = 0; k <= j; k++) {
794 from = lt_packed_index(pickVec[j],pickVec[k]);
795 to = lt_packed_index(j,k);
796 if (to < 0 || to > new_ncov*(new_ncov+1)/2-1) {
797 qCritical("Wrong destination index in pick_chs_omit : %d %d %d",j,k,to);
798 return nullptr;
799 }
800 if (from < 0 || from > ncov*(ncov+1)/2-1) {
801 qCritical("Wrong source index in pick_chs_omit : %d %d %d",pickVec[j],pickVec[k],from);
802 return nullptr;
803 }
804 cov_local[to] = cov[from];
805 if (omit_meg_eeg)
806 if (isMegVec[j] != isMegVec[k])
807 cov_local[to] = 0.0;
808 }
809 }
810 }
811
812 res = MNECovMatrix::create(kind,new_ncov,picked_names,cov_local,cov_diag_local);
813
814 res->bads = bads;
815 res->nbad = nbad;
816 res->proj = proj ? proj->dup() : nullptr;
817 res->sss = sss ? std::make_unique<MNESssData>(*sss) : nullptr;
818
819 if (ch_class.size() > 0) {
820 res->ch_class.resize(res->ncov);
821 for (k = 0; k < res->ncov; k++)
822 res->ch_class[k] = ch_class[pickVec[k]];
823 }
824 return res;
825}
constexpr int FAIL
constexpr int OK
#define FIFFT_DOUBLE
Definition fiff_file.h:233
#define FIFFT_FLOAT
Definition fiff_file.h:232
#define FIFFV_SSS_JOB_NOTHING
Definition fiff_file.h:538
FiffStream class declaration.
FiffTag class declaration, which provides fiff tag I/O and processing methods.
FiffSparseMatrix class declaration.
Fiff constants.
#define FIFF_MNE_COV_KIND
#define FIFFV_EEG_CH
#define FIFF_MNE_COV
#define FIFF_MNE_COV_DIAG
#define FIFF_MNE_ROW_NAMES
#define FIFF_MNE_COV_EIGENVALUES
#define FIFFV_MEG_CH
#define FIFFB_MNE_COV
#define FIFF_MNE_COV_NFREE
#define FIFF_MNE_COV_EIGENVECTORS
#define FIFF_UNIT_T
#define FIFF_MNE_COV_DIM
FiffChInfo class declaration.
MNECovMatrix class declaration.
#define MNE_COV_CH_MEG_MAG
#define MNE_COV_CH_EEG
#define MNE_COV_CH_MEG_GRAD
#define MNE_COV_CH_UNKNOWN
MNE SSS Data (MNESssData) class declaration.
MNEProjItem class declaration.
MNEProjOp class declaration.
int mne_decompose_eigen(const VectorXd &mat, VectorXd &lambda, Eigen::Matrix< float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &vectors, int dim)
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Channel info descriptor.
QSharedPointer< FiffDirNode > SPtr
FIFF sparse matrix storage backed by Eigen.
static FiffSparseMatrix::UPtr fiff_get_float_sparse_matrix(const FIFFLIB::FiffTag::UPtr &tag)
FIFF File I/O routines.
QSharedPointer< FiffStream > SPtr
static QStringList split_name_list(QString p_sNameList)
std::unique_ptr< FiffTag > UPtr
Definition fiff_tag.h:158
Eigen::VectorXd inv_lambda
std::unique_ptr< FIFFLIB::FiffSparseMatrix > cov_sparse
Eigen::VectorXd lambda
int condition(float rank_threshold, int use_rank)
static std::unique_ptr< MNECovMatrix > create_sparse(int kind, int ncov, const QStringList &names, FIFFLIB::FiffSparseMatrix *cov_sparse)
void regularize(const Eigen::Vector3f &regs)
int decompose_eigen_small(float p_small, int use_rank)
int classify_channels(const QList< FIFFLIB::FiffChInfo > &chs, int nchan)
std::unique_ptr< MNEProjOp > proj
static std::unique_ptr< MNECovMatrix > create_diag(int kind, int ncov, const QStringList &names, const Eigen::VectorXd &cov_diag)
static std::unique_ptr< MNECovMatrix > create_dense(int kind, int ncov, const QStringList &names, const Eigen::VectorXd &cov)
std::unique_ptr< MNECovMatrix > dup() const
MNECovMatrix(int p_kind, int p_ncov, const QStringList &p_names, const Eigen::VectorXd &p_cov, const Eigen::VectorXd &p_cov_diag, FIFFLIB::FiffSparseMatrix *p_cov_sparse)
std::unique_ptr< MNESssData > sss
Eigen::VectorXd cov
Eigen::VectorXd cov_diag
static std::unique_ptr< MNECovMatrix > read(const QString &name, int kind)
static std::unique_ptr< MNECovMatrix > create(int kind, int ncov, const QStringList &names, const Eigen::VectorXd &cov, const Eigen::VectorXd &cov_diag)
Eigen::Matrix< float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > eigen
static int lt_packed_index(int j, int k)
std::unique_ptr< MNECovMatrix > pick_chs_omit(const QStringList &new_names, int new_ncov, int omit_meg_eeg, const QList< FIFFLIB::FiffChInfo > &chs) const
Eigen::VectorXi ch_class
int whiten_vector(Eigen::Ref< Eigen::VectorXf > data, Eigen::Ref< Eigen::VectorXf > whitened_data, int nchan) const
static std::unique_ptr< MNEProjOp > read_from_node(FIFFLIB::FiffStream::SPtr &stream, const FIFFLIB::FiffDirNode::SPtr &start)
Read all linear projection items from a FIFF tree node.
static std::unique_ptr< MNESssData > read_from_node(QSharedPointer< FIFFLIB::FiffStream > &stream, const QSharedPointer< FIFFLIB::FiffDirNode > &start)