v2.0.0
Loading...
Searching...
No Matches
sensorset.h
Go to the documentation of this file.
1//=============================================================================================================
34
35#ifndef SENSORSET_H
36#define SENSORSET_H
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
42#include "../inverse_global.h"
43
44//=============================================================================================================
45// QT INCLUDES
46//=============================================================================================================
47
48#include <QObject>
49#include <QSharedPointer>
50
51//=============================================================================================================
52// EIGEN INCLUDES
53//=============================================================================================================
54
55#include <Eigen/Core>
56
57//=============================================================================================================
58// FORWARD DECLARATIONS
59//=============================================================================================================
60
61namespace FWDLIB{
62 class FwdCoil;
63 class FwdCoilSet;
64}
65
66namespace FIFFLIB{
67 class FiffCoordTrans;
68 class FiffDigPointSet;
69 class FiffChInfo;
70}
71
72//=============================================================================================================
73// DEFINE NAMESPACE INVERSELIB
74//=============================================================================================================
75
76namespace INVERSELIB
77{
78enum class Accuracy : int{high = 2, medium = 1, low = 0};
79
84
85public:
86 typedef QSharedPointer<SensorSet> SPtr;
87 typedef QSharedPointer<const SensorSet> ConstSPtr;
88
89 //=========================================================================================================
93 explicit SensorSet() = default;
94
95 //=========================================================================================================
102 explicit SensorSet(const QSharedPointer<FWDLIB::FwdCoilSet> pFwdCoilSet);
103
104 inline int np() const;
105 inline int ncoils() const;
106 inline Eigen::VectorXd ez(int iSensor) const;
107 inline Eigen::MatrixXd ez() const;
108
109 inline Eigen::VectorXd r0(int iSensor) const;
110 inline Eigen::MatrixXd r0() const;
111
112 inline Eigen::MatrixXd rmag(int iSensor) const;
113 inline Eigen::MatrixXd rmag() const;
114
115 inline Eigen::MatrixXd cosmag(int iSensor) const;
116 inline Eigen::MatrixXd cosmag() const;
117
118 inline Eigen::MatrixXd tra(int iSensor) const;
119 inline Eigen::MatrixXd tra() const;
120
121 inline Eigen::RowVectorXd w(int iSensor) const;
122 inline Eigen::RowVectorXd w() const;
123
124 inline bool operator== (const SensorSet &b) const;
125 inline bool operator!= (const SensorSet &b) const;
126
127private:
128 //=========================================================================================================
134 void initFromFwdCoilSet(const QSharedPointer<FWDLIB::FwdCoilSet> pFwdCoilSet);
135
136 //=========================================================================================================
142 void initMatrices(int iNchan, int iNp);
143
144 Eigen::MatrixXd m_ez{Eigen::MatrixXd(0,0)};
145 Eigen::MatrixXd m_r0{Eigen::MatrixXd(0,0)};
146 Eigen::MatrixXd m_rmag{Eigen::MatrixXd(0,0)};
147 Eigen::MatrixXd m_cosmag{Eigen::MatrixXd(0,0)};
148 Eigen::MatrixXd m_tra{Eigen::MatrixXd(0,0)};
149 Eigen::RowVectorXd m_w{Eigen::RowVectorXd(0)};
150 int m_ncoils{0};
151 int m_np{0};
152};
153
154//=============================================================================================================
155// INLINE DEFINITIONS
156//=============================================================================================================
157
158inline int SensorSet::np() const
159{
160 return m_np;
161}
162
163inline int SensorSet::ncoils() const
164{
165 return m_ncoils;
166}
167
168inline Eigen::VectorXd SensorSet::ez(int iSensor) const
169{
170 return m_ez.row(iSensor);
171}
172
173inline Eigen::MatrixXd SensorSet::ez() const
174{
175 return m_ez;
176}
177
178inline Eigen::VectorXd SensorSet::r0(int iSensor) const
179{
180 return m_r0.row(iSensor);
181}
182
183inline Eigen::MatrixXd SensorSet::r0() const
184{
185 return m_r0;
186}
187
188inline Eigen::RowVectorXd SensorSet::w(int iSensor) const
189{
190 return m_w.segment(iSensor*m_np,m_np);
191}
192
193inline Eigen::RowVectorXd SensorSet::w() const
194{
195 return m_w;
196}
197
198inline Eigen::MatrixXd SensorSet::rmag(int iSensor) const
199{
200 return m_rmag.block(iSensor*m_np,0,m_np,3);
201}
202
203inline Eigen::MatrixXd SensorSet::rmag() const
204{
205 return m_rmag;
206}
207
208inline Eigen::MatrixXd SensorSet::cosmag(int iSensor) const
209{
210 return m_cosmag.block(iSensor*m_np,0,m_np,3);
211}
212
213inline Eigen::MatrixXd SensorSet::cosmag() const
214{
215 return m_cosmag;
216}
217
218inline Eigen::MatrixXd SensorSet::tra() const
219{
220 return m_tra;
221}
222//=============================================================================================================
223
224inline bool SensorSet::operator== (const SensorSet &b) const
225{
226 return (this->ez() == b.ez() &&
227 this->r0() == b.r0() &&
228 this->rmag() == b.rmag() &&
229 this->cosmag() == b.cosmag() &&
230 this->tra() == b.tra() &&
231 this->w() == b.w() &&
232 this->np() == b.np() &&
233 this->ncoils() == b.ncoils());
234}
235
236//=============================================================================================================
237
238inline bool SensorSet::operator!= (const SensorSet &b) const
239{
240 bool equal = this==&b;
241 return !(equal);
242}
243
244//=============================================================================================================
251{
252
253public:
254 typedef QSharedPointer<SensorSetCreator> SPtr;
255 typedef QSharedPointer<const SensorSetCreator> ConstSPtr;
256
257 //=========================================================================================================
261 explicit SensorSetCreator();
262
263 //=========================================================================================================
271 SensorSet updateSensorSet(const QList<FIFFLIB::FiffChInfo>& channelList,
272 const Accuracy& accuracy);
273
274private:
275 QSharedPointer<FWDLIB::FwdCoilSet> m_pCoilDefinitions{nullptr}; // the coil definitions as template
276};
277
278//=============================================================================================================
279// INLINE DEFINITIONS
280//=============================================================================================================
281
282} // namespace INVERSELIB
283
284#endif // SENSORSET_H
285
inverse library export/import macros.
#define INVERSESHARED_EXPORT
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Forward modelling (BEM, MEG/EEG lead fields).
Definition compute_fwd.h:95
Channel info descriptor.
Coordinate transformation description.
Holds a set of digitizer points.
Single MEG or EEG sensor coil with integration points, weights, and coordinate frame.
Definition fwd_coil.h:89
Collection of FwdCoil objects representing a full MEG or EEG sensor array.
Stores MEG sensor geometry (positions, orientations, weights, coil count) for a single sensor type.
Definition sensorset.h:83
Eigen::MatrixXd rmag(int iSensor) const
Definition sensorset.h:198
Eigen::MatrixXd tra() const
Definition sensorset.h:218
Eigen::MatrixXd cosmag(int iSensor) const
Definition sensorset.h:208
Eigen::MatrixXd rmag() const
Definition sensorset.h:203
QSharedPointer< const SensorSet > ConstSPtr
Definition sensorset.h:87
Eigen::MatrixXd cosmag() const
Definition sensorset.h:213
bool operator!=(const SensorSet &b) const
Definition sensorset.h:238
bool operator==(const SensorSet &b) const
Definition sensorset.h:224
Eigen::VectorXd r0(int iSensor) const
Definition sensorset.h:178
QSharedPointer< SensorSet > SPtr
Definition sensorset.h:86
Eigen::MatrixXd ez() const
Definition sensorset.h:173
SensorSet(const QSharedPointer< FWDLIB::FwdCoilSet > pFwdCoilSet)
Eigen::VectorXd ez(int iSensor) const
Definition sensorset.h:168
Eigen::MatrixXd r0() const
Definition sensorset.h:183
Eigen::RowVectorXd w(int iSensor) const
Definition sensorset.h:188
Eigen::MatrixXd tra(int iSensor) const
Eigen::RowVectorXd w() const
Definition sensorset.h:193
QSharedPointer< const SensorSetCreator > ConstSPtr
Definition sensorset.h:255
QSharedPointer< SensorSetCreator > SPtr
Definition sensorset.h:254
SensorSet updateSensorSet(const QList< FIFFLIB::FiffChInfo > &channelList, const Accuracy &accuracy)