v2.0.0
Loading...
Searching...
No Matches
mne_raw_data.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "mne_raw_data.h"
42
43#include <QFile>
44#include <QDebug>
45#include <QTextStream>
46
47#include <Eigen/Core>
48
49#define _USE_MATH_DEFINES
50#include <math.h>
51
52//=============================================================================================================
53// USED NAMESPACES
54//=============================================================================================================
55
56using namespace Eigen;
57using namespace FIFFLIB;
58using namespace MNELIB;
59
60constexpr int FAIL = -1;
61constexpr int OK = 0;
62
63#if defined(_WIN32) || defined(_WIN64)
64#define snprintf _snprintf
65#define vsnprintf _vsnprintf
66#define strcasecmp _stricmp
67#define strncasecmp _strnicmp
68#endif
69
70namespace MNELIB
71{
72
74using RowMajorMatrixXf = Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
75
84{
85 struct Entry {
87 };
88
89 std::vector<Entry> entries;
90 int next = 0;
91
92 explicit RingBuffer(int nslots)
93 : entries(static_cast<size_t>(nslots))
94 , next(0)
95 {}
96
98 void allocate(int nrow, int ncol, RowMajorMatrixXf *res)
99 {
100 if (next >= static_cast<int>(entries.size()))
101 next = 0;
102 Entry &e = entries[static_cast<size_t>(next++)];
103 if (e.user) // evict old occupant
104 e.user->resize(0, 0);
105 res->resize(nrow, ncol);
106 e.user = res;
107 }
108};
109
110}
111
112//============================= misc_util.c =============================
113
114//============================= mne_apply_filter.c =============================
115
116namespace MNELIB
117{
118
126{
127 std::vector<float> freq_resp;
128 std::vector<float> eog_freq_resp;
129 std::vector<float> precalc;
130
131 explicit FilterData(int resp_size)
132 : freq_resp(static_cast<size_t>(resp_size), 1.0f)
133 , eog_freq_resp(static_cast<size_t>(resp_size), 1.0f)
134 {}
135};
136
137}
138
140 const MNEFilterDef& f2)
141/*
142 * Return 0 if the two filter definitions are same, 1 otherwise
143 */
144{
145 if (f1.filter_on != f2.filter_on ||
146 std::fabs(f1.lowpass - f2.lowpass) > 0.1 ||
147 std::fabs(f1.lowpass_width - f2.lowpass_width) > 0.1 ||
148 std::fabs(f1.highpass - f2.highpass) > 0.1 ||
149 std::fabs(f1.highpass_width - f2.highpass_width) > 0.1 ||
150 std::fabs(f1.eog_lowpass - f2.eog_lowpass) > 0.1 ||
151 std::fabs(f1.eog_lowpass_width - f2.eog_lowpass_width) > 0.1 ||
152 std::fabs(f1.eog_highpass - f2.eog_highpass) > 0.1 ||
153 std::fabs(f1.eog_highpass_width - f2.eog_highpass_width) > 0.1)
154 return 1;
155 else
156 return 0;
157}
158
159//============================= mne_fft.c =============================
160
161void mne_fft_ana(float *data, int np, std::vector<float>& /*precalc*/)
162/*
163 * FFT analysis for real data
164 */
165{
166 Q_UNUSED(data);
167 Q_UNUSED(np);
168 qCritical("##################### DEBUG Error: FFT analysis needs to be implemented");
169 return;
170}
171
172void mne_fft_syn(float *data, int np, std::vector<float>& /*precalc*/)
173/*
174 * FFT synthesis for real data
175 */
176{
177 Q_UNUSED(data);
178 Q_UNUSED(np);
179 qCritical("##################### DEBUG Error: FFT synthesis needs to be implemented");
180 return;
181}
182
183int mne_apply_filter(const MNEFilterDef& filter, FilterData *d, float *data, int ns, int zero_pad, float dc_offset, int kind)
184/*
185 * Do the magick trick
186 */
187{
188 int k,p,n;
189 float *freq_resp;
190
191 if (ns != filter.size + 2*filter.taper_size) {
192 qCritical("Incorrect data length in apply_filter");
193 return FAIL;
194 }
195 /*
196 * Zero padding
197 */
198 if (zero_pad) {
199 for (k = 0; k < filter.taper_size; k++)
200 data[k] = 0.0;
201 for (k = filter.taper_size + filter.size; k < ns; k++)
202 data[k] = 0.0;
203 }
204 if (!filter.filter_on) /* Nothing else to do */
205 return OK;
206 /*
207 * Make things nice by compensating for the dc offset
208 */
209 if (dc_offset != 0.0) {
210 for (k = filter.taper_size; k < filter.taper_size + filter.size; k++)
211 data[k] = data[k] - dc_offset;
212 }
213 if (!d)
214 return OK;
215 if (d->freq_resp.empty())
216 return OK;
217 /*
218 * Next comes the FFT
219 */
220 mne_fft_ana(data,ns,d->precalc);
221 /*
222 * Multiply with the frequency response
223 * See FFTpack doc for details of the arrangement
224 */
225 n = ns % 2 == 0 ? ns/2 : (ns+1)/2;
226 p = 0;
227 /*
228 * No imaginary part for the DC component
229 */
230 if (kind == FIFFV_EOG_CH)
231 freq_resp = d->eog_freq_resp.data();
232 else
233 freq_resp = d->freq_resp.data();
234 data[p] = data[p]*freq_resp[0]; p++;
235 /*
236 * The other components
237 */
238 for (k = 1 ; k < n ; k++) {
239 data[p] = data[p]*freq_resp[k]; p++;
240 data[p] = data[p]*freq_resp[k]; p++;
241 }
242 /*
243 * Then the last value
244 */
245 if (ns % 2 == 0)
246 data[p] = data[p]*freq_resp[k];
247
248 mne_fft_syn(data,ns,d->precalc);
249
250 return OK;
251}
252
253std::unique_ptr<FilterData> mne_create_filter_response(const MNEFilterDef& filter,
254 float sfreq,
255 int *highpass_effective)
256/*
257 * Create a frequency response
258 */
259{
260 int resp_size;
261 int k,s,w,f;
262 int highpasss,lowpasss;
263 int highpass_widths,lowpass_widths;
264 float lowpass,highpass,lowpass_width,highpass_width;
265 float *freq_resp;
266 float pi4 = M_PI/4.0;
267 float mult,add,c;
268
269 resp_size = (filter.size + 2*filter.taper_size)/2 + 1;
270
271 auto filter_data = std::make_unique<FilterData>(resp_size);
272 *highpass_effective = false;
273
274 for (f = 0; f < 2; f++) {
275 highpass = f == 0 ? filter.highpass : filter.eog_highpass;
276 highpass_width = f == 0 ? filter.highpass_width : filter.eog_highpass_width;
277 lowpass = f == 0 ? filter.lowpass : filter.eog_lowpass;
278 lowpass_width = f == 0 ? filter.lowpass_width : filter.eog_lowpass_width;
279 freq_resp = f == 0 ? filter_data->freq_resp.data() : filter_data->eog_freq_resp.data();
280 /*
281 * Start simple first
282 */
283 highpasss = ((resp_size-1)*highpass)/(0.5*sfreq);
284 lowpasss = ((resp_size-1)*lowpass)/(0.5*sfreq);
285
286 lowpass_widths = ((resp_size-1)*lowpass_width)/(0.5*sfreq);
287 lowpass_widths = (lowpass_widths+1)/2; /* What user specified */
288
289 if (filter.highpass_width > 0.0) {
290 highpass_widths = ((resp_size-1)*highpass_width)/(0.5*sfreq);
291 highpass_widths = (highpass_widths+1)/2; /* What user specified */
292 }
293 else
294 highpass_widths = 3; /* Minimal */
295
296 if (filter.filter_on) {
297 qInfo("filter : %7.3f ... %6.1f Hz bins : %d ... %d of %d hpw : %d lpw : %d\n",
298 highpass,
299 lowpass,
300 highpasss,
301 lowpasss,
302 resp_size,
303 highpass_widths,
304 lowpass_widths);
305 }
306 if (highpasss > highpass_widths + 1) {
307 w = highpass_widths;
308 mult = 1.0/w;
309 add = 3.0;
310 for (k = 0; k < highpasss-w+1; k++)
311 freq_resp[k] = 0.0;
312 for (k = -w+1, s = highpasss-w+1; k < w; k++, s++) {
313 if (s >= 0 && s < resp_size) {
314 c = cos(pi4*(k*mult+add));
315 freq_resp[s] = freq_resp[s]*c*c;
316 }
317 }
318 *highpass_effective = true;
319 }
320 else
321 *highpass_effective = *highpass_effective || (filter.highpass == 0.0);
322
323 if (lowpass_widths > 0) {
324 w = lowpass_widths;
325 mult = 1.0/w;
326 add = 1.0;
327 for (k = -w+1, s = lowpasss-w+1; k < w; k++, s++) {
328 if (s >= 0 && s < resp_size) {
329 c = cos(pi4*(k*mult+add));
330 freq_resp[s] = freq_resp[s]*c*c;
331 }
332 }
333 for (k = s; k < resp_size; k++)
334 freq_resp[k] = 0.0;
335 }
336 else {
337 for (k = lowpasss; k < resp_size; k++)
338 freq_resp[k] = 0.0;
339 }
340 if (filter.filter_on) {
341 if (*highpass_effective)
342 qInfo("Highpass filter will work as specified.\n");
343 else
344 qWarning("NOTE: Highpass filter omitted due to a too low corner frequency.\n");
345 }
346 else
347 qWarning("NOTE: Filter is presently switched off.\n");
348 }
349 return filter_data;
350}
351
352//============================= mne_raw_routines.c =============================
353
354int mne_read_raw_buffer_t(//fiffFile in, /* Input file */
355 FiffStream::SPtr& stream,
356 const FiffDirEntry::SPtr& ent, /* The directory entry to read */
357 RowMajorMatrixXf& data, /* Matrix [npick x nsamp] to fill */
358 int nchan, /* Number of channels in the data */
359 int nsamp, /* Expected number of samples */
360 const QList<FIFFLIB::FiffChInfo>& chs, /* Channel info for ALL channels */
361 int *pickno, /* Which channels to pick */
362 int npick) /* How many */
363
364{
365 FiffTag::UPtr t_pTag;
366// fiffTagRec tag;
367 fiff_short_t *this_samples;
368 const fiff_float_t *this_samplef;
369 fiff_int_t *this_sample;
370
371 int s,c;
372 int do_all;
373
374// tag.data = NULL;
375
376 Eigen::VectorXi pickno_vec;
377 if (npick == 0) {
378 pickno_vec = Eigen::VectorXi::LinSpaced(nchan, 0, nchan - 1);
379 pickno = pickno_vec.data();
380 do_all = true;
381 npick = nchan;
382 }
383 else
384 do_all = false;
385
386 Eigen::VectorXf mult(npick);
387 for (c = 0; c < npick; c++)
388 mult[c] = chs[pickno[c]].cal*chs[pickno[c]].range;
389
390// if (fiff_read_this_tag(in->fd,ent->pos,&tag) == FIFF_FAIL)
391// goto bad;
392 if (!stream->read_tag(t_pTag,ent->pos))
393 return FAIL;
394
395 if (ent->type == FIFFT_FLOAT) {
396 if (static_cast<int>(t_pTag->size()/(sizeof(fiff_float_t)*nchan)) != nsamp) {
397 qCritical("Incorrect number of samples in buffer.");
398 return FAIL;
399 }
400 qDebug() << "ToDo: Check whether this_samplef contains the right stuff!!! - use VectorXf instead";
401 this_samplef = t_pTag->toFloat();
402 for (s = 0; s < nsamp; s++, this_samplef += nchan) {
403 for (c = 0; c < npick; c++)
404 data(c,s) = mult[c]*this_samplef[pickno[c]];
405 }
406 }
407 else if (ent->type == FIFFT_SHORT || ent->type == FIFFT_DAU_PACK16) {
408 if (static_cast<int>(t_pTag->size()/(sizeof(fiff_short_t)*nchan)) != nsamp) {
409 qCritical("Incorrect number of samples in buffer.");
410 return FAIL;
411 }
412 qDebug() << "ToDo: Check whether this_samples contains the right stuff!!! - use VectorXi instead";
413 this_samples = (fiff_short_t *)t_pTag->data();
414 for (s = 0; s < nsamp; s++, this_samples += nchan) {
415 for (c = 0; c < npick; c++)
416 data(c,s) = mult[c]*this_samples[pickno[c]];
417 }
418 }
419 else if (ent->type == FIFFT_INT) {
420 if (static_cast<int>(t_pTag->size()/(sizeof(fiff_int_t)*nchan)) != nsamp) {
421 qCritical("Incorrect number of samples in buffer.");
422 return FAIL;
423 }
424 qDebug() << "ToDo: Check whether this_sample contains the right stuff!!! - use VectorXi instead";
425 this_sample = t_pTag->toInt();
426 for (s = 0; s < nsamp; s++, this_sample += nchan) {
427 for (c = 0; c < npick; c++)
428 data(c,s) = mult[c]*this_sample[pickno[c]];
429 }
430 }
431 else {
432 qCritical("We are not prepared to handle raw data type: %d",ent->type);
433 return FAIL;
434 }
435 return OK;
436}
437
438//============================= mne_process_bads.c =============================
439
441 const FiffDirNode::SPtr& pNode, QStringList& listp, int& nlistp)
442{
443 FiffDirNode::SPtr node,bad;
444 QList<FiffDirNode::SPtr> temp;
445 QStringList list;
446 int nlist = 0;
447 FiffTag::UPtr t_pTag;
448 QString names;
449
450 if (pNode->isEmpty())
451 node = stream->dirtree();
452 else
453 node = pNode;
454
455 temp = node->dir_tree_find(FIFFB_MNE_BAD_CHANNELS);
456 if (temp.size() > 0) {
457 bad = temp[0];
458
459 bad->find_tag(stream, FIFF_MNE_CH_NAME_LIST, t_pTag);
460 if (t_pTag) {
461 names = t_pTag->toString();
462 list = FiffStream::split_name_list(names);
463 nlist = list.size();
464 }
465 }
466 listp = list;
467 nlistp = nlist;
468 return OK;
469}
470
471int mne_read_bad_channel_list(const QString& name, QStringList& listp, int& nlistp)
472
473{
474 QFile file(name);
475 FiffStream::SPtr stream(new FiffStream(&file));
476
477 int res;
478
479 if(!stream->open())
480 return FAIL;
481
482 res = mne_read_bad_channel_list_from_node(stream,stream->dirtree(),listp,nlistp);
483
484 stream->close();
485
486 return res;
487}
488
489int mne_sparse_vec_mult2(FiffSparseMatrix* mat, /* The sparse matrix */
490 float *vector, /* Vector to be multiplied */
491 float *res) /* Result of the multiplication */
492/*
493 * Multiply a vector by a sparse matrix using Eigen.
494 */
495{
496 Eigen::Map<const Eigen::VectorXf> vecIn(vector, mat->cols());
497 Eigen::Map<Eigen::VectorXf> vecOut(res, mat->rows());
498 vecOut = mat->eigen() * vecIn;
499 return 0;
500}
501
502int mne_sparse_mat_mult2(FiffSparseMatrix* mat, /* The sparse matrix */
503 const RowMajorMatrixXf& mult, /* Matrix to be multiplied */
504 int ncol, /* How many columns in the above */
505 RowMajorMatrixXf& res) /* Result of the multiplication */
506/*
507 * Multiply a dense matrix by a sparse matrix using Eigen.
508 */
509{
510 Q_UNUSED(ncol);
511 // mat->eigen() is column-major sparse, mult is row-major dense
512 // Result: res = sparse * mult (rows: mat->rows(), cols: mult.cols())
513 res = mat->eigen() * mult;
514 return 0;
515}
516
517#define APPROX_RING_BUF_SIZE (600*1024*1024)
518
519static int approx_ring_buf_size = APPROX_RING_BUF_SIZE;
520
521//=============================================================================================================
522// DEFINE MEMBER METHODS
523//=============================================================================================================
524
526:info(nullptr)
527,nbad(0)
528,first_samp(0)
529,omit_samp(0)
532,nsamp(0)
533,proj(nullptr)
534,sss(nullptr)
535,comp(nullptr)
538,max_event(0)
540,deriv(nullptr)
541,deriv_matched(nullptr)
542{
543}
544
545//=============================================================================================================
546
548{
549// fiff_close(this->file);
550 if (this->stream)
551 this->stream->close();
552 this->filename.clear();
553 this->ch_names.clear();
554
555 this->badlist.clear();
556
557 this->dig_trigger.clear();
558 this->event_list.reset();
559}
560
561//=============================================================================================================
562
563void MNERawData::add_filter_response(int *highpass_effective)
564/*
565 * Add the standard filter frequency response function
566 */
567{
568 /*
569 * Free the previous filter definition
570 */
571 filter_data.reset();
572 /*
573 * Nothing more to do if there is no filter
574 */
575 if (!filter)
576 return;
577 /*
578 * Create a new one
579 */
581 info->sfreq,
582 highpass_effective);
583}
584
585//=============================================================================================================
586
588/*
589 * These will hold the filtered data
590 */
591{
592 MNEFilterDef* filter = this->filter.get();
593 int nfilt_buf;
594 int k;
595 int firstsamp;
596 int nring_buf;
597 int highpass_effective;
598
599 this->filt_bufs.clear();
600 this->filt_ring.reset();
601
602 if (!this->filter || filter->size <= 0)
603 return;
604
605 for (nfilt_buf = 0, firstsamp = this->first_samp-filter->taper_size;
606 firstsamp < this->nsamp + this->first_samp;
607 firstsamp = firstsamp + filter->size)
608 nfilt_buf++;
609#ifdef DEBUG
610 qInfo("%d filter buffers needed\n",nfilt_buf);
611#endif
612 this->filt_bufs.resize(nfilt_buf);
613 for (k = 0, firstsamp = this->first_samp-filter->taper_size; k < nfilt_buf; k++,
614 firstsamp = firstsamp + filter->size) {
615 filt_bufs[k].ns = filter->size + 2*filter->taper_size;
616 filt_bufs[k].firsts = firstsamp;
617 filt_bufs[k].lasts = firstsamp + filt_bufs[k].ns - 1;
618 // bufs[k].ent = NULL;
619 filt_bufs[k].nchan = this->info->nchan;
620 filt_bufs[k].is_skip = false;
621 filt_bufs[k].valid = false;
622 filt_bufs[k].ch_filtered = Eigen::VectorXi::Zero(this->info->nchan);
623 filt_bufs[k].comp_status = MNE_CTFV_NOGRAD;
624 }
625 nring_buf = approx_ring_buf_size/((2*filter->taper_size+filter->size)*
626 this->info->nchan*sizeof(float));
627 this->filt_ring = std::make_unique<RingBuffer>(nring_buf);
628 add_filter_response(&highpass_effective);
629
630 return;
631}
632
633//=============================================================================================================
634
636/*
637 * load just one
638 */
639{
640 if (buf->ent->kind == FIFF_DATA_SKIP) {
641 qCritical("Cannot load a skip");
642 return FAIL;
643 }
644 if (buf->vals.size() == 0) { /* The data space may have been reused */
645 buf->valid = false;
646 ring->allocate(buf->nchan,buf->ns,&buf->vals);
647 }
648 if (buf->valid)
649 return OK;
650
651#ifdef DEBUG
652 qDebug("Read buffer %d .. %d\n",buf->firsts,buf->lasts);
653#endif
654
656 buf->ent,
657 buf->vals,
658 buf->nchan,
659 buf->ns,
660 info->chInfo,
661 nullptr,0) != OK) {
662 buf->valid = false;
663 return FAIL;
664 }
665 buf->valid = true;
666 buf->comp_status = comp_file;
667 return OK;
668}
669
670//=============================================================================================================
671
673/*
674 * Apply compensation channels
675 */
676{
677
678 if (!comp)
679 return OK;
680 if (!comp->undo && !comp->current)
681 return OK;
682 if (buf->comp_status == comp_now)
683 return OK;
684 if (buf->vals.size() == 0)
685 return OK;
686 /*
687 * vals is now a RowMajorMatrixXf — wrap in a column-major MatrixXf for compensation
688 */
689 {
690 Eigen::MatrixXf dataMat = buf->vals; /* implicit copy/conversion */
691
692 if (comp->undo) {
693 std::swap(comp->current, comp->undo);
694 /*
695 * Undo the previous compensation
696 */
697 if (comp->apply_transpose(false, dataMat) != OK) {
698 std::swap(comp->current, comp->undo);
699 return FAIL;
700 }
701 std::swap(comp->current, comp->undo);
702 }
703 if (comp->current) {
704 /*
705 * Apply new compensation
706 */
707 if (comp->apply_transpose(true, dataMat) != OK)
708 return FAIL;
709 }
710 /*
711 * Copy result back to buf->vals
712 */
713 buf->vals = dataMat;
714 }
715 buf->comp_status = comp_now;
716 return OK;
717}
718
719//=============================================================================================================
720
721int MNERawData::pick_data(mneChSelection sel, int firsts, int ns, float **picked)
722/*
723 * Data from a selection of channels
724 */
725{
726 int k,s,p,start,c,fills;
727 int ns2,s2;
728 MNERawBufDef* this_buf;
729 float *values;
730 int need_some;
731
732 RowMajorMatrixXf deriv_vals;
733 int deriv_ns = 0;
734 int nderiv = 0;
735
736 if (firsts < first_samp) {
737 for (s = 0, p = firsts; p < first_samp; s++, p++) {
738 if (sel)
739 for (c = 0; c < sel->nchan; c++)
740 picked[c][s] = 0.0;
741 else
742 for (c = 0; c < info->nchan; c++)
743 picked[c][s] = 0.0;
744 }
745 ns = ns - s;
746 firsts = first_samp;
747 }
748 else
749 s = 0;
750 /*
751 * There is possibly nothing to do
752 */
753 if (sel) {
754 for (c = 0, need_some = false; c < sel->nchan; c++) {
755 if (sel->pick[c] >= 0 || sel->pick_deriv[c] >= 0) {
756 need_some = true;
757 break;
758 }
759 }
760 if (!need_some)
761 return OK;
762 }
763 /*
764 * Have to to the hard work
765 */
766 for (k = 0, this_buf = bufs.data(), s = 0; k < static_cast<int>(bufs.size()); k++, this_buf++) {
767 if (this_buf->lasts >= firsts) {
768 start = firsts - this_buf->firsts;
769 if (start < 0)
770 start = 0;
771 if (this_buf->is_skip) {
772 for (p = start; p < this_buf->ns && ns > 0; p++, ns--, s++) {
773 if (sel) {
774 for (c = 0; c < sel->nchan; c++)
775 if (sel->pick[c] >= 0)
776 picked[c][s] = 0.0;
777 }
778 else {
779 for (c = 0; c < info->nchan; c++)
780 picked[c][s] = 0.0;
781 }
782 }
783 }
784 else {
785 /*
786 * Load the buffer
787 */
788 if (load_one_buffer(this_buf) != OK)
789 return FAIL;
790 /*
791 * Apply compensation
792 */
793 if (compensate_buffer(this_buf) != OK)
794 return FAIL;
795 ns2 = s2 = 0;
796 if (sel) {
797 /*
798 * Do we need the derived channels?
799 */
800 if (sel->nderiv > 0 && deriv_matched) {
801 if (deriv_ns < this_buf->ns || nderiv != deriv_matched->deriv_data->nrow) {
802 deriv_vals.resize(deriv_matched->deriv_data->nrow, this_buf->ns);
803 nderiv = deriv_matched->deriv_data->nrow;
804 deriv_ns = this_buf->ns;
805 }
806 if (mne_sparse_mat_mult2(deriv_matched->deriv_data->data.get(),this_buf->vals,this_buf->ns,deriv_vals) == FAIL) {
807 return FAIL;
808 }
809 }
810 for (c = 0; c < sel->nchan; c++) {
811 /*
812 * First pick the ordinary channels...
813 */
814 if (sel->pick[c] >= 0) {
815 for (p = start, s2 = s, ns2 = ns; p < this_buf->ns && ns2 > 0; p++, ns2--, s2++)
816 picked[c][s2] = this_buf->vals(sel->pick[c], p);
817 }
818 /*
819 * ...then the derived ones
820 */
821 else if (sel->pick_deriv[c] >= 0 && deriv_matched) {
822 for (p = start, s2 = s, ns2 = ns; p < this_buf->ns && ns2 > 0; p++, ns2--, s2++)
823 picked[c][s2] = deriv_vals(sel->pick_deriv[c], p);
824 }
825 }
826 }
827 else {
828 for (c = 0; c < info->nchan; c++)
829 for (p = start, s2 = s, ns2 = ns; p < this_buf->ns && ns2 > 0; p++, ns2--, s2++)
830 picked[c][s2] = this_buf->vals(c, p);
831 }
832 s = s2;
833 ns = ns2;
834 }
835 if (ns == 0)
836 break;
837 }
838 }
839 /*
840 * Extend with the last available sample or zero if the request is beyond the data
841 */
842 if (s > 0) {
843 fills = s-1;
844 for (; ns > 0; ns--, s++) {
845 if (sel)
846 for (c = 0; c < sel->nchan; c++)
847 picked[c][s] = picked[c][fills];
848 else
849 for (c = 0; c < info->nchan; c++)
850 picked[c][s] = picked[c][fills];
851 }
852 }
853 else {
854 for (; ns > 0; ns--, s++) {
855 if (sel)
856 for (c = 0; c < sel->nchan; c++)
857 picked[c][s] = 0;
858 else
859 for (c = 0; c < info->nchan; c++)
860 picked[c][s] = 0;
861 }
862 }
863 return OK;
864}
865
866//=============================================================================================================
867
868int MNERawData::pick_data_proj(mneChSelection sel, int firsts, int ns, float **picked)
869/*
870 * Data from a set of channels, apply projection
871 */
872{
873 int k,s,p,start,c,fills;
874 MNERawBufDef* this_buf;
875 RowMajorMatrixXf *values;
876 Eigen::VectorXf deriv_pvalues_vec;
877
878 if (!proj || (sel && !proj->affect(sel->chspick,sel->nchan) && !proj->affect(sel->chspick_nospace,sel->nchan)))
879 return pick_data(sel,firsts,ns,picked);
880
881 if (firsts < first_samp) {
882 for (s = 0, p = firsts; p < first_samp; s++, p++) {
883 if (sel)
884 for (c = 0; c < sel->nchan; c++)
885 picked[c][s] = 0.0;
886 else
887 for (c = 0; c < info->nchan; c++)
888 picked[c][s] = 0.0;
889 }
890 ns = ns - s;
891 firsts = first_samp;
892 }
893 else
894 s = 0;
895 Eigen::VectorXf pvalues(info->nchan);
896 for (k = 0, this_buf = bufs.data(); k < static_cast<int>(bufs.size()); k++, this_buf++) {
897 if (this_buf->lasts >= firsts) {
898 start = firsts - this_buf->firsts;
899 if (start < 0)
900 start = 0;
901 if (this_buf->is_skip) {
902 for (p = start; p < this_buf->ns && ns > 0; p++, ns--, s++) {
903 if (sel) {
904 for (c = 0; c < sel->nchan; c++)
905 if (sel->pick[c] >= 0)
906 picked[c][s] = 0.0;
907 }
908 else {
909 for (c = 0; c < info->nchan; c++)
910 picked[c][s] = 0.0;
911 }
912 }
913 }
914 else {
915 /*
916 * Load the buffer
917 */
918 if (load_one_buffer(this_buf) != OK)
919 return FAIL;
920 /*
921 * Apply compensation
922 */
923 if (compensate_buffer(this_buf) != OK)
924 return FAIL;
925 /*
926 * Apply projection
927 */
928 values = &this_buf->vals;
929 if (sel && sel->nderiv > 0 && deriv_matched) {
930 deriv_pvalues_vec.resize(deriv_matched->deriv_data->nrow);
931 }
932 for (p = start; p < this_buf->ns && ns > 0; p++, ns--, s++) {
933 for (c = 0; c < info->nchan; c++)
934 pvalues[c] = (*values)(c,p);
935 if (proj->project_vector(pvalues,true) != OK)
936 qWarning()<<"Error";
937 if (sel) {
938 if (sel->nderiv > 0 && deriv_matched) {
939 if (mne_sparse_vec_mult2(deriv_matched->deriv_data->data.get(),pvalues.data(),deriv_pvalues_vec.data()) == FAIL)
940 return FAIL;
941 }
942 for (c = 0; c < sel->nchan; c++) {
943 /*
944 * First try the ordinary channels...
945 */
946 if (sel->pick[c] >= 0)
947 picked[c][s] = pvalues[sel->pick[c]];
948 /*
949 * ...then the derived ones
950 */
951 else if (sel->pick_deriv[c] >= 0 && deriv_matched)
952 picked[c][s] = deriv_pvalues_vec[sel->pick_deriv[c]];
953 }
954 }
955 else {
956 for (c = 0; c < info->nchan; c++) {
957 picked[c][s] = pvalues[c];
958 }
959 }
960 }
961 }
962 if (ns == 0)
963 break;
964 }
965 }
966
967 /*
968 * Extend with the last available sample or zero if the request is beyond the data
969 */
970 if (s > 0) {
971 fills = s-1;
972 for (; ns > 0; ns--, s++) {
973 if (sel)
974 for (c = 0; c < sel->nchan; c++)
975 picked[c][s] = picked[c][fills];
976 else
977 for (c = 0; c < info->nchan; c++)
978 picked[c][s] = picked[c][fills];
979 }
980 }
981 else {
982 for (; ns > 0; ns--, s++) {
983 if (sel)
984 for (c = 0; c < sel->nchan; c++)
985 picked[c][s] = 0;
986 else
987 for (c = 0; c < info->nchan; c++)
988 picked[c][s] = 0;
989 }
990 }
991 return OK;
992}
993
994//=============================================================================================================
995
997/*
998 * Load and filter one buffer
999 */
1000{
1001 int k;
1002 int res;
1003
1004 if (buf->vals.size() == 0) {
1005 buf->valid = false;
1006 filt_ring->allocate(buf->nchan, buf->ns,&buf->vals);
1007 }
1008 if (buf->valid)
1009 return OK;
1010
1011 std::vector<float*> vals_storage(buf->nchan);
1012 float **vals = vals_storage.data();
1013 for (k = 0; k < buf->nchan; k++) {
1014 buf->ch_filtered[k] = false;
1015 vals[k] = buf->vals.row(k).data() + filter->taper_size;
1016 }
1017
1018 res = pick_data_proj(nullptr,buf->firsts + filter->taper_size,buf->ns - 2*filter->taper_size,vals);
1019
1020#ifdef DEBUG
1021 if (res == OK)
1022 qDebug("Loaded filtered buffer %d...%d %d %d last = %d\n",
1023 buf->firsts,buf->lasts,buf->lasts-buf->firsts+1,buf->ns,first_samp + nsamp);
1024#endif
1025 buf->valid = res == OK;
1026 return res;
1027}
1028
1029//=============================================================================================================
1030
1031int MNERawData::pick_data_filt(mneChSelection sel, int firsts, int ns, float **picked)
1032/*
1033 * Data for a selection (filtered and picked)
1034 */
1035{
1036 int k,s,bs,c;
1037 int bs1,bs2,s1,s2,lasts;
1038 MNERawBufDef* this_buf;
1039 float *values;
1040 RowMajorMatrixXf deriv_vals;
1041 Eigen::VectorXf dc;
1042 float dc_offset;
1043 int deriv_ns = 0;
1044 int nderiv = 0;
1045 int filter_was;
1046
1047 if (!filter->filter_on)
1048 return pick_data_proj(sel,firsts,ns,picked);
1049
1050 if (sel) {
1051 for (s = 0; s < ns; s++)
1052 for (c = 0; c < sel->nchan; c++)
1053 picked[c][s] = 0.0;
1054 }
1055 else {
1056 for (s = 0; s < ns; s++)
1057 for (c = 0; c < info->nchan; c++)
1058 picked[c][s] = 0.0;
1059 }
1060 lasts = firsts + ns - 1;
1061 /*
1062 * Take into account the initial dc offset (compensate and project)
1063 */
1064 if (first_sample_val.size() > 0) {
1065 dc = first_sample_val;
1066 /*
1067 * Is this correct??
1068 */
1069 if (comp && comp->current)
1070 if (comp->apply(true,dc) != OK)
1071 return FAIL;
1072 if (proj)
1073 if (proj->project_vector(dc,true) != OK)
1074 return FAIL;
1075 }
1076 filter_was = filter->filter_on;
1077 /*
1078 * Find the first buffer to consider
1079 */
1080 for (k = 0, this_buf = filt_bufs.data(); k < static_cast<int>(filt_bufs.size()); k++, this_buf++) {
1081 if (this_buf->lasts >= firsts)
1082 break;
1083 }
1084 for (; k < static_cast<int>(filt_bufs.size()) && this_buf->firsts <= lasts; k++, this_buf++) {
1085#ifdef DEBUG
1086 qDebug("this_buf (%d): %d..%d\n",k,this_buf->firsts,this_buf->lasts);
1087#endif
1088 /*
1089 * Load the buffer first and apply projection
1090 */
1091 if (load_one_filt_buf(this_buf) != OK)
1092 return FAIL;
1093 /*
1094 * Then filter all relevant channels (not stimuli)
1095 */
1096 if (sel) {
1097 for (c = 0; c < sel->nchan; c++) {
1098 if (sel->pick[c] >= 0) {
1099 if (!this_buf->ch_filtered[sel->pick[c]]) {
1100 /*
1101 * Do not filter stimulus channels
1102 */
1103 dc_offset = 0.0;
1104 if (info->chInfo[sel->pick[c]].kind == FIFFV_STIM_CH)
1105 filter->filter_on = false;
1106 else if (dc.size() > 0)
1107 dc_offset = dc[sel->pick[c]];
1108 if (mne_apply_filter(*filter,filter_data.get(),this_buf->vals.row(sel->pick[c]).data(),this_buf->ns,true,
1109 dc_offset,info->chInfo[sel->pick[c]].kind) != OK) {
1110 filter->filter_on = filter_was;
1111 return FAIL;
1112 }
1113 this_buf->ch_filtered[sel->pick[c]] = true;
1114 filter->filter_on = filter_was;
1115 }
1116 }
1117 }
1118 /*
1119 * Also check channels included in derivations if they are used
1120 */
1121 if (sel->nderiv > 0 && deriv_matched) {
1122 MNEDeriv* der = deriv_matched.get();
1123 for (c = 0; c < der->deriv_data->ncol; c++) {
1124 if (der->in_use[c] > 0 &&
1125 !this_buf->ch_filtered[c]) {
1126 /*
1127 * Do not filter stimulus channels
1128 */
1129 dc_offset = 0.0;
1130 if (info->chInfo[c].kind == FIFFV_STIM_CH)
1131 filter->filter_on = false;
1132 else if (dc.size() > 0)
1133 dc_offset = dc[c];
1134 if (mne_apply_filter(*filter,filter_data.get(),this_buf->vals.row(c).data(),this_buf->ns,true,
1135 dc_offset,info->chInfo[c].kind) != OK) {
1136 filter->filter_on = filter_was;
1137 return FAIL;
1138 }
1139 this_buf->ch_filtered[c] = true;
1140 filter->filter_on = filter_was;
1141 }
1142 }
1143 }
1144 }
1145 else {
1146 /*
1147 * Simply filter all channels if there is no selection
1148 */
1149 for (c = 0; c < info->nchan; c++) {
1150 if (!this_buf->ch_filtered[c]) {
1151 /*
1152 * Do not filter stimulus channels
1153 */
1154 dc_offset = 0.0;
1155 if (info->chInfo[c].kind == FIFFV_STIM_CH)
1156 filter->filter_on = false;
1157 else if (dc.size() > 0)
1158 dc_offset = dc[c];
1159 if (mne_apply_filter(*filter,filter_data.get(),this_buf->vals.row(c).data(),this_buf->ns,true,
1160 dc_offset,info->chInfo[c].kind) != OK) {
1161 filter->filter_on = filter_was;
1162 return FAIL;
1163 }
1164 this_buf->ch_filtered[c] = true;
1165 filter->filter_on = filter_was;
1166 }
1167 }
1168 }
1169 /*
1170 * Decide the picking limits
1171 */
1172 if (firsts >= this_buf->firsts) {
1173 bs1 = firsts - this_buf->firsts;
1174 s1 = 0;
1175 }
1176 else {
1177 bs1 = 0;
1178 s1 = this_buf->firsts - firsts;
1179 }
1180 if (lasts >= this_buf->lasts) {
1181 bs2 = this_buf->ns;
1182 s2 = this_buf->lasts - lasts + ns;
1183 }
1184 else {
1185 bs2 = lasts - this_buf->lasts + this_buf->ns;
1186 s2 = ns;
1187 }
1188#ifdef DEBUG
1189 qDebug("buf : %d..%d %d\n",bs1,bs2,bs2-bs1);
1190 qDebug("dest : %d..%d %d\n",s1,s2,s2-s1);
1191#endif
1192 /*
1193 * Then pick data from all relevant channels
1194 */
1195 if (sel) {
1196 if (sel->nderiv > 0 && deriv_matched) {
1197 /*
1198 * Compute derived data if we need it
1199 */
1200 if (deriv_ns < this_buf->ns || nderiv != deriv_matched->deriv_data->nrow) {
1201 deriv_vals.resize(deriv_matched->deriv_data->nrow, this_buf->ns);
1202 nderiv = deriv_matched->deriv_data->nrow;
1203 deriv_ns = this_buf->ns;
1204 }
1205 if (mne_sparse_mat_mult2(deriv_matched->deriv_data->data.get(),this_buf->vals,this_buf->ns,deriv_vals) == FAIL)
1206 return FAIL;
1207 }
1208 for (c = 0; c < sel->nchan; c++) {
1209 /*
1210 * First the ordinary channels
1211 */
1212 if (sel->pick[c] >= 0) {
1213 values = this_buf->vals.row(sel->pick[c]).data();
1214 for (s = s1, bs = bs1; s < s2; s++, bs++)
1215 picked[c][s] += values[bs];
1216 }
1217 else if (sel->pick_deriv[c] >= 0 && deriv_matched) {
1218 for (s = s1, bs = bs1; s < s2; s++, bs++)
1219 picked[c][s] += deriv_vals(sel->pick_deriv[c], bs);
1220 }
1221 }
1222 }
1223 else {
1224 for (c = 0; c < info->nchan; c++) {
1225 values = this_buf->vals.row(c).data();
1226 for (s = s1, bs = bs1; s < s2; s++, bs++)
1227 picked[c][s] += values[bs];
1228 }
1229 }
1230 }
1231 (void)bs2; // squash compiler warning, this is unused
1232 return OK;
1233}
1234
1235//=============================================================================================================
1236
1238 int omit_skip,
1239 int allow_maxshield,
1240 const MNEFilterDef& filter,
1241 int comp_set)
1242/*
1243 * Open a raw data file
1244 */
1245{
1246 std::unique_ptr<MNERawInfo> info;
1247 std::unique_ptr<MNERawData> data;
1248
1249 auto filePtr = std::make_unique<QFile>(name);
1250 FiffStream::SPtr stream(new FiffStream(filePtr.get()));
1251 // fiffFile in = NULL;
1252
1254 QList<FiffDirEntry::SPtr> dir0;
1255 // fiffTagRec tag;
1256 FiffTag::UPtr t_pTag;
1257 FiffChInfo ch;
1258 int k, b, nbuf, ndir;
1259 int current_dir0 = 0;
1260
1261 // tag.data = NULL;
1262
1263 if (MNERawInfo::load(name,allow_maxshield,info) == FAIL)
1264 return nullptr;
1265
1266 for (k = 0; k < info->nchan; k++) {
1267 ch = info->chInfo.at(k);
1268 if (QString::compare(ch.ch_name,MNE_DEFAULT_TRIGGER_CH) == 0) {
1269 if (std::fabs(1.0 - ch.range) > 1e-5) {
1270 ch.range = 1.0;
1271 qInfo("%s range set to %f\n",MNE_DEFAULT_TRIGGER_CH,ch.range);
1272 }
1273 }
1274 /*
1275 * Take care of the nonzero unit multiplier
1276 */
1277 if (ch.unit_mul != 0) {
1278 ch.cal = pow(10.0,static_cast<double>(ch.unit_mul))*ch.cal;
1279 qInfo("Ch %s unit multiplier %d -> 0\n",ch.ch_name.toLatin1().data(),ch.unit_mul);
1280 ch.unit_mul = 0;
1281 }
1282 }
1283 // if ((in = fiff_open(name)) == NULL)
1284 // goto bad;
1285 if(!stream->open())
1286 return nullptr;
1287
1288 data = std::make_unique<MNERawData>();
1289 data->filename = name;
1290 data->file = std::move(filePtr);
1291 data->stream = stream;
1292 data->info = std::move(info);
1293 /*
1294 * Add the channel name list
1295 */
1296 data->ch_names.clear();
1297 for (int i = 0; i < data->info->nchan; i++)
1298 data->ch_names.append(data->info->chInfo[i].ch_name);
1299 if (data->ch_names.size() != data->info->nchan) {
1300 qCritical("Channel names were not translated correctly into a name list");
1301 return nullptr;
1302 }
1303 /*
1304 * Compensation data
1305 */
1306 data->comp = MNECTFCompDataSet::read(data->filename);
1307 if (data->comp) {
1308 if (data->comp->ncomp > 0)
1309 qInfo("Read %d compensation data sets from %s\n",data->comp->ncomp,data->filename.toUtf8().constData());
1310 else
1311 qInfo("No compensation data in %s\n",data->filename.toUtf8().constData());
1312 }
1313 else
1314 qWarning() << "err_print_error()";
1315 if ((data->comp_file = MNECTFCompDataSet::get_comp(data->info->chInfo,data->info->nchan)) == FAIL)
1316 return nullptr;
1317 qInfo("Compensation in file : %s\n",MNECTFCompDataSet::explain_comp(MNECTFCompDataSet::map_comp_kind(data->comp_file)).toUtf8().constData());
1318 if (comp_set < 0)
1319 data->comp_now = data->comp_file;
1320 else
1321 data->comp_now = comp_set;
1322
1323 if (!data->comp) {
1324 if (data->comp_now != MNE_CTFV_NOGRAD) {
1325 qCritical("Cannot do compensation because compensation data are missing");
1326 return nullptr;
1327 }
1328 } else if (data->comp->set_compensation(data->comp_now,
1329 data->info->chInfo,
1330 data->info->nchan,
1331 QList<FIFFLIB::FiffChInfo>(),
1332 0) == FAIL)
1333 return nullptr;
1334 /*
1335 * SSS data
1336 */
1337 data->sss = MNESssData::read(data->filename);
1338 if (data->sss && data->sss->job != FIFFV_SSS_JOB_NOTHING && data->sss->comp_info.size() > 0) {
1339 qInfo("SSS data read from %s :\n",data->filename.toUtf8().constData());
1340 QTextStream errStream(stderr);
1341 data->sss->print(errStream);
1342 }
1343 else {
1344 qInfo("No SSS data in %s\n",data->filename.toUtf8().constData());
1345 data->sss.reset();
1346 }
1347 /*
1348 * Buffers
1349 */
1350 dir0 = data->info->rawDir;
1351 ndir = data->info->ndir;
1352 /*
1353 * Take into account the first sample
1354 */
1355 if (dir0[current_dir0]->kind == FIFF_FIRST_SAMPLE) {
1356 // if (fiff_read_this_tag(in->fd,dir0->pos,&tag) == FIFF_FAIL)
1357 // goto bad;
1358 if (!stream->read_tag(t_pTag,dir0[current_dir0]->pos))
1359 return nullptr;
1360 data->first_samp = *t_pTag->toInt();
1361 current_dir0++;
1362 ndir--;
1363 }
1364 if (dir0[current_dir0]->kind == FIFF_DATA_SKIP) {
1365 int nsamp_skip;
1366 // if (fiff_read_this_tag(in->fd,dir0->pos,&tag) == FIFF_FAIL)
1367 // goto bad;
1368 if (!stream->read_tag(t_pTag,dir0[current_dir0]->pos))
1369 return nullptr;
1370 nsamp_skip = data->info->buf_size*(*t_pTag->toInt());
1371 qInfo("Data skip of %d samples in the beginning\n",nsamp_skip);
1372 current_dir0++;
1373 ndir--;
1374 if (dir0[current_dir0]->kind == FIFF_FIRST_SAMPLE) {
1375 // if (fiff_read_this_tag(in->fd,dir0->pos,&tag) == FIFF_FAIL)
1376 // goto bad;
1377 if (!stream->read_tag(t_pTag,dir0[current_dir0]->pos))
1378 return nullptr;
1379 data->first_samp += *t_pTag->toInt();
1380 current_dir0++;
1381 ndir--;
1382 }
1383 if (omit_skip) {
1384 data->omit_samp = data->first_samp + nsamp_skip;
1385 data->omit_samp_old = nsamp_skip;
1386 data->first_samp = 0;
1387 }
1388 else {
1389 data->first_samp = data->first_samp + nsamp_skip;
1390 }
1391 }
1392 else if (omit_skip) {
1393 data->omit_samp = data->first_samp;
1394 data->first_samp = 0;
1395 }
1396#ifdef DEBUG
1397 qInfo("data->first_samp = %d\n",data->first_samp);
1398#endif
1399 /*
1400 * Figure out the buffers
1401 */
1402 // for (k = 0, dir = dir0, nbuf = 0; k < ndir; k++, dir++)
1403 for (k = 0, nbuf = 0; k < ndir; k++)
1404 if (dir0[k]->kind == FIFF_DATA_BUFFER ||
1405 dir0[k]->kind == FIFF_DATA_SKIP)
1406 nbuf++;
1407 data->bufs.resize(nbuf);
1408
1409 // for (k = 0, nbuf = 0, dir = dir0; k < ndir; k++, dir++)
1410 for (k = 0, nbuf = 0; k < ndir; k++)
1411 if (dir0[k]->kind == FIFF_DATA_BUFFER ||
1412 dir0[k]->kind == FIFF_DATA_SKIP) {
1413 data->bufs[nbuf].ns = 0;
1414 data->bufs[nbuf].ent = dir0[k];
1415 data->bufs[nbuf].nchan = data->info->nchan;
1416 data->bufs[nbuf].is_skip = dir0[k]->kind == FIFF_DATA_SKIP;
1417 data->bufs[nbuf].valid = false;
1418 data->bufs[nbuf].comp_status = data->comp_file;
1419 nbuf++;
1420 }
1421 data->nsamp = 0;
1422 for (k = 0; k < nbuf; k++) {
1423 dir = data->bufs[k].ent;
1424 if (dir->kind == FIFF_DATA_BUFFER) {
1425 if (dir->type == FIFFT_DAU_PACK16 || dir->type == FIFFT_SHORT)
1426 data->bufs[k].ns = dir->size/(data->info->nchan*sizeof(fiff_dau_pack16_t));
1427 else if (dir->type == FIFFT_FLOAT)
1428 data->bufs[k].ns = dir->size/(data->info->nchan*sizeof(fiff_float_t));
1429 else if (dir->type == FIFFT_INT)
1430 data->bufs[k].ns = dir->size/(data->info->nchan*sizeof(fiff_int_t));
1431 else {
1432 qCritical("We are not prepared to handle raw data type: %d",dir->type);
1433 return nullptr;
1434 }
1435 }
1436 else if (dir->kind == FIFF_DATA_SKIP) {
1437 // if (fiff_read_this_tag(in->fd,dir->pos,&tag) == FIFF_FAIL)
1438 // goto bad;
1439 if (!stream->read_tag(t_pTag,dir->pos))
1440 return nullptr;
1441 data->bufs[k].ns = data->info->buf_size*(*t_pTag->toInt());
1442 }
1443 data->bufs[k].firsts = k == 0 ? data->first_samp : data->bufs[k-1].lasts + 1;
1444 data->bufs[k].lasts = data->bufs[k].firsts + data->bufs[k].ns - 1;
1445 data->nsamp += data->bufs[k].ns;
1446 }
1447 // FREE_36(tag.data);
1448 /*
1449 * Set up the first sample values
1450 */
1451 data->bad = Eigen::VectorXi::Zero(data->info->nchan);
1452 data->offsets = Eigen::VectorXf::Zero(data->info->nchan);
1453 /*
1454 * Th bad channel stuff
1455 */
1456 {
1457 if (mne_read_bad_channel_list(name,data->badlist,data->nbad) == OK) {
1458 for (b = 0; b < data->nbad; b++) {
1459 for (k = 0; k < data->info->nchan; k++) {
1460 if (QString::compare(data->info->chInfo[k].ch_name,data->badlist[b],Qt::CaseInsensitive) == 0) {
1461 data->bad[k] = 1;
1462 break;
1463 }
1464 }
1465 }
1466 qInfo("%d bad channels read from %s%s",data->nbad,name.toUtf8().constData(),data->nbad > 0 ? ":\n" : "\n");
1467 if (data->nbad > 0) {
1468 qInfo("\t");
1469 for (k = 0; k < data->nbad; k++)
1470 qInfo("%s%c",data->badlist[k].toUtf8().constData(),k < data->nbad-1 ? ' ' : '\n');
1471 }
1472 }
1473 }
1474 /*
1475 * Initialize the raw data buffers
1476 */
1477 nbuf = approx_ring_buf_size/(data->info->buf_size*data->info->nchan*sizeof(float));
1478 data->ring = std::make_unique<RingBuffer>(nbuf);
1479 /*
1480 * Initialize the filter buffers
1481 */
1482 data->filter = std::make_unique<MNEFilterDef>(filter);
1483 data->setup_filter_bufs();
1484
1485 {
1486 std::vector<float> vals_storage(data->info->nchan, 0.0f);
1487 std::vector<float*> vals_rows(data->info->nchan);
1488 for (int i = 0; i < data->info->nchan; i++)
1489 vals_rows[i] = &vals_storage[i];
1490 float **vals = vals_rows.data();
1491
1492 if (data->pick_data(nullptr,data->first_samp,1,vals) == FAIL)
1493 return nullptr;
1494 data->first_sample_val.resize(data->info->nchan);
1495 for (k = 0; k < data->info->nchan; k++)
1496 data->first_sample_val[k] = vals[k][0];
1497 qInfo("Initial dc offsets determined\n");
1498 }
1499 qInfo("Raw data file %s:\n",name.toUtf8().constData());
1500 qInfo("\tnchan = %d\n",data->info->nchan);
1501 qInfo("\tnsamp = %d\n",data->nsamp);
1502 qInfo("\tsfreq = %-8.3f Hz\n",data->info->sfreq);
1503 qInfo("\tlength = %-8.3f sec\n",data->nsamp/data->info->sfreq);
1504
1505 return data.release();
1506}
1507
1508//=============================================================================================================
1509
1510MNERawData *MNERawData::open_file(const QString& name, int omit_skip, int allow_maxshield, const MNEFilterDef& filter)
1511/*
1512 * Wrapper for open_file to work as before
1513 */
1514{
1515 return open_file_comp(name,omit_skip,allow_maxshield,filter,-1);
1516}
#define M_PI
constexpr int FAIL
constexpr int OK
#define FIFF_DATA_BUFFER
Definition fiff_file.h:556
#define FIFFT_INT
Definition fiff_file.h:231
#define FIFFT_SHORT
Definition fiff_file.h:230
#define FIFF_FIRST_SAMPLE
Definition fiff_file.h:461
#define FIFFT_DAU_PACK16
Definition fiff_file.h:243
#define FIFFT_FLOAT
Definition fiff_file.h:232
#define FIFFV_SSS_JOB_NOTHING
Definition fiff_file.h:538
#define FIFF_DATA_SKIP
Definition fiff_file.h:557
#define FIFFV_EOG_CH
#define FIFF_MNE_CH_NAME_LIST
#define FIFFV_STIM_CH
#define FIFFB_MNE_BAD_CHANNELS
#define MNE_DEFAULT_TRIGGER_CH
Default digital trigger channel name.
Definition mne_types.h:149
#define MNE_CTFV_NOGRAD
Definition mne_types.h:131
MNERawData class declaration.
void mne_fft_syn(float *data, int np, std::vector< float > &)
int mne_sparse_vec_mult2(FiffSparseMatrix *mat, float *vector, float *res)
int mne_read_raw_buffer_t(FiffStream::SPtr &stream, const FiffDirEntry::SPtr &ent, RowMajorMatrixXf &data, int nchan, int nsamp, const QList< FIFFLIB::FiffChInfo > &chs, int *pickno, int npick)
int mne_sparse_mat_mult2(FiffSparseMatrix *mat, const RowMajorMatrixXf &mult, int ncol, RowMajorMatrixXf &res)
int mne_apply_filter(const MNEFilterDef &filter, FilterData *d, float *data, int ns, int zero_pad, float dc_offset, int kind)
int mne_read_bad_channel_list(const QString &name, QStringList &listp, int &nlistp)
#define APPROX_RING_BUF_SIZE
int mne_read_bad_channel_list_from_node(FiffStream::SPtr &stream, const FiffDirNode::SPtr &pNode, QStringList &listp, int &nlistp)
int mne_compare_filters(const MNEFilterDef &f1, const MNEFilterDef &f2)
void mne_fft_ana(float *data, int np, std::vector< float > &)
std::unique_ptr< FilterData > mne_create_filter_response(const MNEFilterDef &filter, float sfreq, int *highpass_effective)
Core MNE data structures (source spaces, source estimates, hemispheres).
MNEChSelection * mneChSelection
Eigen::Matrix< float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > RowMajorMatrixXf
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
float fiff_float_t
Definition fiff_types.h:93
qint16 fiff_dau_pack16_t
Definition fiff_types.h:97
qint16 fiff_short_t
Definition fiff_types.h:87
Channel info descriptor.
QSharedPointer< FiffDirEntry > SPtr
QSharedPointer< FiffDirNode > SPtr
FIFF sparse matrix storage backed by Eigen.
Eigen::SparseMatrix< float > & eigen()
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::VectorXi pick_deriv
static std::unique_ptr< MNECTFCompDataSet > read(const QString &name)
static QString explain_comp(int kind)
static int get_comp(const QList< FIFFLIB::FiffChInfo > &chs, int nch)
One item in a derivation data set.
Definition mne_deriv.h:79
std::unique_ptr< MNESparseNamedMatrix > deriv_data
Definition mne_deriv.h:99
Eigen::VectorXi in_use
Definition mne_deriv.h:100
Definition of one raw data buffer within a FIFF file.
FIFFLIB::FiffDirEntry::SPtr ent
Eigen::VectorXi ch_filtered
Eigen::Matrix< float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > vals
void allocate(int nrow, int ncol, RowMajorMatrixXf *res)
std::vector< Entry > entries
RingBuffer(int nslots)
RowMajorMatrixXf * user
Pre-computed frequency-domain filter state used for FFT-based raw data filtering.
FilterData(int resp_size)
std::vector< float > eog_freq_resp
std::vector< float > freq_resp
std::vector< float > precalc
int load_one_filt_buf(MNERawBufDef *buf)
std::unique_ptr< MNELIB::MNEProjOp > proj
std::unique_ptr< MNELIB::MNECTFCompDataSet > comp
std::unique_ptr< FilterData > filter_data
void add_filter_response(int *highpass_effective)
int compensate_buffer(MNERawBufDef *buf)
int pick_data_proj(mneChSelection sel, int firsts, int ns, float **picked)
QStringList ch_names
std::unique_ptr< RingBuffer > filt_ring
std::unique_ptr< MNELIB::MNERawInfo > info
std::unique_ptr< MNELIB::MNEDeriv > deriv_matched
unsigned int dig_trigger_mask
QStringList badlist
std::vector< MNELIB::MNERawBufDef > filt_bufs
static MNERawData * open_file_comp(const QString &name, int omit_skip, int allow_maxshield, const MNEFilterDef &filter, int comp_set)
std::unique_ptr< MNEFilterDef > filter
int load_one_buffer(MNERawBufDef *buf)
Eigen::VectorXf first_sample_val
std::unique_ptr< MNEEventList > event_list
std::unique_ptr< MNELIB::MNEDerivSet > deriv
unsigned int max_event
std::vector< MNELIB::MNERawBufDef > bufs
int pick_data(mneChSelection sel, int firsts, int ns, float **picked)
FIFFLIB::FiffStream::SPtr stream
std::unique_ptr< MNELIB::MNESssData > sss
static MNERawData * open_file(const QString &name, int omit_skip, int allow_maxshield, const MNEFilterDef &filter)
int pick_data_filt(mneChSelection sel, int firsts, int ns, float **picked)
std::unique_ptr< RingBuffer > ring
static int load(const QString &name, int allow_maxshield, std::unique_ptr< MNERawInfo > &infop)
static std::unique_ptr< MNESssData > read(const QString &name)