MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
sensorset.h
Go to the documentation of this file.
1//=============================================================================================================
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 FiffCoordTransOld;
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
81
82public:
83 typedef QSharedPointer<SensorSet> SPtr;
84 typedef QSharedPointer<const SensorSet> ConstSPtr;
86 //=========================================================================================================
90 explicit SensorSet() = default;
91
92 //=========================================================================================================
99 explicit SensorSet(const QSharedPointer<FWDLIB::FwdCoilSet> pFwdCoilSet);
100
101 inline int np() const;
102 inline int ncoils() const;
103 inline Eigen::VectorXd ez(int iSensor) const;
104 inline Eigen::MatrixXd ez() const;
105
106 inline Eigen::VectorXd r0(int iSensor) const;
107 inline Eigen::MatrixXd r0() const;
108
109 inline Eigen::MatrixXd rmag(int iSensor) const;
110 inline Eigen::MatrixXd rmag() const;
111
112 inline Eigen::MatrixXd cosmag(int iSensor) const;
113 inline Eigen::MatrixXd cosmag() const;
114
115 inline Eigen::MatrixXd tra(int iSensor) const;
116 inline Eigen::MatrixXd tra() const;
117
118 inline Eigen::RowVectorXd w(int iSensor) const;
119 inline Eigen::RowVectorXd w() const;
120
121 inline bool operator== (const SensorSet &b) const;
122 inline bool operator!= (const SensorSet &b) const;
123
124private:
125 //=========================================================================================================
131 void initFromFwdCoilSet(const QSharedPointer<FWDLIB::FwdCoilSet> pFwdCoilSet);
132
133 //=========================================================================================================
139 void initMatrices(int iNchan, int iNp);
140
141 Eigen::MatrixXd m_ez{Eigen::MatrixXd(0,0)};
142 Eigen::MatrixXd m_r0{Eigen::MatrixXd(0,0)};
143 Eigen::MatrixXd m_rmag{Eigen::MatrixXd(0,0)};
144 Eigen::MatrixXd m_cosmag{Eigen::MatrixXd(0,0)};
145 Eigen::MatrixXd m_tra{Eigen::MatrixXd(0,0)};
146 Eigen::RowVectorXd m_w{Eigen::RowVectorXd(0)};
147 int m_ncoils{0};
148 int m_np{0};
149};
150
151//=============================================================================================================
152// INLINE DEFINITIONS
153//=============================================================================================================
154
155inline int SensorSet::np() const
156{
157 return m_np;
158}
159
160inline int SensorSet::ncoils() const
161{
162 return m_ncoils;
163}
164
165inline Eigen::VectorXd SensorSet::ez(int iSensor) const
166{
167 return m_ez.row(iSensor);
168}
169
170inline Eigen::MatrixXd SensorSet::ez() const
171{
172 return m_ez;
173}
174
175inline Eigen::VectorXd SensorSet::r0(int iSensor) const
176{
177 return m_r0.row(iSensor);
178}
179
180inline Eigen::MatrixXd SensorSet::r0() const
181{
182 return m_r0;
183}
184
185inline Eigen::RowVectorXd SensorSet::w(int iSensor) const
186{
187 return m_w.segment(iSensor*m_np,m_np);
188}
189
190inline Eigen::RowVectorXd SensorSet::w() const
191{
192 return m_w;
193}
194
195inline Eigen::MatrixXd SensorSet::rmag(int iSensor) const
196{
197 return m_rmag.block(iSensor*m_np,0,m_np,3);
198}
199
200inline Eigen::MatrixXd SensorSet::rmag() const
201{
202 return m_rmag;
203}
204
205inline Eigen::MatrixXd SensorSet::cosmag(int iSensor) const
206{
207 return m_cosmag.block(iSensor*m_np,0,m_np,3);
208}
209
210inline Eigen::MatrixXd SensorSet::cosmag() const
211{
212 return m_cosmag;
213}
214
215inline Eigen::MatrixXd SensorSet::tra() const
216{
217 return m_tra;
218}
219//=============================================================================================================
220
221inline bool SensorSet::operator== (const SensorSet &b) const
222{
223 return (this->ez() == b.ez() &&
224 this->r0() == b.r0() &&
225 this->rmag() == b.rmag() &&
226 this->cosmag() == b.cosmag() &&
227 this->tra() == b.tra() &&
228 this->w() == b.w() &&
229 this->np() == b.np() &&
230 this->ncoils() == b.ncoils());
231}
232
233//=============================================================================================================
234
235inline bool SensorSet::operator!= (const SensorSet &b) const
236{
237 bool equal = this==&b;
238 return !(equal);
239}
240
241//=============================================================================================================
248{
249
250public:
251 typedef QSharedPointer<SensorSetCreator> SPtr;
252 typedef QSharedPointer<const SensorSetCreator> ConstSPtr;
254 //=========================================================================================================
258 explicit SensorSetCreator();
259
260 //=========================================================================================================
268 SensorSet updateSensorSet(const QList<FIFFLIB::FiffChInfo>& channelList,
269 const Accuracy& accuracy);
270
271private:
272 QSharedPointer<FWDLIB::FwdCoilSet> m_pCoilDefinitions{nullptr}; // the coil definitions as template
273};
274
275//=============================================================================================================
276// INLINE DEFINITIONS
277//=============================================================================================================
278
279} // namespace INVERSELIB
280
281#endif // SENSORSET_H
282
#define INVERSESHARED_EXPORT
QSharedPointer< const SensorSet > ConstSPtr
Definition sensorset.h:84
QSharedPointer< SensorSet > SPtr
Definition sensorset.h:83
SensorSet(const QSharedPointer< FWDLIB::FwdCoilSet > pFwdCoilSet)
Brief description of this class.
Definition sensorset.h:248
QSharedPointer< const SensorSetCreator > ConstSPtr
Definition sensorset.h:252
QSharedPointer< SensorSetCreator > SPtr
Definition sensorset.h:251