v2.0.0
Loading...
Searching...
No Matches
fiff_raw_data.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "fiff_raw_data.h"
42#include "fiff_events.h"
43#include "fiff_tag.h"
44#include "fiff_stream.h"
45#include "cstdlib"
46
47//=============================================================================================================
48// USED NAMESPACES
49//=============================================================================================================
50
51using namespace FIFFLIB;
52using namespace Eigen;
53
54//=============================================================================================================
55// DEFINE MEMBER METHODS
56//=============================================================================================================
57
63
64//=============================================================================================================
65
66FiffRawData::FiffRawData(QIODevice &p_IODevice)
67: first_samp(-1)
68, last_samp(-1)
69{
70 //setup FiffRawData object
71 if(!FiffStream::setup_read_raw(p_IODevice, *this))
72 {
73 printf("\tError during fiff setup raw read.\n");
74 //exit(EXIT_FAILURE); //ToDo Throw here, e.g.: throw std::runtime_error("IO Error! File not found");
75 return;
76 }
77}
78
79//=============================================================================================================
80
81FiffRawData::FiffRawData(QIODevice &p_IODevice, bool b_littleEndian)
82: first_samp(-1)
83, last_samp(-1)
84{
85 //setup FiffRawData object
86 if(!FiffStream::setup_read_raw(p_IODevice, *this, false, b_littleEndian))
87 {
88 printf("\tError during fiff setup raw read.\n");
89 //exit(EXIT_FAILURE); //ToDo Throw here, e.g.: throw std::runtime_error("IO Error! File not found");
90 return;
91 }
92}
93
94//=============================================================================================================
95
97: file(p_FiffRawData.file)
98, info(p_FiffRawData.info)
99, first_samp(p_FiffRawData.first_samp)
100, last_samp(p_FiffRawData.last_samp)
101, cals(p_FiffRawData.cals)
102, rawdir(p_FiffRawData.rawdir)
103, proj(p_FiffRawData.proj)
104, comp(p_FiffRawData.comp)
105{
106}
107
108//=============================================================================================================
109
113
114//=============================================================================================================
115
117{
118 info.clear();
119 first_samp = -1;
120 last_samp = -1;
121 cals = RowVectorXd();
122 rawdir.clear();
123 proj = MatrixXd();
124 comp.clear();
125}
126
127//=============================================================================================================
128#include <QElapsedTimer>
129#include <QDebug>
130bool FiffRawData::read_raw_segment(MatrixXd& data,
131 MatrixXd& times,
132 fiff_int_t from,
133 fiff_int_t to,
134 const RowVectorXi& sel,
135 bool do_debug) const
136{
137 bool projAvailable = true;
138
139 if (this->proj.size() == 0) {
140 //qDebug() << "FiffRawData::read_raw_segment - No projectors setup. Consider calling MNE::setup_compensators.";
141 projAvailable = false;
142 }
143
144 if(from == -1)
145 from = this->first_samp;
146 if(to == -1)
147 to = this->last_samp;
148 //
149 // Initial checks
150 //
151 if(from < this->first_samp)
152 from = this->first_samp;
153 if(to > this->last_samp)
154 to = this->last_samp;
155 //
156 if(from > to)
157 {
158 printf("No data in this range %d ... %d = %9.3f ... %9.3f secs...", from, to, ((float)from)/this->info.sfreq, ((float)to)/this->info.sfreq);
159 return false;
160 }
161 //printf("Reading %d ... %d = %9.3f ... %9.3f secs...", from, to, ((float)from)/this->info.sfreq, ((float)to)/this->info.sfreq);
162 //
163 // Initialize the data and calibration vector
164 //
165 qint32 nchan = this->info.nchan;
166 qint32 dest = 0;//1;
167 qint32 i, k, r;
168
169 typedef Eigen::Triplet<double> T;
170 std::vector<T> tripletList;
171 tripletList.reserve(nchan);
172 for(i = 0; i < nchan; ++i)
173 tripletList.push_back(T(i, i, this->cals[i]));
174
175 SparseMatrix<double> cal(nchan, nchan);
176 cal.setFromTriplets(tripletList.begin(), tripletList.end());
177// cal.makeCompressed();
178
179 MatrixXd mult_full;
180 //
181 if (sel.size() == 0)
182 {
183 data = MatrixXd(nchan, to-from+1);
184// data->setZero();
185 if (projAvailable || this->comp.kind != -1)
186 {
187 if (!projAvailable)
188 mult_full = this->comp.data->data*cal;
189 else if (this->comp.kind == -1)
190 mult_full = this->proj*cal;
191 else
192 mult_full = this->proj*this->comp.data->data*cal;
193 }
194 }
195 else
196 {
197 data = MatrixXd(sel.size(),to-from+1);
198// data->setZero();
199
200 MatrixXd selVect(sel.size(), nchan);
201
202 selVect.setZero();
203
204 if (!projAvailable && this->comp.kind == -1)
205 {
206 tripletList.clear();
207 tripletList.reserve(sel.size());
208 for(i = 0; i < sel.size(); ++i)
209 tripletList.push_back(T(i, i, this->cals[sel[i]]));
210 cal = SparseMatrix<double>(sel.size(), sel.size());
211 cal.setFromTriplets(tripletList.begin(), tripletList.end());
212 }
213 else
214 {
215 if (!projAvailable)
216 {
217 qDebug() << "This has to be debugged! #1";
218 for( i = 0; i < sel.size(); ++i)
219 selVect.row(i) = this->comp.data->data.block(sel[i],0,1,nchan);
220 mult_full = selVect*cal;
221 }
222 else if (this->comp.kind == -1)
223 {
224 for( i = 0; i < sel.size(); ++i)
225 selVect.row(i) = this->proj.block(sel[i],0,1,nchan);
226
227 mult_full = selVect*cal;
228 }
229 else
230 {
231 qDebug() << "This has to be debugged! #3";
232 for( i = 0; i < sel.size(); ++i)
233 selVect.row(i) = this->proj.block(sel[i],0,1,nchan);
234
235 mult_full = selVect*this->comp.data->data*cal;
236 }
237 }
238 }
239
240 //
241 // Make mult sparse
242 //
243 tripletList.clear();
244 tripletList.reserve(mult_full.rows()*mult_full.cols());
245 for(i = 0; i < mult_full.rows(); ++i)
246 for(k = 0; k < mult_full.cols(); ++k)
247 if(mult_full(i,k) != 0)
248 tripletList.push_back(T(i, k, mult_full(i,k)));
249
250 SparseMatrix<double> mult(mult_full.rows(),mult_full.cols());
251 if(tripletList.size() > 0)
252 mult.setFromTriplets(tripletList.begin(), tripletList.end());
253// mult.makeCompressed();
254
256 if (!this->file->device()->isOpen())
257 {
258 if (!this->file->device()->open(QIODevice::ReadOnly))
259 {
260 printf("Cannot open file %s",this->info.filename.toUtf8().constData());
261 }
262 fid = this->file;
263 }
264 else
265 {
266 fid = this->file;
267 }
268
269 MatrixXd one, newData, tmp_data;
270 FiffRawDir thisRawDir;
271 FiffTag::SPtr t_pTag;
272 fiff_int_t first_pick, last_pick, picksamp;
273 for(k = 0; k < this->rawdir.size(); ++k)
274 {
275 thisRawDir = this->rawdir[k];
276 //
277 // Do we need this buffer
278 //
279 if (thisRawDir.last > from)
280 {
281 if (thisRawDir.ent->kind == -1)
282 {
283 //
284 // Take the easy route: skip is translated to zeros
285 //
286 if(do_debug)
287 printf("S");
288 if (sel.cols() <= 0)
289 one.resize(nchan,thisRawDir.nsamp);
290 else
291 one.resize(sel.cols(),thisRawDir.nsamp);
292
293 one.setZero();
294 }
295 else
296 {
297 fid->read_tag(t_pTag, thisRawDir.ent->pos);
298 //
299 // Depending on the state of the projection and selection
300 // we proceed a little bit differently
301 //
302 if (mult.cols() == 0)
303 {
304 if (sel.cols() == 0)
305 {
306 if (t_pTag->type == FIFFT_DAU_PACK16)
307 one = cal*(Map< MatrixDau16 >( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
308 else if(t_pTag->type == FIFFT_INT)
309 one = cal*(Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
310 else if(t_pTag->type == FIFFT_FLOAT)
311 one = cal*(Map< MatrixXf >( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
312 else if(t_pTag->type == FIFFT_SHORT)
313 one = cal*(Map< MatrixShort >( t_pTag->toShort(),nchan, thisRawDir.nsamp)).cast<double>();
314 else
315 printf("Data Storage Format not known yet [1]!! Type: %d\n", t_pTag->type);
316 }
317 else
318 {
319 //ToDo find a faster solution for this!! --> make cal and mul sparse like in MATLAB
320 newData.resize(sel.cols(), thisRawDir.nsamp); //ToDo this can be done much faster, without newData
321
322 if (t_pTag->type == FIFFT_DAU_PACK16)
323 {
324 tmp_data = (Map< MatrixDau16 > ( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
325
326 for(r = 0; r < sel.size(); ++r)
327 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
328 }
329 else if(t_pTag->type == FIFFT_INT)
330 {
331 tmp_data = (Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
332
333 for(r = 0; r < sel.size(); ++r)
334 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
335 }
336 else if(t_pTag->type == FIFFT_FLOAT)
337 {
338 tmp_data = (Map< MatrixXf > ( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
339
340 for(r = 0; r < sel.size(); ++r)
341 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
342 }
343 else if(t_pTag->type == FIFFT_SHORT)
344 {
345 tmp_data = (Map< MatrixShort > ( t_pTag->toShort(),nchan, thisRawDir.nsamp)).cast<double>();
346
347 for(r = 0; r < sel.size(); ++r)
348 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
349 }
350 else
351 {
352 printf("Data Storage Format not known yet [2]!! Type: %d\n", t_pTag->type);
353 }
354
355 one = cal*newData;
356 }
357 }
358 else
359 {
360 if (t_pTag->type == FIFFT_DAU_PACK16)
361 one = mult*(Map< MatrixDau16 >( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
362 else if(t_pTag->type == FIFFT_INT)
363 one = mult*(Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
364 else if(t_pTag->type == FIFFT_FLOAT)
365 one = mult*(Map< MatrixXf >( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
366 else
367 printf("Data Storage Format not known yet [3]!! Type: %d\n", t_pTag->type);
368 }
369 }
370 //
371 // The picking logic is a bit complicated
372 //
373 if (to >= thisRawDir.last && from <= thisRawDir.first)
374 {
375 //
376 // We need the whole buffer
377 //
378 first_pick = 0;//1;
379 last_pick = thisRawDir.nsamp - 1;
380 if (do_debug)
381 printf("W");
382 }
383 else if (from > thisRawDir.first)
384 {
385 first_pick = from - thisRawDir.first;// + 1;
386 if(to < thisRawDir.last)
387 {
388 //
389 // Something from the middle
390 //
391// qDebug() << "This needs to be debugged!";
392 last_pick = thisRawDir.nsamp + to - thisRawDir.last - 1;//is this alright?
393 if (do_debug)
394 printf("M");
395 }
396 else
397 {
398 //
399 // From the middle to the end
400 //
401 last_pick = thisRawDir.nsamp - 1;
402 if (do_debug)
403 printf("E");
404 }
405 }
406 else
407 {
408 //
409 // From the beginning to the middle
410 //
411 first_pick = 0;//1;
412 last_pick = to - thisRawDir.first;// + 1;
413 if (do_debug)
414 printf("B");
415 }
416 //
417 // Now we are ready to pick
418 //
419 picksamp = last_pick - first_pick + 1;
420
421 if(do_debug)
422 {
423 qDebug() << "first_pick: " << first_pick;
424 qDebug() << "last_pick: " << last_pick;
425 qDebug() << "picksamp: " << picksamp;
426 }
427
428 if (picksamp > 0)
429 {
430// for(r = 0; r < data->rows(); ++r)
431// for(c = 0; c < picksamp; ++c)
432// (*data)(r,dest + c) = one(r,first_pick + c);
433 data.block(0,dest,data.rows(),picksamp) = one.block(0, first_pick, data.rows(), picksamp);
434
435 dest += picksamp;
436 }
437 }
438 //
439 // Done?
440 //
441 if (thisRawDir.last >= to)
442 {
443 //printf(" [done]\n");
444 break;
445 }
446 }
447
448 if (!this->file->device()->isOpen()) {
449 this->file->device()->close();
450 }
451
452 times = MatrixXd(1, to-from+1);
453
454 for (i = 0; i < times.cols(); ++i)
455 times(0, i) = ((float)(from+i)) / this->info.sfreq;
456
457 return true;
458}
459
460//=============================================================================================================
461
462bool FiffRawData::read_raw_segment(MatrixXd& data,
463 MatrixXd& times,
464 SparseMatrix<double>& multSegment,
465 fiff_int_t from,
466 fiff_int_t to,
467 const RowVectorXi& sel,
468 bool do_debug) const
469{
470 bool projAvailable = true;
471
472 if (this->proj.size() == 0) {
473 //qInfo() << "FiffRawData::read_raw_segment - No projectors setup. Consider calling MNE::setup_compensators.";
474 projAvailable = false;
475 }
476
477 if(from == -1)
478 from = this->first_samp;
479 if(to == -1)
480 to = this->last_samp;
481 //
482 // Initial checks
483 //
484 if(from < this->first_samp)
485 from = this->first_samp;
486 if(to > this->last_samp)
487 to = this->last_samp;
488 //
489 if(from > to)
490 {
491 printf("No data in this range\n");
492 return false;
493 }
494 //printf("Reading %d ... %d = %9.3f ... %9.3f secs...", from, to, ((float)from)/this->info.sfreq, ((float)to)/this->info.sfreq);
495 //
496 // Initialize the data and calibration vector
497 //
498 qint32 nchan = this->info.nchan;
499 qint32 dest = 0;//1;
500 qint32 i, k, r;
501
502 typedef Eigen::Triplet<double> T;
503 std::vector<T> tripletList;
504 tripletList.reserve(nchan);
505 for(i = 0; i < nchan; ++i)
506 tripletList.push_back(T(i, i, this->cals[i]));
507
508 SparseMatrix<double> cal(nchan, nchan);
509 cal.setFromTriplets(tripletList.begin(), tripletList.end());
510// cal.makeCompressed();
511
512 MatrixXd mult_full;
513 //
514 if (sel.size() == 0)
515 {
516 data = MatrixXd(nchan, to-from+1);
517// data->setZero();
518 if (projAvailable || this->comp.kind != -1)
519 {
520 if (!projAvailable)
521 mult_full = this->comp.data->data*cal;
522 else if (this->comp.kind == -1)
523 mult_full = this->proj*cal;
524 else
525 mult_full = this->proj*this->comp.data->data*cal;
526 }
527 }
528 else
529 {
530 data = MatrixXd(sel.size(),to-from+1);
531// data->setZero();
532
533 MatrixXd selVect(sel.size(), nchan);
534
535 selVect.setZero();
536
537 if (!projAvailable && this->comp.kind == -1)
538 {
539 tripletList.clear();
540 tripletList.reserve(sel.size());
541 for(i = 0; i < sel.size(); ++i)
542 tripletList.push_back(T(i, i, this->cals[sel[i]]));
543 cal = SparseMatrix<double>(sel.size(), sel.size());
544 cal.setFromTriplets(tripletList.begin(), tripletList.end());
545 }
546 else
547 {
548 if (!projAvailable)
549 {
550 qDebug() << "This has to be debugged! #1";
551 for( i = 0; i < sel.size(); ++i)
552 selVect.row(i) = this->comp.data->data.block(sel[i],0,1,nchan);
553 mult_full = selVect*cal;
554 }
555 else if (this->comp.kind == -1)
556 {
557 for( i = 0; i < sel.size(); ++i)
558 selVect.row(i) = this->proj.block(sel[i],0,1,nchan);
559
560 mult_full = selVect*cal;
561 }
562 else
563 {
564 qDebug() << "This has to be debugged! #3";
565 for( i = 0; i < sel.size(); ++i)
566 selVect.row(i) = this->proj.block(sel[i],0,1,nchan);
567
568 mult_full = selVect*this->comp.data->data*cal;
569 }
570 }
571 }
572
573 //
574 // Make mult sparse
575 //
576 tripletList.clear();
577 tripletList.reserve(mult_full.rows()*mult_full.cols());
578 for(i = 0; i < mult_full.rows(); ++i)
579 for(k = 0; k < mult_full.cols(); ++k)
580 if(mult_full(i,k) != 0)
581 tripletList.push_back(T(i, k, mult_full(i,k)));
582
583 SparseMatrix<double> mult(mult_full.rows(),mult_full.cols());
584 if(tripletList.size() > 0)
585 mult.setFromTriplets(tripletList.begin(), tripletList.end());
586// mult.makeCompressed();
587
588 //
589
591 if (!this->file->device()->isOpen())
592 {
593 if (!this->file->device()->open(QIODevice::ReadOnly))
594 {
595 printf("Cannot open file %s",this->info.filename.toUtf8().constData());
596 }
597 fid = this->file;
598 }
599 else
600 {
601 fid = this->file;
602 }
603
604 MatrixXd one;
605 fiff_int_t first_pick, last_pick, picksamp;
606 for(k = 0; k < this->rawdir.size(); ++k)
607 {
608 FiffRawDir thisRawDir = this->rawdir[k];
609 //
610 // Do we need this buffer
611 //
612 if (thisRawDir.last > from)
613 {
614 if (thisRawDir.ent->kind == -1)
615 {
616 //
617 // Take the easy route: skip is translated to zeros
618 //
619 if(do_debug)
620 printf("S");
621 if (sel.cols() <= 0)
622 one.resize(nchan,thisRawDir.nsamp);
623 else
624 one.resize(sel.cols(),thisRawDir.nsamp);
625
626 one.setZero();
627 }
628 else
629 {
630 FiffTag::SPtr t_pTag;
631 fid->read_tag(t_pTag, thisRawDir.ent->pos);
632 //
633 // Depending on the state of the projection and selection
634 // we proceed a little bit differently
635 //
636 if (mult.cols() == 0)
637 {
638 if (sel.cols() == 0)
639 {
640 if (t_pTag->type == FIFFT_DAU_PACK16)
641 one = cal*(Map< MatrixDau16 >( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
642 else if(t_pTag->type == FIFFT_INT)
643 one = cal*(Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
644 else if(t_pTag->type == FIFFT_FLOAT)
645 one = cal*(Map< MatrixXf >( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
646 else if(t_pTag->type == FIFFT_SHORT)
647 one = cal*(Map< MatrixShort >( t_pTag->toShort(),nchan, thisRawDir.nsamp)).cast<double>();
648 else
649 printf("Data Storage Format not known yet [1]!! Type: %d\n", t_pTag->type);
650 }
651 else
652 {
653
654 //ToDo find a faster solution for this!! --> make cal and mul sparse like in MATLAB
655 MatrixXd newData(sel.cols(), thisRawDir.nsamp); //ToDo this can be done much faster, without newData
656
657 if (t_pTag->type == FIFFT_DAU_PACK16)
658 {
659 MatrixXd tmp_data = (Map< MatrixDau16 > ( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
660
661 for(r = 0; r < sel.size(); ++r)
662 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
663 }
664 else if(t_pTag->type == FIFFT_INT)
665 {
666 MatrixXd tmp_data = (Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
667
668 for(r = 0; r < sel.size(); ++r)
669 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
670 }
671 else if(t_pTag->type == FIFFT_FLOAT)
672 {
673 MatrixXd tmp_data = (Map< MatrixXf > ( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
674
675 for(r = 0; r < sel.size(); ++r)
676 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
677 }
678 else if(t_pTag->type == FIFFT_SHORT)
679 {
680 MatrixXd tmp_data = (Map< MatrixShort > ( t_pTag->toShort(),nchan, thisRawDir.nsamp)).cast<double>();
681
682 for(r = 0; r < sel.size(); ++r)
683 newData.block(r,0,1,thisRawDir.nsamp) = tmp_data.block(sel[r],0,1,thisRawDir.nsamp);
684 }
685 else
686 {
687 printf("Data Storage Format not known yet [2]!! Type: %d\n", t_pTag->type);
688 }
689
690 one = cal*newData;
691 }
692 }
693 else
694 {
695 if (t_pTag->type == FIFFT_DAU_PACK16)
696 one = mult*(Map< MatrixDau16 >( t_pTag->toDauPack16(),nchan, thisRawDir.nsamp)).cast<double>();
697 else if(t_pTag->type == FIFFT_INT)
698 one = mult*(Map< MatrixXi >( t_pTag->toInt(),nchan, thisRawDir.nsamp)).cast<double>();
699 else if(t_pTag->type == FIFFT_FLOAT)
700 one = mult*(Map< MatrixXf >( t_pTag->toFloat(),nchan, thisRawDir.nsamp)).cast<double>();
701 else
702 printf("Data Storage Format not known yet [3]!! Type: %d\n", t_pTag->type);
703 }
704 }
705 //
706 // The picking logic is a bit complicated
707 //
708 if (to >= thisRawDir.last && from <= thisRawDir.first)
709 {
710 //
711 // We need the whole buffer
712 //
713 first_pick = 0;//1;
714 last_pick = thisRawDir.nsamp - 1;
715 if (do_debug)
716 printf("W");
717 }
718 else if (from > thisRawDir.first)
719 {
720 first_pick = from - thisRawDir.first;// + 1;
721 if(to < thisRawDir.last)
722 {
723 //
724 // Something from the middle
725 //
726// qDebug() << "This needs to be debugged!";
727 last_pick = thisRawDir.nsamp + to - thisRawDir.last - 1;//is this alright?
728 if (do_debug)
729 printf("M");
730 }
731 else
732 {
733 //
734 // From the middle to the end
735 //
736 last_pick = thisRawDir.nsamp - 1;
737 if (do_debug)
738 printf("E");
739 }
740 }
741 else
742 {
743 //
744 // From the beginning to the middle
745 //
746 first_pick = 0;//1;
747 last_pick = to - thisRawDir.first;// + 1;
748 if (do_debug)
749 printf("B");
750 }
751 //
752 // Now we are ready to pick
753 //
754 picksamp = last_pick - first_pick + 1;
755
756 if(do_debug)
757 {
758 qDebug() << "first_pick: " << first_pick;
759 qDebug() << "last_pick: " << last_pick;
760 qDebug() << "picksamp: " << picksamp;
761 }
762
763 if (picksamp > 0)
764 {
765// for(r = 0; r < data->rows(); ++r)
766// for(c = 0; c < picksamp; ++c)
767// (*data)(r,dest + c) = one(r,first_pick + c);
768 data.block(0,dest,data.rows(),picksamp) = one.block(0, first_pick, data.rows(), picksamp);
769
770 dest += picksamp;
771 }
772 }
773 //
774 // Done?
775 //
776 if (thisRawDir.last >= to)
777 {
778 //printf(" [done]\n");
779 break;
780 }
781 }
782
783 if(mult.cols()==0)
784 multSegment = cal;
785 else
786 multSegment = mult;
787
788 if (!this->file->device()->isOpen()) {
789 this->file->device()->close();
790 }
791
792 times = MatrixXd(1, to-from+1);
793
794 for (i = 0; i < times.cols(); ++i)
795 times(0, i) = ((float)(from+i)) / this->info.sfreq;
796
797 return true;
798}
799
800//=============================================================================================================
801
803 MatrixXd& times,
804 float from,
805 float to,
806 const RowVectorXi& sel) const
807{
808 //
809 // Convert to samples
810 //
811 from = floor(from*this->info.sfreq);
812 to = ceil(to*this->info.sfreq);
813 //
814 // Read it
815 //
816 return this->read_raw_segment(data, times, (qint32)from, (qint32)to, sel);
817}
818
819//=============================================================================================================
820
821bool FiffRawData::save(QIODevice &p_IODevice,
822 const RowVectorXi &picks,
823 int decim,
824 int from,
825 int to) const
826{
827 if (decim < 1) decim = 1;
828
829 int firstSamp = (from >= 0) ? from : first_samp;
830 int lastSamp = (to >= 0) ? to : last_samp;
831
832 if (firstSamp > lastSamp) {
833 qWarning() << "[FiffRawData::save] Invalid sample range.";
834 return false;
835 }
836
837 // Prepare output info
838 FiffInfo outInfo;
839 if (picks.size() > 0) {
840 outInfo = info.pick_info(picks);
841 } else {
842 outInfo = info;
843 }
844
845 // Adjust sampling frequency for decimation
846 if (decim > 1) {
847 outInfo.sfreq = info.sfreq / static_cast<float>(decim);
848 }
849
850 // Use the standard start_writing_raw pipeline
851 RowVectorXd calsOut;
852 FiffStream::SPtr pStream = FiffStream::start_writing_raw(p_IODevice, outInfo, calsOut, picks);
853 if (!pStream) {
854 qWarning() << "[FiffRawData::save] Cannot start writing raw file.";
855 return false;
856 }
857
858 // Write data in blocks
859 const int blockSize = 2000;
860 int blockSamples = decim * blockSize;
861
862 for (int samp = firstSamp; samp <= lastSamp; samp += blockSamples) {
863 int nsamp = qMin(blockSamples, lastSamp - samp + 1);
864
865 MatrixXd segData;
866 MatrixXd segTimes;
867 if (!read_raw_segment(segData, segTimes, samp, samp + nsamp - 1, picks)) {
868 qWarning() << "[FiffRawData::save] Error reading data at sample" << samp;
869 pStream->finish_writing_raw();
870 return false;
871 }
872
873 // Decimate if needed
874 if (decim > 1) {
875 int nOut = (nsamp + decim - 1) / decim;
876 MatrixXd decimData(segData.rows(), nOut);
877 for (int s = 0, idx = 0; s < nsamp && idx < nOut; s += decim, ++idx) {
878 decimData.col(idx) = segData.col(s);
879 }
880 segData = decimData;
881 }
882
883 pStream->write_raw_buffer(segData, calsOut);
884 }
885
886 pStream->finish_writing_raw();
887
888 qInfo() << "[FiffRawData::save] Saved raw data from sample" << firstSamp
889 << "to" << lastSamp << "(decim=" << decim << ")";
890 return true;
891}
FiffTag class declaration, which provides fiff tag I/O and processing methods.
FiffStream class declaration.
FiffRawData class declaration.
#define FIFFT_INT
Definition fiff_file.h:231
#define FIFFT_SHORT
Definition fiff_file.h:230
#define FIFFT_DAU_PACK16
Definition fiff_file.h:243
#define FIFFT_FLOAT
Definition fiff_file.h:232
FiffEvents class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
FIFF measurement file information.
Definition fiff_info.h:85
Eigen::RowVectorXd cals
bool read_raw_segment(Eigen::MatrixXd &data, Eigen::MatrixXd &times, fiff_int_t from=-1, fiff_int_t to=-1, const Eigen::RowVectorXi &sel=defaultRowVectorXi, bool do_debug=false) const
FiffStream::SPtr file
bool save(QIODevice &p_IODevice, const Eigen::RowVectorXi &picks=Eigen::RowVectorXi(), int decim=1, int from=-1, int to=-1) const
Eigen::MatrixXd proj
bool read_raw_segment_times(Eigen::MatrixXd &data, Eigen::MatrixXd &times, float from, float to, const Eigen::RowVectorXi &sel=defaultRowVectorXi) const
QList< FiffRawDir > rawdir
FiffDirEntry::SPtr ent
static bool setup_read_raw(QIODevice &p_IODevice, FiffRawData &data, bool allow_maxshield=true, bool is_littleEndian=false)
static FiffStream::SPtr start_writing_raw(QIODevice &p_IODevice, const FiffInfo &info, Eigen::RowVectorXd &cals, Eigen::MatrixXi sel=defaultMatrixXi, bool bResetRange=true)
QSharedPointer< FiffStream > SPtr
QSharedPointer< FiffTag > SPtr
Definition fiff_tag.h:155