v2.0.0
Loading...
Searching...
No Matches
mne_ctf_comp_data_set.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
42#include "mne_ctf_comp_data.h"
43
44#include "mne_types.h"
45
47
48#include <fiff/fiff_types.h>
49
50#include <Eigen/Core>
51
52//=============================================================================================================
53// QT INCLUDES
54//=============================================================================================================
55
56#include <QFile>
57
58#ifndef TRUE
59#define TRUE 1
60#endif
61
62#ifndef FALSE
63#define FALSE 0
64#endif
65
66#ifndef FAIL
67#define FAIL -1
68#endif
69
70#ifndef OK
71#define OK 0
72#endif
73
74
75
76//=============================================================================================================
77// USED NAMESPACES
78//=============================================================================================================
79
80using namespace Eigen;
81using namespace FIFFLIB;
82using namespace MNELIB;
83
84//============================= mne_read_forward_solution.c =============================
85
86int mne_read_meg_comp_eeg_ch_info_32(const QString& name,
87 QList<FIFFLIB::FiffChInfo>& megp, /* MEG channels */
88 int *nmegp,
89 QList<FIFFLIB::FiffChInfo>& meg_compp,
90 int *nmeg_compp,
91 QList<FIFFLIB::FiffChInfo>& eegp, /* EEG channels */
92 int *neegp,
93 FiffCoordTrans *meg_head_t,
94 fiffId *idp) /* The measurement ID */
95/*
96 * Read the channel information and split it into three arrays,
97 * one for MEG, one for MEG compensation channels, and one for EEG
98 */
99{
100 QFile file(name);
101 FiffStream::SPtr stream(new FiffStream(&file));
102
103 QList<FIFFLIB::FiffChInfo> chs;
104 int nchan = 0;
105 QList<FIFFLIB::FiffChInfo> meg;
106 int nmeg = 0;
107 QList<FIFFLIB::FiffChInfo> meg_comp;
108 int nmeg_comp = 0;
109 QList<FIFFLIB::FiffChInfo> eeg;
110 int neeg = 0;
111 std::unique_ptr<FiffId> id;
112 QList<FiffDirNode::SPtr> nodes;
114 FiffTag::UPtr t_pTag;
115 FIFFLIB::FiffChInfo this_ch;
117 fiff_int_t kind, pos;
118 int j,k,to_find;
119
120 if(!stream->open())
121 goto bad;
122
123 nodes = stream->dirtree()->dir_tree_find(FIFFB_MNE_PARENT_MEAS_FILE);
124
125 if (nodes.size() == 0) {
126 nodes = stream->dirtree()->dir_tree_find(FIFFB_MEAS_INFO);
127 if (nodes.size() == 0) {
128 qCritical ("Could not find the channel information.");
129 goto bad;
130 }
131 }
132 info = nodes[0];
133 to_find = 0;
134 for (k = 0; k < info->nent(); k++) {
135 kind = info->dir[k]->kind;
136 pos = info->dir[k]->pos;
137 switch (kind) {
138 case FIFF_NCHAN :
139 if (!stream->read_tag(t_pTag,pos))
140 goto bad;
141 nchan = *t_pTag->toInt();
142
143 for (j = 0; j < nchan; j++) {
144 chs.append(FiffChInfo());
145 chs[j].scanNo = -1;
146 }
147 to_find = nchan;
148 break;
149
151 if(!stream->read_tag(t_pTag, pos))
152 goto bad;
153 id = std::make_unique<FiffId>(*(FiffId*)t_pTag->data());
154 break;
155
156 case FIFF_COORD_TRANS :
157 if(!stream->read_tag(t_pTag, pos))
158 goto bad;
159// t = t_pTag->toCoordTrans();
160 t = FiffCoordTrans::readFromTag( t_pTag );
162 t = FiffCoordTrans();
163 break;
164
165 case FIFF_CH_INFO : /* Information about one channel */
166 if(!stream->read_tag(t_pTag, pos))
167 goto bad;
168
169 this_ch = t_pTag->toChInfo();
170 if (this_ch.scanNo <= 0 || this_ch.scanNo > nchan) {
171 printf ("FIFF_CH_INFO : scan # out of range %d (%d)!",this_ch.scanNo,nchan);
172 goto bad;
173 }
174 else
175 chs[this_ch.scanNo-1] = this_ch;
176 to_find--;
177 break;
178 }
179 }
180 if (to_find != 0) {
181 qCritical("Some of the channel information was missing.");
182 goto bad;
183 }
184 if (t.isEmpty() && meg_head_t != NULL) {
185 /*
186 * Try again in a more general fashion
187 */
189 if (t.isEmpty()) {
190 qCritical("MEG -> head coordinate transformation not found.");
191 goto bad;
192 }
193 }
194 /*
195 * Sort out the channels
196 */
197 for (k = 0; k < nchan; k++) {
198 if (chs[k].kind == FIFFV_MEG_CH) {
199 meg.append(chs[k]);
200 nmeg++;
201 } else if (chs[k].kind == FIFFV_REF_MEG_CH) {
202 meg_comp.append(chs[k]);
203 nmeg_comp++;
204 } else if (chs[k].kind == FIFFV_EEG_CH) {
205 eeg.append(chs[k]);
206 neeg++;
207 }
208 }
209// fiff_close(in);
210 stream->close();
211
212 megp = meg;
213 if(nmegp) {
214 *nmegp = nmeg;
215 }
216
217 meg_compp = meg_comp;
218 if(nmeg_compp) {
219 *nmeg_compp = nmeg_comp;
220 }
221
222 eegp = eeg;
223 if(neegp) {
224 *neegp = neeg;
225 }
226
227 if (idp == NULL) {
228 /* id is auto-deleted by unique_ptr */
229 }
230 else
231 *idp = id.release();
232 if (meg_head_t == NULL) {
233 }
234 else
235 *meg_head_t = t;
236
237 return FIFF_OK;
238
239bad : {
240 stream->close();
241 return FIFF_FAIL;
242 }
243}
244
245#define MNE_CTFV_COMP_UNKNOWN -1
246#define MNE_CTFV_COMP_NONE 0
247#define MNE_CTFV_COMP_G1BR 0x47314252
248#define MNE_CTFV_COMP_G2BR 0x47324252
249#define MNE_CTFV_COMP_G3BR 0x47334252
250#define MNE_CTFV_COMP_G2OI 0x47324f49
251#define MNE_CTFV_COMP_G3OI 0x47334f49
252
253static struct {
256} compMap[] = { { MNE_CTFV_NOGRAD, MNE_CTFV_COMP_NONE },
260{ MNE_4DV_COMP1, MNE_4DV_COMP1 }, /* One-to-one mapping for 4D data */
262
264
265{
266 int k;
267
268 for (k = 0; compMap[k].grad_comp >= 0; k++)
269 if (ctf_comp == compMap[k].ctf_comp)
270 return compMap[k].grad_comp;
271 return ctf_comp;
272}
273
274FiffSparseMatrix* mne_convert_to_sparse(const Eigen::MatrixXf& dense, /* The dense matrix to be converted */
275 int stor_type, /* Either FIFFTS_MC_CCS or FIFFTS_MC_RCS */
276 float small) /* How small elements should be ignored? */
277/*
278 * Create the compressed row or column storage sparse matrix representation
279 * including a vector containing the nonzero matrix element values,
280 * the row or column pointer vector and the appropriate index vector(s).
281 */
282{
283 int j,k;
284 int nz;
285 int ptr;
286 int nrow = static_cast<int>(dense.rows());
287 int ncol = static_cast<int>(dense.cols());
288 FiffSparseMatrix* sparse = NULL;
289
290 if (small < 0) { /* Automatic scaling */
291 float maxval = dense.cwiseAbs().maxCoeff();
292 if (maxval > 0)
293 small = maxval*std::fabs(small);
294 else
295 small = std::fabs(small);
296 }
297 for (j = 0, nz = 0; j < nrow; j++)
298 for (k = 0; k < ncol; k++) {
299 if (std::fabs(dense(j,k)) > small)
300 nz++;
301 }
302
303 if (nz <= 0) {
304 qWarning("No nonzero elements found.");
305 return NULL;
306 }
307 if (stor_type != FIFFTS_MC_CCS && stor_type != FIFFTS_MC_RCS) {
308 qWarning("Unknown sparse matrix storage type: %d",stor_type);
309 return NULL;
310 }
311 sparse = new FiffSparseMatrix;
312 sparse->coding = stor_type;
313 sparse->m = nrow;
314 sparse->n = ncol;
315 sparse->nz = nz;
316 sparse->data.resize(nz);
317 sparse->inds.resize(nz);
318
319 if (stor_type == FIFFTS_MC_RCS) {
320 sparse->ptrs.resize(nrow + 1);
321 for (j = 0, nz = 0; j < nrow; j++) {
322 ptr = -1;
323 for (k = 0; k < ncol; k++)
324 if (std::fabs(dense(j,k)) > small) {
325 sparse->data[nz] = dense(j,k);
326 if (ptr < 0)
327 ptr = nz;
328 sparse->inds[nz++] = k;
329 }
330 sparse->ptrs[j] = ptr;
331 }
332 sparse->ptrs[nrow] = nz;
333 for (j = nrow - 1; j >= 0; j--) /* Take care of the empty rows */
334 if (sparse->ptrs[j] < 0)
335 sparse->ptrs[j] = sparse->ptrs[j+1];
336 }
337 else if (stor_type == FIFFTS_MC_CCS) {
338 sparse->ptrs.resize(ncol + 1);
339 for (k = 0, nz = 0; k < ncol; k++) {
340 ptr = -1;
341 for (j = 0; j < nrow; j++)
342 if (std::fabs(dense(j,k)) > small) {
343 sparse->data[nz] = dense(j,k);
344 if (ptr < 0)
345 ptr = nz;
346 sparse->inds[nz++] = j;
347 }
348 sparse->ptrs[k] = ptr;
349 }
350 sparse->ptrs[ncol] = nz;
351 for (k = ncol-1; k >= 0; k--) /* Take care of the empty columns */
352 if (sparse->ptrs[k] < 0)
353 sparse->ptrs[k] = sparse->ptrs[k+1];
354 }
355 return sparse;
356}
357
358int mne_sparse_vec_mult2_32(FiffSparseMatrix* mat, /* The sparse matrix */
359 float *vector, /* Vector to be multiplied */
360 float *res) /* Result of the multiplication */
361/*
362 * Multiply a vector by a sparse matrix.
363 */
364{
365 int i,j;
366
367 if (mat->coding == FIFFTS_MC_RCS) {
368 for (i = 0; i < mat->m; i++) {
369 res[i] = 0.0;
370 for (j = mat->ptrs[i]; j < mat->ptrs[i+1]; j++)
371 res[i] += mat->data[j]*vector[mat->inds[j]];
372 }
373 return 0;
374 }
375 else if (mat->coding == FIFFTS_MC_CCS) {
376 for (i = 0; i < mat->m; i++)
377 res[i] = 0.0;
378 for (i = 0; i < mat->n; i++)
379 for (j = mat->ptrs[i]; j < mat->ptrs[i+1]; j++)
380 res[mat->inds[j]] += mat->data[j]*vector[i];
381 return 0;
382 }
383 else {
384 qWarning("mne_sparse_vec_mult2: unknown sparse matrix storage type: %d",mat->coding);
385 return -1;
386 }
387}
388
389float mne_dot_vectors_32 (float *v1,
390 float *v2,
391 int nn)
392
393{
394#ifdef BLAS
395 int one = 1;
396 float res = sdot(&nn,v1,&one,v2,&one);
397 return res;
398#else
399 float res = 0.0;
400 int k;
401
402 for (k = 0; k < nn; k++)
403 res = res + v1[k]*v2[k];
404 return res;
405#endif
406}
407
408void mne_mat_vec_mult2_32 (float **m,float *v,float *result, int d1,int d2)
409/*
410 * Matrix multiplication
411 * result(d1) = m(d1 x d2) * v(d2)
412 */
413
414{
415 int j;
416
417 for (j = 0; j < d1; j++)
418 result[j] = mne_dot_vectors_32 (m[j],v,d2);
419 return;
420}
421
422//=============================================================================================================
423// DEFINE MEMBER METHODS
424//=============================================================================================================
425
427:ncomp(0)
428,nch(0)
429,undo(nullptr)
430,current(nullptr)
431{
432}
433
434//=============================================================================================================
435
437:ncomp(0)
438,nch(set.nch)
439,undo(nullptr)
440,current(nullptr)
441{
442 if (set.ncomp > 0) {
443 for (int k = 0; k < set.ncomp; k++)
444 if(set.comps[k])
445 this->comps.append(new MNECTFCompData(*set.comps[k]));
446 this->ncomp = this->comps.size();
447 }
448
449 this->chs = set.chs;
450
451 if(set.undo)
452 this->undo = std::make_unique<MNECTFCompData>(*set.undo);
453
454 if(set.current)
455 this->current = std::make_unique<MNECTFCompData>(*set.current);
456}
457
458//=============================================================================================================
459
461{
462
463 for (int k = 0; k < comps.size(); k++)
464 if(comps[k])
465 delete comps[k];
466}
467
468//=============================================================================================================
469
470std::unique_ptr<MNECTFCompDataSet> MNECTFCompDataSet::read(const QString &name)
471/*
472 * Read all CTF compensation data from a given file
473 */
474{
475 QFile file(name);
476 FiffStream::SPtr stream(new FiffStream(&file));
477
478 std::unique_ptr<MNECTFCompDataSet> set;
479 MNECTFCompData* one;
480 QList<FiffDirNode::SPtr> nodes;
481 QList<FiffDirNode::SPtr> comps;
482 int ncomp;
483 int kind,k;
484 FiffTag::UPtr t_pTag;
485 QList<FiffChInfo> chs;
486 int nch = 0;
487 int calibrated;
488 /*
489 * Read the channel information
490 */
491 {
492 QList<FiffChInfo> comp_chs, temp;
493 int ncompch = 0;
494
496 chs,
497 &nch,
498 comp_chs,
499 &ncompch,
500 temp,
501 NULL,
502 NULL,
503 NULL) == FAIL)
504 goto bad;
505 if (ncompch > 0) {
506 for (k = 0; k < ncompch; k++)
507 chs.append(comp_chs[k]);
508 nch = nch + ncompch;
509 }
510 }
511 /*
512 * Read the rest of the stuff
513 */
514 if(!stream->open())
515 goto bad;
516 set = std::make_unique<MNECTFCompDataSet>();
517 /*
518 * Locate the compensation data sets
519 */
520 nodes = stream->dirtree()->dir_tree_find(FIFFB_MNE_CTF_COMP);
521 if (nodes.size() == 0)
522 goto good; /* Nothing more to do */
523 comps = nodes[0]->dir_tree_find(FIFFB_MNE_CTF_COMP_DATA);
524 if (comps.size() == 0)
525 goto good;
526 ncomp = comps.size();
527 /*
528 * Set the channel info
529 */
530 set->chs = chs;
531 set->nch = nch;
532 /*
533 * Read each data set
534 */
535 for (k = 0; k < ncomp; k++) {
537 if (!mat)
538 goto bad;
539 comps[k]->find_tag(stream, FIFF_MNE_CTF_COMP_KIND, t_pTag);
540 if (t_pTag) {
541 kind = *t_pTag->toInt();
542 }
543 else
544 goto bad;
545 comps[k]->find_tag(stream, FIFF_MNE_CTF_COMP_CALIBRATED, t_pTag);
546 if (t_pTag) {
547 calibrated = *t_pTag->toInt();
548 }
549 else
550 calibrated = FALSE;
551 /*
552 * Add these data to the set
553 */
554 one = new MNECTFCompData;
555 one->data = std::move(mat);
556 one->kind = kind;
558 one->calibrated = calibrated;
559
560 if (one->calibrate(set->chs,set->nch,TRUE) == FAIL) {
561 printf("Warning: Compensation data for '%s' omitted\n", explain_comp(one->kind).toUtf8().constData());//,err_get_error(),explain_comp(one->kind));
562 if(one)
563 delete one;
564 }
565 else {
566 set->comps.append(one);
567 set->ncomp++;
568 }
569 }
570#ifdef DEBUG
571 printf("%d CTF compensation data sets read from %s\n",set->ncomp,name);
572#endif
573 goto good;
574
575bad : {
576 stream->close();
577 return nullptr;
578 }
579
580good : {
581 stream->close();
582 return set;
583 }
584}
585
586//=============================================================================================================
587
588int MNECTFCompDataSet::make_comp(const QList<FiffChInfo>& chs,
589 int nch,
590 QList<FiffChInfo> compchs,
591 int ncomp) /* How many of these */
592/*
593 * Make compensation data to apply to a set of channels to yield (or uncompensated) compensated data
594 */
595{
596 Eigen::VectorXi comps;
597 int need_comp;
598 int first_comp;
599 MNECTFCompData* this_comp;
600 Eigen::VectorXi comp_sel;
601 QStringList names;
602 QString name;
603 int j,k,p;
604
605 std::unique_ptr<FiffSparseMatrix> presel;
606 std::unique_ptr<FiffSparseMatrix> postsel;
607 std::unique_ptr<MNENamedMatrix> data;
608
609 QStringList emptyList;
610
611 if (compchs.isEmpty()) {
612 compchs = chs;
613 ncomp = nch;
614 }
615 printf("Setting up compensation data...\n");
616 if (nch == 0)
617 return OK;
618 current.reset();
619 comps.resize(nch);
620 for (k = 0, need_comp = 0, first_comp = MNE_CTFV_COMP_NONE; k < nch; k++) {
621 if (chs[k].kind == FIFFV_MEG_CH) {
622 comps[k] = chs[k].chpos.coil_type >> 16;
623 if (comps[k] != MNE_CTFV_COMP_NONE) {
624 if (first_comp == MNE_CTFV_COMP_NONE)
625 first_comp = comps[k];
626 else {
627 if (comps[k] != first_comp) {
628 printf("We do not support nonuniform compensation yet.");
629 return FAIL;
630 }
631 }
632 need_comp++;
633 }
634 }
635 else
637 }
638 if (need_comp == 0) {
639 printf("\tNo compensation set. Nothing more to do.\n");
640 return OK;
641 }
642 printf("\t%d out of %d channels have the compensation set.\n",need_comp,nch);
643 /*
644 * Find the desired compensation data matrix
645 */
646 for (k = 0, this_comp = NULL; k < this->ncomp; k++) {
647 if (this->comps[k]->mne_kind == first_comp) {
648 this_comp = this->comps[k];
649 break;
650 }
651 }
652 if (!this_comp) {
653 printf("Did not find the desired compensation data : %s",
654 explain_comp(map_comp_kind(first_comp)).toUtf8().constData());
655 return FAIL;
656 }
657 printf("\tDesired compensation data (%s) found.\n",explain_comp(map_comp_kind(first_comp)).toUtf8().constData());
658 /*
659 * Find the compensation channels
660 */
661 comp_sel.resize(this_comp->data->ncol);
662 for (k = 0; k < this_comp->data->ncol; k++) {
663 comp_sel[k] = -1;
664 name = this_comp->data->collist[k];
665 for (p = 0; p < ncomp; p++)
666 if (QString::compare(name,compchs[p].ch_name) == 0) {
667 comp_sel[k] = p;
668 break;
669 }
670 if (comp_sel[k] < 0) {
671 printf("Compensation channel %s not found",name.toUtf8().constData());
672 return FAIL;
673 }
674 }
675 printf("\tAll compensation channels found.\n");
676 /*
677 * Create the preselector
678 */
679 {
680 Eigen::MatrixXf sel = Eigen::MatrixXf::Zero(this_comp->data->ncol, ncomp);
681 for (j = 0; j < this_comp->data->ncol; j++)
682 sel(j, comp_sel[j]) = 1.0f;
684 if (!ps)
685 return FAIL;
686 presel.reset(ps);
687 printf("\tPreselector created.\n");
688 }
689 /*
690 * Pick the desired channels
691 */
692 for (k = 0; k < nch; k++) {
693 if (comps[k] != MNE_CTFV_COMP_NONE)
694 names.append(chs[k].ch_name);
695 }
696
697 {
698 auto d = this_comp->data->pick(names, need_comp, emptyList, 0);
699 if (!d)
700 return FAIL;
701 data = std::move(d);
702 }
703 printf("\tCompensation data matrix created.\n");
704 /*
705 * Create the postselector
706 */
707 {
708 Eigen::MatrixXf sel = Eigen::MatrixXf::Zero(nch, data->nrow);
709 for (j = 0, p = 0; j < nch; j++) {
710 if (comps[j] != MNE_CTFV_COMP_NONE)
711 sel(j, p++) = 1.0f;
712 }
714 if (!ps)
715 return FAIL;
716 postsel.reset(ps);
717 printf("\tPostselector created.\n");
718 }
719 current = std::make_unique<MNECTFCompData>();
720 current->kind = this_comp->kind;
721 current->mne_kind = this_comp->mne_kind;
722 current->data = std::move(data);
723 current->presel = std::move(presel);
724 current->postsel = std::move(postsel);
725
726 printf("\tCompensation set up.\n");
727 return OK;
728}
729
730//=============================================================================================================
731
732int MNECTFCompDataSet::set_comp(QList<FIFFLIB::FiffChInfo>& chs,
733 int nch,
734 int comp)
735/*
736 * Set the compensation bits to the desired value
737 */
738{
739 int k;
740 int nset;
741 for (k = 0, nset = 0; k < nch; k++) {
742 if (chs[k].kind == FIFFV_MEG_CH) {
743 chs[k].chpos.coil_type = (chs[k].chpos.coil_type & 0xFFFF) | (comp << 16);
744 nset++;
745 }
746 }
747 printf("A new compensation value (%s) was assigned to %d MEG channels.\n",
748 explain_comp(map_comp_kind(comp)).toUtf8().constData(),nset);
749 return nset;
750}
751
752//=============================================================================================================
753
754int MNECTFCompDataSet::apply(int do_it, Eigen::Ref<Eigen::VectorXf> data)
755{
756 return apply(do_it, data, data);
757}
758
759//=============================================================================================================
760
761int MNECTFCompDataSet::apply(int do_it, Eigen::Ref<Eigen::VectorXf> data, Eigen::Ref<const Eigen::VectorXf> compdata)
762/*
763 * Apply compensation or revert to uncompensated data
764 */
765{
766 MNECTFCompData* this_comp;
767 int ndata = static_cast<int>(data.size());
768 int ncompdata = static_cast<int>(compdata.size());
769
770 if (!current)
771 return OK;
772 this_comp = current.get();
773 /*
774 * Dimension checks
775 */
776 if (this_comp->presel) {
777 if (this_comp->presel->n != ncompdata) {
778 printf("Compensation data dimension mismatch. Expected %d, got %d channels.",
779 this_comp->presel->n,ncompdata);
780 return FAIL;
781 }
782 }
783 else if (this_comp->data->ncol != ncompdata) {
784 printf("Compensation data dimension mismatch. Expected %d, got %d channels.",
785 this_comp->data->ncol,ncompdata);
786 return FAIL;
787 }
788 if (this_comp->postsel) {
789 if (this_comp->postsel->m != ndata) {
790 printf("Data dimension mismatch. Expected %d, got %d channels.",
791 this_comp->postsel->m,ndata);
792 return FAIL;
793 }
794 }
795 else if (this_comp->data->nrow != ndata) {
796 printf("Data dimension mismatch. Expected %d, got %d channels.",
797 this_comp->data->nrow,ndata);
798 return FAIL;
799 }
800 /*
801 * Preselection is optional
802 */
803 const float *presel;
804 if (this_comp->presel) {
805 if (this_comp->presel_data.size() == 0)
806 this_comp->presel_data.resize(this_comp->presel->m);
807 if (mne_sparse_vec_mult2_32(this_comp->presel.get(),const_cast<float*>(compdata.data()),this_comp->presel_data.data()) != OK)
808 return FAIL;
809 presel = this_comp->presel_data.data();
810 }
811 else
812 presel = compdata.data();
813 /*
814 * This always happens
815 */
816 if (this_comp->comp_data.size() == 0)
817 this_comp->comp_data.resize(this_comp->data->nrow);
818 {
819 Eigen::Map<const Eigen::VectorXf> preselVec(presel, this_comp->data->ncol);
820 Eigen::Map<Eigen::VectorXf> compVec(this_comp->comp_data.data(), this_comp->data->nrow);
821 compVec = this_comp->data->data * preselVec;
822 }
823 /*
824 * Optional postselection
825 */
826 const float *comp;
827 if (!this_comp->postsel)
828 comp = this_comp->comp_data.data();
829 else {
830 if (this_comp->postsel_data.size() == 0)
831 this_comp->postsel_data.resize(this_comp->postsel->m);
832 if (mne_sparse_vec_mult2_32(this_comp->postsel.get(),this_comp->comp_data.data(),this_comp->postsel_data.data()) != OK)
833 return FAIL;
834 comp = this_comp->postsel_data.data();
835 }
836 /*
837 * Compensate or revert compensation?
838 */
839 Eigen::Map<const Eigen::VectorXf> compVec(comp, ndata);
840 if (do_it)
841 data -= compVec;
842 else
843 data += compVec;
844 return OK;
845}
846
847//=============================================================================================================
848
849int MNECTFCompDataSet::apply_transpose(int do_it, Eigen::MatrixXf& data)
850/*
851 * Apply compensation or revert to uncompensated data
852 */
853{
854 using RowMatrixXf = Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
855
856 MNECTFCompData* this_comp;
857 int ndata = static_cast<int>(data.rows());
858 int ns = static_cast<int>(data.cols());
859 int ncompdata = ndata;
860
861 if (!current)
862 return OK;
863 this_comp = current.get();
864 /*
865 * Dimension checks
866 */
867 if (this_comp->presel) {
868 if (this_comp->presel->n != ncompdata) {
869 printf("Compensation data dimension mismatch. Expected %d, got %d channels.",
870 this_comp->presel->n,ncompdata);
871 return FAIL;
872 }
873 }
874 else if (this_comp->data->ncol != ncompdata) {
875 printf("Compensation data dimension mismatch. Expected %d, got %d channels.",
876 this_comp->data->ncol,ncompdata);
877 return FAIL;
878 }
879 if (this_comp->postsel) {
880 if (this_comp->postsel->m != ndata) {
881 printf("Data dimension mismatch. Expected %d, got %d channels.",
882 this_comp->postsel->m,ndata);
883 return FAIL;
884 }
885 }
886 else if (this_comp->data->nrow != ndata) {
887 printf("Data dimension mismatch. Expected %d, got %d channels.",
888 this_comp->data->nrow,ndata);
889 return FAIL;
890 }
891 /*
892 * Preselection is optional — sparse matrix * data
893 */
894 Eigen::MatrixXf preselMat;
895 if (this_comp->presel) {
896 preselMat.resize(this_comp->presel->m, ns);
897 FiffSparseMatrix* sp = this_comp->presel.get();
898 if (sp->coding == FIFFTS_MC_RCS) {
899 for (int i = 0; i < sp->m; i++)
900 for (int c = 0; c < ns; c++) {
901 float val = 0.0f;
902 for (int j = sp->ptrs[i]; j < sp->ptrs[i+1]; j++)
903 val += sp->data[j] * data(sp->inds[j], c);
904 preselMat(i, c) = val;
905 }
906 } else if (sp->coding == FIFFTS_MC_CCS) {
907 preselMat.setZero();
908 for (int c = 0; c < ns; c++)
909 for (int i = 0; i < sp->n; i++)
910 for (int j = sp->ptrs[i]; j < sp->ptrs[i+1]; j++)
911 preselMat(sp->inds[j], c) += sp->data[j] * data(i, c);
912 } else {
913 printf("Unknown sparse matrix storage type: %d", sp->coding);
914 return FAIL;
915 }
916 }
917 else {
918 preselMat = data;
919 }
920 /*
921 * Compensation: comp = data_matrix * preselMat
922 * data_matrix is float** (contiguous row-major via mne_cmatrix)
923 */
924 Eigen::MatrixXf comp = this_comp->data->data * preselMat;
925 /*
926 * Optional postselection — sparse matrix * comp
927 */
928 if (this_comp->postsel) {
929 Eigen::MatrixXf postselMat(this_comp->postsel->m, ns);
930 FiffSparseMatrix* sp = this_comp->postsel.get();
931 if (sp->coding == FIFFTS_MC_RCS) {
932 for (int i = 0; i < sp->m; i++)
933 for (int c = 0; c < ns; c++) {
934 float val = 0.0f;
935 for (int j = sp->ptrs[i]; j < sp->ptrs[i+1]; j++)
936 val += sp->data[j] * comp(sp->inds[j], c);
937 postselMat(i, c) = val;
938 }
939 } else if (sp->coding == FIFFTS_MC_CCS) {
940 postselMat.setZero();
941 for (int c = 0; c < ns; c++)
942 for (int i = 0; i < sp->n; i++)
943 for (int j = sp->ptrs[i]; j < sp->ptrs[i+1]; j++)
944 postselMat(sp->inds[j], c) += sp->data[j] * comp(i, c);
945 } else {
946 printf("Unknown sparse matrix storage type: %d", sp->coding);
947 return FAIL;
948 }
949 comp = std::move(postselMat);
950 }
951 /*
952 * Compensate or revert compensation?
953 */
954 if (do_it)
955 data -= comp;
956 else
957 data += comp;
958 return OK;
959}
960
961//=============================================================================================================
962
963int MNECTFCompDataSet::get_comp(const QList<FIFFLIB::FiffChInfo> &chs, int nch)
964{
965 int res = MNE_CTFV_NOGRAD;
966 int first_comp,comp;
967 int k;
968
969 for (k = 0, first_comp = -1; k < nch; k++) {
970 if (chs[k].kind == FIFFV_MEG_CH) {
971 comp = chs[k].chpos.coil_type >> 16;
972 if (first_comp < 0)
973 first_comp = comp;
974 else if (first_comp != comp) {
975 printf("Non uniform compensation not supported.");
976 return FAIL;
977 }
978 }
979 }
980 if (first_comp >= 0)
981 res = first_comp;
982 return res;
983}
984
985//=============================================================================================================
986
988/*
989 * Simple mapping
990 */
991{
992 int k;
993
994 for (k = 0; compMap[k].grad_comp >= 0; k++)
995 if (grad == compMap[k].grad_comp)
996 return compMap[k].ctf_comp;
997 return grad;
998}
999
1000//=============================================================================================================
1001
1003{
1004 static const struct {
1005 int kind;
1006 const char *expl;
1007 } explain[] = { { MNE_CTFV_COMP_NONE, "uncompensated" },
1008 { MNE_CTFV_COMP_G1BR, "first order gradiometer" },
1009 { MNE_CTFV_COMP_G2BR, "second order gradiometer" },
1010 { MNE_CTFV_COMP_G3BR, "third order gradiometer" },
1011 { MNE_4DV_COMP1, "4D comp 1" },
1012 { MNE_CTFV_COMP_UNKNOWN, "unknown" } };
1013 int k;
1014
1015 for (k = 0; explain[k].kind != MNE_CTFV_COMP_UNKNOWN; k++)
1016 if (explain[k].kind == kind)
1017 return QString::fromLatin1(explain[k].expl);
1018 return QString::fromLatin1(explain[k].expl);
1019}
1020
1021//=============================================================================================================
1022
1024 QList<FiffChInfo>& chs,
1025 int nchan,
1026 QList<FiffChInfo> comp_chs,
1027 int ncomp_chan) /* How many */
1028/*
1029 * Make data which has the third-order gradient compensation applied
1030 */
1031{
1032 int k;
1033 int have_comp_chs;
1034 int comp_was = MNE_CTFV_COMP_UNKNOWN;
1035
1036 if (comp_chs.isEmpty()) {
1037 comp_chs = chs;
1038 ncomp_chan = nchan;
1039 }
1040 undo.reset();
1041 current.reset();
1042 for (k = 0, have_comp_chs = 0; k < ncomp_chan; k++)
1043 if (comp_chs[k].kind == FIFFV_REF_MEG_CH)
1044 have_comp_chs++;
1045 if (have_comp_chs == 0 && compensate_to != MNE_CTFV_NOGRAD) {
1046 printf("No compensation channels in these data.");
1047 return FAIL;
1048 }
1049 /*
1050 * Update the 'current' field to reflect the compensation possibly present in the data now
1051 */
1052 if (make_comp(chs,nchan,comp_chs,ncomp_chan) == FAIL)
1053 goto bad;
1054 /*
1055 * Are we there already?
1056 */
1057 if (current && current->mne_kind == compensate_to) {
1058 printf("No further compensation necessary (comp = %s)\n",explain_comp(current->kind).toUtf8().constData());
1059 current.reset();
1060 return OK;
1061 }
1062 undo = std::move(current);
1063 if (compensate_to == MNE_CTFV_NOGRAD) {
1064 printf("No compensation was requested.\n");
1065 set_comp(chs,nchan,compensate_to);
1066 return OK;
1067 }
1068 if (set_comp(chs,nchan,compensate_to) > 0) {
1069 if (undo)
1070 comp_was = undo->mne_kind;
1071 else
1072 comp_was = MNE_CTFV_NOGRAD;
1073 if (make_comp(chs,nchan,comp_chs,ncomp_chan) == FAIL)
1074 goto bad;
1075 printf("Compensation set up as requested (%s -> %s).\n",
1076 explain_comp(map_comp_kind(comp_was)).toUtf8().constData(),
1077 explain_comp(current->kind).toUtf8().constData());
1078 }
1079 return OK;
1080
1081bad : {
1082 if (comp_was != MNE_CTFV_COMP_UNKNOWN)
1083 set_comp(chs,nchan,comp_was);
1084 return FAIL;
1085 }
1086}
#define TRUE
#define FALSE
#define OK
#define FAIL
#define FIFFV_EEG_CH
#define FIFF_OK
#define FIFF_MNE_CTF_COMP_KIND
#define FIFFV_COORD_DEVICE
#define FIFFB_MNE_CTF_COMP_DATA
#define FIFF_MNE_CTF_COMP_CALIBRATED
#define FIFF_FAIL
#define FIFFV_REF_MEG_CH
#define FIFFV_MEG_CH
#define FIFF_MNE_CTF_COMP_DATA
#define FIFFV_COORD_HEAD
#define FIFFB_MNE_CTF_COMP
#define FIFFB_MNE_PARENT_MEAS_FILE
#define FIFF_PARENT_BLOCK_ID
Definition fiff_file.h:333
#define FIFF_NCHAN
Definition fiff_file.h:453
#define FIFFTS_MC_RCS
Definition fiff_file.h:270
#define FIFFTS_MC_CCS
Definition fiff_file.h:269
#define FIFF_COORD_TRANS
Definition fiff_file.h:475
#define FIFF_CH_INFO
Definition fiff_file.h:456
#define FIFFB_MEAS_INFO
Definition fiff_file.h:363
Old fiff_type declarations - replace them.
FiffCoordTrans class declaration.
#define MNE_CTFV_COMP_G2BR
#define MNE_CTFV_COMP_NONE
#define MNE_CTFV_COMP_UNKNOWN
#define MNE_CTFV_COMP_G3BR
#define MNE_CTFV_COMP_G1BR
MNECTFCompDataSet class declaration.
Legacy MNE-C constants and common typedefs.
#define MNE_CTFV_GRAD1
Definition mne_types.h:132
#define MNE_CTFV_GRAD2
Definition mne_types.h:133
#define MNE_CTFV_GRAD3
Definition mne_types.h:134
#define MNE_4DV_COMP1
Definition mne_types.h:141
#define MNE_CTFV_NOGRAD
Definition mne_types.h:131
MNECTFCompData class declaration.
float mne_dot_vectors_32(float *v1, float *v2, int nn)
void mne_mat_vec_mult2_32(float **m, float *v, float *result, int d1, int d2)
FiffSparseMatrix * mne_convert_to_sparse(const Eigen::MatrixXf &dense, int stor_type, float small)
int mne_unmap_ctf_comp_kind(int ctf_comp)
int mne_sparse_vec_mult2_32(FiffSparseMatrix *mat, float *vector, float *res)
int mne_read_meg_comp_eeg_ch_info_32(const QString &name, QList< FIFFLIB::FiffChInfo > &megp, int *nmegp, QList< FIFFLIB::FiffChInfo > &meg_compp, int *nmeg_compp, QList< FIFFLIB::FiffChInfo > &eegp, int *neegp, FiffCoordTrans *meg_head_t, fiffId *idp)
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
FiffId * fiffId
Backward-compatible pointer typedef for the old fiffId pointer.
Definition fiff_types.h:130
Channel info descriptor.
Coordinate transformation description.
QSharedPointer< FiffDirNode > SPtr
Universally unique identifier.
Definition fiff_id.h:68
FIFF sparse matrix storage.
FIFF File I/O routines.
QSharedPointer< FiffStream > SPtr
std::unique_ptr< FiffTag > UPtr
Definition fiff_tag.h:158
Represents a single CTF compensation data element.
std::unique_ptr< FIFFLIB::FiffSparseMatrix > postsel
std::unique_ptr< MNENamedMatrix > data
int calibrate(const QList< FIFFLIB::FiffChInfo > &chs, int nch, int do_it)
std::unique_ptr< FIFFLIB::FiffSparseMatrix > presel
std::unique_ptr< MNECTFCompData > undo
int make_comp(const QList< FIFFLIB::FiffChInfo > &chs, int nch, QList< FIFFLIB::FiffChInfo > compchs, int ncomp)
std::unique_ptr< MNECTFCompData > current
QList< FIFFLIB::FiffChInfo > chs
static int set_comp(QList< FIFFLIB::FiffChInfo > &chs, int nch, int comp)
static std::unique_ptr< MNECTFCompDataSet > read(const QString &name)
int apply_transpose(int do_it, Eigen::MatrixXf &data)
static QString explain_comp(int kind)
QList< MNECTFCompData * > comps
static int get_comp(const QList< FIFFLIB::FiffChInfo > &chs, int nch)
int set_compensation(int compensate_to, QList< FIFFLIB::FiffChInfo > &chs, int nchan, QList< FIFFLIB::FiffChInfo > comp_chs, int ncomp_chan)
int apply(int do_it, Eigen::Ref< Eigen::VectorXf > data, Eigen::Ref< const Eigen::VectorXf > compdata)
static std::unique_ptr< MNENamedMatrix > read(QSharedPointer< FIFFLIB::FiffStream > &stream, const QSharedPointer< FIFFLIB::FiffDirNode > &node, int kind)
Factory: read a named matrix from a FIFF file.
static FiffCoordTrans readMeasTransform(const QString &name)
static FiffCoordTrans readFromTag(const std::unique_ptr< FiffTag > &tag)