MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
sensorset.cpp
Go to the documentation of this file.
1//=============================================================================================================
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
40#include <iostream>
41#include <fwd/fwd_coil_set.h>
42#include <fiff/fiff_ch_info.h>
43
44//=============================================================================================================
45// QT INCLUDES
46//=============================================================================================================
47
48#include <QCoreApplication>
49
50//=============================================================================================================
51// EIGEN INCLUDES
52//=============================================================================================================
53
54#include <Eigen/Dense>
55
56//=============================================================================================================
57// USED NAMESPACES
58//=============================================================================================================
59
60using namespace INVERSELIB;
61using namespace FIFFLIB;
62using namespace FWDLIB;
63using namespace Eigen;
64
65//=============================================================================================================
66// DEFINE GLOBAL METHODS
67//=============================================================================================================
68
69//=============================================================================================================
70// DEFINE MEMBER METHODS
71//=============================================================================================================
72
73SensorSet::SensorSet(const FwdCoilSet::SPtr pFwdSensorSet)
74{
75 if(pFwdSensorSet!=nullptr) {
76 initFromFwdCoilSet(pFwdSensorSet);
77 }
78}
79
80//=============================================================================================================
81
82void SensorSet::initFromFwdCoilSet(const QSharedPointer<FWDLIB::FwdCoilSet> pFwdSensorSet)
83{
84 m_ncoils = pFwdSensorSet->ncoil;
85 m_np = pFwdSensorSet->coils[0]->np;
86 initMatrices(m_ncoils,m_np);
87
88 // get data froms Fwd Coilset
89 for(int i = 0; i < m_ncoils; i++){
90 const FwdCoil* coil = (pFwdSensorSet->coils[i]);
91 MatrixXd matRmag = MatrixXd::Zero(m_np,3);
92 MatrixXd matCosmag = MatrixXd::Zero(m_np,3);
93 RowVectorXd vecW(m_np);
94
95 m_r0(i,0) = coil->r0[0];
96 m_r0(i,1) = coil->r0[1];
97 m_r0(i,2) = coil->r0[2];
98
99 m_ez(i,0) = coil->ez[0];
100 m_ez(i,1) = coil->ez[1];
101 m_ez(i,2) = coil->ez[2];
102
103 for (int p = 0; p < m_np; p++){
104 m_w(i*m_np+p) = coil->w[p];
105 for (int c = 0; c < 3; c++) {
106 matRmag(p,c) = coil->rmag[p][c];
107 matCosmag(p,c) = coil->cosmag[p][c];
108 }
109 }
110
111 m_cosmag.block(i*m_np,0,m_np,3) = matCosmag;
112 m_rmag.block(i*m_np,0,m_np,3) = matRmag;
113 }
114 m_tra = MatrixXd::Identity(m_ncoils,m_ncoils);
115}
116
117//=============================================================================================================
118
119void SensorSet::initMatrices(int ncoils, int np)
120{
121 m_ez = MatrixXd(ncoils,3);
122 m_r0 = MatrixXd(ncoils,3);
123 m_rmag = MatrixXd(ncoils*np,3);
124 m_cosmag = MatrixXd(ncoils*np,3);
125 m_cosmag = MatrixXd(ncoils*np,3);
126 m_tra = MatrixXd(ncoils,ncoils);
127 m_w = RowVectorXd(ncoils*np);
128}
129
130//=============================================================================================================
131
133{
134 const QString qPath = QString(QCoreApplication::applicationDirPath() + "/../resources/general/coilDefinitions/coil_def.dat");
135 m_pCoilDefinitions = FwdCoilSet::SPtr(FwdCoilSet::read_coil_defs(qPath));
136}
137
138//=============================================================================================================
139
140SensorSet SensorSetCreator::updateSensorSet(const QList<FIFFLIB::FiffChInfo>& channelList,
141 const Accuracy& accuracy)
142{
143 if(channelList.isEmpty()) {
144 return SensorSet();
145 } else {
146 auto pCoilMeg = FwdCoilSet::SPtr(m_pCoilDefinitions->create_meg_coils(channelList, channelList.size(), static_cast<int>(accuracy), nullptr));
147 return SensorSet(pCoilMeg);
148 }
149}
150
151//=============================================================================================================
FiffChInfo class declaration.
FwdCoilSet class declaration.
SensorSet class declaration.
FwdCoil description.
Definition fwd_coil.h:89
float ** rmag
Definition fwd_coil.h:180
float ez[3]
Definition fwd_coil.h:178
float ** cosmag
Definition fwd_coil.h:181
float r0[3]
Definition fwd_coil.h:175
static FwdCoilSet * read_coil_defs(const QString &name)
QSharedPointer< FwdCoilSet > SPtr
SensorSet updateSensorSet(const QList< FIFFLIB::FiffChInfo > &channelList, const Accuracy &accuracy)