v2.0.0
Loading...
Searching...
No Matches
fiff_info.cpp
Go to the documentation of this file.
1//=============================================================================================================
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
42#include "fiff_info.h"
43#include "fiff_stream.h"
44#include "fiff_file.h"
45
46#include <utils/ioutils.h>
47#include <iostream>
48#include <QDebug>
49
50//=============================================================================================================
51// USED NAMESPACES
52//=============================================================================================================
53
54using namespace FIFFLIB;
55using namespace UTILSLIB;
56using namespace Eigen;
57
58//=============================================================================================================
59// DEFINE MEMBER METHODS
60//=============================================================================================================
61
63: FiffInfoBase()//nchan(-1)
64, sfreq(-1.0)
65, linefreq(-1.0)
66, highpass(-1.0)
67, lowpass(-1.0)
68, gantry_angle(-1)
69, acq_pars("")
70, acq_stim("")
71{
72 meas_date[0] = -1;
73}
74
75//=============================================================================================================
76
77FiffInfo::FiffInfo(const FiffInfo& p_FiffInfo)
78: FiffInfoBase(p_FiffInfo)
79, file_id(p_FiffInfo.file_id)
80, sfreq(p_FiffInfo.sfreq)
81, linefreq(p_FiffInfo.linefreq)
82, highpass(p_FiffInfo.highpass)
83, lowpass(p_FiffInfo.lowpass)
84, proj_id(p_FiffInfo.proj_id)
85, proj_name(p_FiffInfo.proj_name)
87, experimenter(p_FiffInfo.experimenter)
88, description(p_FiffInfo.description)
89, utc_offset(p_FiffInfo.utc_offset)
90, gantry_angle(p_FiffInfo.gantry_angle)
91, dev_ctf_t(p_FiffInfo.dev_ctf_t)
92, dig(p_FiffInfo.dig)
93, dig_trans(p_FiffInfo.dig_trans)
94, projs(p_FiffInfo.projs)
95, comps(p_FiffInfo.comps)
96, acq_pars(p_FiffInfo.acq_pars)
97, acq_stim(p_FiffInfo.acq_stim)
98{
99 meas_date[0] = p_FiffInfo.meas_date[0];
100 meas_date[1] = p_FiffInfo.meas_date[1];
101}
102
103//=============================================================================================================
104
108
109//=============================================================================================================
110
112{
114 file_id = FiffId();
115 meas_date[0] = -1;
116 sfreq = -1.0;
117 linefreq = -1.0;
118 highpass = -1.0;
119 lowpass = -1.0;
120 proj_id = -1;
121 proj_name = "";
122 xplotter_layout = "";
123 experimenter = "";
124 description = "";
125 utc_offset = "";
126 gantry_angle = -1;
127 dev_ctf_t.clear();
128 dig.clear();
129 dig_trans.clear();
130 projs.clear();
131 comps.clear();
132 acq_pars = "";
133 acq_stim = "";
134}
135
136//=============================================================================================================
137
139{
140 qint32 comp = 0;
141 qint32 first_comp = -1;
142
143 qint32 k = 0;
144 for (k = 0; k < this->nchan; ++k)
145 {
146 if (this->chs[k].kind == FIFFV_MEG_CH)
147 {
148 comp = this->chs[k].chpos.coil_type >> 16;
149 if (first_comp < 0)
150 first_comp = comp;
151 else if (comp != first_comp)
152 qWarning("Compensation is not set equally on all MEG channels");
153 }
154 }
155 return comp;
156}
157
158//=============================================================================================================
159
160bool FiffInfo::make_compensator(fiff_int_t from, fiff_int_t to, FiffCtfComp& ctf_comp, bool exclude_comp_chs) const
161{
162 MatrixXd C1, C2, comp_tmp;
163
164// if(ctf_comp.data)
165// delete ctf_comp.data;
166 ctf_comp.data->clear();
167
168 if (from == to)
169 {
170 ctf_comp.data->data = MatrixXd::Identity(this->nchan, this->nchan);
171 return false;
172 }
173
174 if (from == 0)
175 C1 = MatrixXd::Zero(this->nchan,this->nchan);
176 else
177 {
178 if (!this->make_compensator(from, C1))
179 {
180 qWarning("Cannot create compensator C1\n");
181 qWarning("Desired compensation matrix (kind = %d) not found\n",from);
182 return false;
183 }
184 }
185
186 if (to == 0)
187 C2 = MatrixXd::Zero(this->nchan,this->nchan);
188 else
189 {
190 if (!this->make_compensator(to, C2))
191 {
192 qWarning("Cannot create compensator C2\n");
193 qWarning("Desired compensation matrix (kind = %d) not found\n",to);
194 return false;
195 }
196 }
197 //
198 // s_orig = s_from + C1*s_from = (I + C1)*s_from
199 // s_to = s_orig - C2*s_orig = (I - C2)*s_orig
200 // s_to = (I - C2)*(I + C1)*s_from = (I + C1 - C2 - C2*C1)*s_from
201 //
202 comp_tmp = MatrixXd::Identity(this->nchan,this->nchan) + C1 - C2 - C2*C1;
203
204 qint32 k;
205 if (exclude_comp_chs)
206 {
207 VectorXi pick = MatrixXi::Zero(1,this->nchan);
208 qint32 npick = 0;
209 for (k = 0; k < this->nchan; ++k)
210 {
211 if (this->chs[k].kind != FIFFV_REF_MEG_CH)
212 {
213 pick(npick) = k;
214 ++npick;
215 }
216 }
217 if (npick == 0)
218 {
219 qWarning("Nothing remains after excluding the compensation channels\n");
220 return false;
221 }
222
223 ctf_comp.data->data.resize(npick,this->nchan);
224 for (k = 0; k < npick; ++k)
225 ctf_comp.data->data.row(k) = comp_tmp.block(pick(k), 0, 1, this->nchan);
226 } else {
227 ctf_comp.data->data = comp_tmp;
228 }
229
230 return true;
231}
232
233//=============================================================================================================
234
235bool FiffInfo::make_compensator(fiff_int_t kind, MatrixXd& this_comp) const//private method
236{
237 FiffNamedMatrix::SDPtr this_data;
238 MatrixXd presel, postsel;
239 qint32 k, col, c, ch=0, row, row_ch=0, channelAvailable;
240
241 for (k = 0; k < this->comps.size(); ++k)
242 {
243 if (this->comps[k].kind == kind)
244 {
245 this_data = this->comps[k].data;
246
247 //
248 // Create the preselector
249 //
250 presel = MatrixXd::Zero(this_data->ncol,this->nchan);
251
252 for(col = 0; col < this_data->ncol; ++col)
253 {
254 channelAvailable = 0;
255 for (c = 0; c < this->ch_names.size(); ++c)
256 {
257 if (QString::compare(this_data->col_names.at(col),this->ch_names.at(c)) == 0)
258 {
259 ++channelAvailable;
260 ch = c;
261 }
262 }
263 if (channelAvailable == 0)
264 {
265 qWarning("Channel %s is not available in data\n",this_data->col_names.at(col).toUtf8().constData());
266 return false;
267 }
268 else if (channelAvailable > 1)
269 {
270 qWarning("Ambiguous channel %s",this_data->col_names.at(col).toUtf8().constData());
271 return false;
272 }
273 presel(col,ch) = 1.0;
274 }
275 //
276 // Create the postselector
277 //
278 postsel = MatrixXd::Zero(this->nchan,this_data->nrow);
279
280 for (c = 0; c < this->nchan; ++c)
281 {
282 channelAvailable = 0;
283 for (row = 0; row < this_data->row_names.size(); ++row)
284 {
285 if (QString::compare(this->ch_names.at(c),this_data->row_names.at(row)) == 0)
286 {
287 ++channelAvailable;
288 row_ch = row;
289 }
290 }
291 if (channelAvailable > 1)
292 {
293 qWarning("Ambiguous channel %s", this->ch_names.at(c).toUtf8().constData());
294 return false;
295 }
296 else if (channelAvailable == 1)
297 {
298 postsel(c,row_ch) = 1.0;
299 }
300 }
301 this_comp = postsel*this_data->data*presel;
302 return true;
303 }
304 }
305 this_comp = defaultMatrixXd;
306 return false;
307}
308
309//=============================================================================================================
310
311FiffInfo FiffInfo::pick_info(const RowVectorXi &sel) const
312{
313 FiffInfo res = *this;//new FiffInfo(this);
314 if (sel.size() == 0)
315 return res;
316
317 //ToDo when pointer List do delation
318 res.chs.clear();
319 res.ch_names.clear();
320
321 qint32 idx;
322 for(qint32 i = 0; i < sel.size(); ++i)
323 {
324 idx = sel[i];
325 res.chs.append(this->chs[idx]);
326 res.ch_names.append(this->ch_names[idx]);
327 }
328 res.nchan = static_cast<int>(sel.size());
329
330 return res;
331}
332
333//=============================================================================================================
334
335QList<FiffChInfo> FiffInfo::set_current_comp(QList<FiffChInfo>& listFiffChInfo, fiff_int_t value)
336{
337 QList<FiffChInfo> newList;
338 qint32 k;
339 fiff_int_t coil_type;
340
341 for(k = 0; k < listFiffChInfo.size(); ++k)
342 newList.append(listFiffChInfo[k]);
343
344 qint32 lower_half = 65535;// hex2dec('FFFF');
345 for (k = 0; k < listFiffChInfo.size(); ++k)
346 {
347 if (listFiffChInfo[k].kind == FIFFV_MEG_CH)
348 {
349 coil_type = listFiffChInfo[k].chpos.coil_type & lower_half;
350 newList[k].chpos.coil_type = (coil_type | (value << 16));
351 }
352 }
353 return newList;
354}
355
356//=============================================================================================================
357
358bool FiffInfo::readMegEegChannels(const QString& name,
359 bool do_meg,
360 bool do_eeg,
361 const QStringList& bads,
362 QList<FiffChInfo>& chsp,
363 int& nmegp,
364 int& neegp)
365{
366 QFile file(name);
367 FiffStream::SPtr stream(new FiffStream(&file));
368
369 if (!stream->open())
370 return false;
371
372 FiffInfo info;
373 FiffDirNode::SPtr infoNode;
374 if (!stream->read_meas_info(stream->dirtree(), info, infoNode)) {
375 qCritical("%s : could not read measurement info", name.toUtf8().data());
376 stream->close();
377 return false;
378 }
379 stream->close();
380
381 QList<FiffChInfo> meg;
382 QList<FiffChInfo> eeg;
383
384 for (int k = 0; k < info.chs.size(); k++) {
385 if (bads.contains(info.chs[k].ch_name))
386 continue;
387 if (do_meg && info.chs[k].kind == FIFFV_MEG_CH)
388 meg.append(info.chs[k]);
389 else if (do_eeg && info.chs[k].isValidEeg())
390 eeg.append(info.chs[k]);
391 }
392
393 chsp.clear();
394 chsp.reserve(meg.size() + eeg.size());
395 chsp.append(meg);
396 chsp.append(eeg);
397
398 nmegp = meg.size();
399 neegp = eeg.size();
400 return true;
401}
402
403//=============================================================================================================
404
406{
407 //
408 // We will always write floats
409 //
410 fiff_int_t data_type = 4;
411 qint32 k;
412 QList<FiffChInfo> chs;
413
414 for(k = 0; k < this->nchan; ++k)
415 chs << this->chs[k];
416
417 fiff_int_t nchan = chs.size();
418
419 //
420 // write the essentials
421 //
422 p_pStream->start_block(FIFFB_MEAS);//4
423 p_pStream->write_id(FIFF_BLOCK_ID);//5
424 if(this->meas_id.version != -1)
425 {
426 p_pStream->write_id(FIFF_PARENT_BLOCK_ID,this->meas_id);//6
427 }
428 //
429 // Measurement info
430 //
431 p_pStream->start_block(FIFFB_MEAS_INFO);//7
432
433 //
434 // Blocks from the original -> skip this
435 //
436// QList<fiff_int_t> blocks;
437// blocks << FIFFB_SUBJECT << FIFFB_HPI_MEAS << FIFFB_HPI_RESULT << FIFFB_ISOTRAK << FIFFB_PROCESSING_HISTORY;
438 bool have_hpi_result = false;
439 bool have_isotrak = false;
440 //
441 // megacq parameters
442 //
443 if (!this->acq_pars.isEmpty() || !this->acq_stim.isEmpty())
444 {
445 p_pStream->start_block(FIFFB_DACQ_PARS);
446 if (!this->acq_pars.isEmpty())
447 p_pStream->write_string(FIFF_DACQ_PARS, this->acq_pars);
448
449 if (!this->acq_stim.isEmpty())
450 p_pStream->write_string(FIFF_DACQ_STIM, this->acq_stim);
451
452 p_pStream->end_block(FIFFB_DACQ_PARS);
453 }
454 //
455 // Coordinate transformations if the HPI result block was not there
456 //
457 if (!have_hpi_result)
458 {
459 if (!this->dev_head_t.isEmpty())
460 p_pStream->write_coord_trans(this->dev_head_t);
461
462 if (!this->ctf_head_t.isEmpty())
463 p_pStream->write_coord_trans(this->ctf_head_t);
464 }
465 //
466 // Polhemus data
467 //
468 if (this->dig.size() > 0 && !have_isotrak)
469 {
470 p_pStream->start_block(FIFFB_ISOTRAK);
471 for (qint32 k = 0; k < this->dig.size(); ++k)
472 p_pStream->write_dig_point(this->dig[k]);
473
474 p_pStream->end_block(FIFFB_ISOTRAK);
475 }
476 //
477 // Projectors
478 //
479 p_pStream->write_proj(this->projs);
480 //
481 // CTF compensation info
482 //
483 p_pStream->write_ctf_comp(this->comps);
484 //
485 // Bad channels
486 //
487 if (this->bads.size() > 0)
488 {
492 }
493 //
494 // General
495 //
496 p_pStream->write_float(FIFF_SFREQ,&this->sfreq);
497 p_pStream->write_float(FIFF_LINE_FREQ,&this->linefreq);
498 p_pStream->write_float(FIFF_HIGHPASS,&this->highpass);
499 p_pStream->write_float(FIFF_LOWPASS,&this->lowpass);
501 p_pStream->write_string(FIFF_DESCRIPTION,this->description);
502 p_pStream->write_string(FIFF_UNIT_C,this->utc_offset);
503 p_pStream->write_string(FIFF_PROJ_NAME,this->proj_name);
504 p_pStream->write_int(FIFF_PROJ_ID,&this->proj_id);
505 p_pStream->write_int(FIFF_GANTRY_ANGLE,&this->gantry_angle);
506 p_pStream->write_int(FIFF_NCHAN,&nchan);
507 p_pStream->write_int(FIFF_DATA_PACK,&data_type);
508 if (this->meas_date[0] != -1)
509 p_pStream->write_int(FIFF_MEAS_DATE,this->meas_date, 2);
510 //
511 // Channel info
512 //
513 MatrixXd cals(1,nchan);
514
515 for(k = 0; k < nchan; ++k)
516 {
517 //
518 // Scan numbers may have been messed up
519 //
520 chs[k].scanNo = k+1;//+1 because
521// chs[k].range = 1.0f;//Why? -> cause its already calibrated through reading
522 cals(0,k) = static_cast<double>(chs[k].cal); //ToDo whats going on with cals?
523 p_pStream->write_ch_info(chs[k]);
524 }
525 //
526 //
527 p_pStream->end_block(FIFFB_MEAS_INFO);
528}
529
530//=============================================================================================================
531
532void FiffInfo::print() const
533{
534 std::cout << "Sample frequency: " << sfreq << "\n";
535 std::cout << "LineFreq: " << linefreq << " | Highpass: " << highpass << " | Lowpass: " << lowpass << "\n";
536 std::cout << "Number of digitizer points: " << dig.size() << "\n";
537 for (auto& point : dig){
538 if (point.kind == FIFFV_POINT_HPI){
539 std::cout << "HPI Point " << point.ident << " - " << point.r[0] << ", " << point.r[1] << ", " << point.r[2] << "\n";
540 }
541 }
542}
FiffInfo class declaration.
#define FIFF_UNIT_C
#define FIFFV_REF_MEG_CH
#define FIFFV_MEG_CH
#define FIFF_MNE_CH_NAME_LIST
#define FIFFB_MNE_BAD_CHANNELS
#define FIFFV_POINT_HPI
FiffStream class declaration.
Header file describing the numerical values used in fif files.
#define FIFF_GANTRY_ANGLE
Definition fiff_file.h:555
#define FIFF_DACQ_STIM
Definition fiff_file.h:351
#define FIFF_PARENT_BLOCK_ID
Definition fiff_file.h:333
#define FIFF_NCHAN
Definition fiff_file.h:453
#define FIFF_PROJ_ID
Definition fiff_file.h:575
#define FIFFB_ISOTRAK
Definition fiff_file.h:370
#define FIFF_HIGHPASS
Definition fiff_file.h:476
#define FIFF_EXPERIMENTER
Definition fiff_file.h:465
#define FIFF_PROJ_NAME
Definition fiff_file.h:576
#define FIFF_DATA_PACK
Definition fiff_file.h:455
#define FIFF_DESCRIPTION
Definition fiff_file.h:486
#define FIFF_LINE_FREQ
Definition fiff_file.h:489
#define FIFFB_MEAS
Definition fiff_file.h:362
#define FIFFB_DACQ_PARS
Definition fiff_file.h:379
#define FIFF_BLOCK_ID
Definition fiff_file.h:326
#define FIFF_DACQ_PARS
Definition fiff_file.h:350
#define FIFF_MEAS_DATE
Definition fiff_file.h:457
#define FIFF_LOWPASS
Definition fiff_file.h:472
#define FIFFB_MEAS_INFO
Definition fiff_file.h:363
#define FIFF_SFREQ
Definition fiff_file.h:454
IOUtils class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
CTF software compensation data.
QSharedPointer< FiffDirNode > SPtr
Universally unique identifier.
Definition fiff_id.h:68
FiffInfo pick_info(const Eigen::RowVectorXi &sel=defaultVectorXi) const
void set_current_comp(fiff_int_t value)
Definition fiff_info.h:299
void print() const
static bool readMegEegChannels(const QString &name, bool do_meg, bool do_eeg, const QStringList &bads, QList< FiffChInfo > &chsp, int &nmegp, int &neegp)
QString description
Definition fiff_info.h:269
FiffCoordTrans dev_ctf_t
Definition fiff_info.h:272
fiff_int_t gantry_angle
Definition fiff_info.h:271
void writeToStream(FiffStream *p_pStream) const
FiffCoordTrans dig_trans
Definition fiff_info.h:274
QList< FiffCtfComp > comps
Definition fiff_info.h:276
qint32 get_current_comp()
QString xplotter_layout
Definition fiff_info.h:267
fiff_int_t meas_date[2]
Definition fiff_info.h:260
bool make_compensator(fiff_int_t from, fiff_int_t to, FiffCtfComp &ctf_comp, bool exclude_comp_chs=false) const
QList< FiffDigPoint > dig
Definition fiff_info.h:273
QString experimenter
Definition fiff_info.h:268
QList< FiffProj > projs
Definition fiff_info.h:275
QList< FiffChInfo > chs
FiffCoordTrans ctf_head_t
FiffCoordTrans dev_head_t
QSharedDataPointer< FiffNamedMatrix > SDPtr
FIFF File I/O routines.
fiff_long_t start_block(fiff_int_t kind)
fiff_long_t write_proj(const QList< FiffProj > &projs)
QSharedPointer< FiffStream > SPtr
fiff_long_t write_dig_point(const FiffDigPoint &dig)
fiff_long_t write_int(fiff_int_t kind, const fiff_int_t *data, fiff_int_t nel=1, fiff_int_t next=FIFFV_NEXT_SEQ)
fiff_long_t write_float(fiff_int_t kind, const float *data, fiff_int_t nel=1)
fiff_long_t write_id(fiff_int_t kind, const FiffId &id=FiffId::getDefault())
fiff_long_t write_coord_trans(const FiffCoordTrans &trans)
fiff_long_t write_name_list(fiff_int_t kind, const QStringList &data)
fiff_long_t write_string(fiff_int_t kind, const QString &data)
fiff_long_t end_block(fiff_int_t kind, fiff_int_t next=FIFFV_NEXT_SEQ)
fiff_long_t write_ch_info(const FiffChInfo &ch)
fiff_long_t write_ctf_comp(const QList< FiffCtfComp > &comps)