v2.0.0
Loading...
Searching...
No Matches
hpifit.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef HPIFIT_H
38#define HPIFIT_H
39
40//=============================================================================================================
41// INCLUDES
42//=============================================================================================================
43
44#include "../inverse_global.h"
45#include "fiff/fiff_ch_info.h"
46#include "sensorset.h"
47#include "signalmodel.h"
49#include <fiff/fiff_dig_point.h>
51
52//=============================================================================================================
53// EIGEN INCLUDES
54//=============================================================================================================
55
56#include <Eigen/Core>
57
58//=============================================================================================================
59// QT INCLUDES
60//=============================================================================================================
61
62#include <QSharedPointer>
63
64//=============================================================================================================
65// FORWARD DECLARATIONS
66//=============================================================================================================
67
68namespace FWDLIB{
69 class FwdCoil;
70 class FwdCoilSet;
71}
72
73namespace FIFFLIB{
74 class FiffInfo;
75 class FiffCoordTrans;
76 class FiffDigPointSet;
77}
78
79//=============================================================================================================
80// DEFINE NAMESPACE INVERSELIB
81//=============================================================================================================
82
83namespace INVERSELIB
84{
85 class SensorSet;
86 class SignalModel;
88//=============================================================================================================
89// Declare all structures to be used
90//=============================================================================================================
91
97struct CoilParam {
98 Eigen::MatrixXd pos;
99 Eigen::MatrixXd mom;
100 Eigen::VectorXd dpfiterror;
101 Eigen::VectorXd dpfitnumitr;
102
103 CoilParam(int iNumCoils)
104 : pos(Eigen::MatrixXd(iNumCoils,3)),
105 mom(Eigen::MatrixXd::Zero(iNumCoils,3)),
106 dpfiterror(Eigen::VectorXd::Zero(iNumCoils)),
107 dpfitnumitr(Eigen::VectorXd::Zero(iNumCoils))
108 {}
109};
110
127
128//=============================================================================================================
129// INVERSELIB FORWARD DECLARATIONS
130//=============================================================================================================
131
132//=============================================================================================================
139{
140
141public:
142 typedef QSharedPointer<HPIFit> SPtr;
143 typedef QSharedPointer<const HPIFit> ConstSPtr;
144
145 //=========================================================================================================
150 explicit HPIFit();
151
152 //=========================================================================================================
158 explicit HPIFit(const SensorSet& sensorSet);
159
160 //=========================================================================================================
166 void checkForUpdate(const SensorSet& sensorSet);
167
168 //=========================================================================================================
179 void fit(const Eigen::MatrixXd& matProjectedData,
180 const Eigen::MatrixXd& matProjectors,
181 const HpiModelParameters& hpiModelParameters,
182 const Eigen::MatrixXd& matCoilsHead,
183 HpiFitResult& hpiFitResult);
184
185 void fit(const Eigen::MatrixXd& matProjectedData,
186 const Eigen::MatrixXd& matProjectors,
187 const HpiModelParameters& hpiModelParameters,
188 const Eigen::MatrixXd& matCoilsHead,
189 const bool bOrderFrequencies,
190 HpiFitResult& hpiFitResult);
191
192 //=========================================================================================================
206 static void storeHeadPosition(float fTime,
207 const Eigen::MatrixXf& matTransDevHead,
208 Eigen::MatrixXd& matPosition,
209 const Eigen::VectorXd& vecGoF,
210 const QVector<double>& vecError);
211
212private:
213
214 //=========================================================================================================
223 Eigen::MatrixXd computeAmplitudes(const Eigen::MatrixXd& matProjectedData,
224 const HpiModelParameters& hpiModelParameters);
225
226 //=========================================================================================================
236 Eigen::MatrixXd computeSeedPoints(const Eigen::MatrixXd& matAmplitudes,
237 const FIFFLIB::FiffCoordTrans& transDevHead,
238 const QVector<double>& vecError,
239 const Eigen::MatrixXd& matCoilsHead);
240
241 //=========================================================================================================
255 CoilParam dipfit(const Eigen::MatrixXd matCoilsSeed,
256 const SensorSet& sensors,
257 const Eigen::MatrixXd &matData,
258 const int iNumCoils,
259 const Eigen::MatrixXd &t_matProjectors,
260 const int iMaxIterations,
261 const float fAbortError);
262
263 //=========================================================================================================
270 Eigen::VectorXd computeGoF(const Eigen::VectorXd& vecDipFitError);
271
272 //=========================================================================================================
283 FIFFLIB::FiffCoordTrans computeDeviceHeadTransformation(const Eigen::MatrixXd& matCoilsDev,
284 const Eigen::MatrixXd& matCoilsHead);
285
286 //=========================================================================================================
295 QVector<double> computeEstimationError(const Eigen::MatrixXd& matCoilsDev,
296 const Eigen::MatrixXd& matCoilsHead,
297 const FIFFLIB::FiffCoordTrans& transDevHead);
298
299 //=========================================================================================================
307 FIFFLIB::FiffDigPointSet getFittedPointSet(const Eigen::MatrixXd& matCoilsDev);
308
309 //=========================================================================================================
318 Eigen::Matrix4d computeTransformation(Eigen::MatrixXd matNH,
319 Eigen::MatrixXd matBT);
320
321 //=========================================================================================================
330 std::vector<int> findCoilOrder(const Eigen::MatrixXd& matCoilsDev,
331 const Eigen::MatrixXd& matCoilsHead);
332
333 //=========================================================================================================
342 Eigen::MatrixXd order(const std::vector<int>& vecOrder,
343 const Eigen::MatrixXd& matToOrder);
344
345 QVector<int> order(const std::vector<int>& vecOrder,
346 const QVector<int>& vecToOrder);
347
348 //=========================================================================================================
358 double objectTrans(const Eigen::MatrixXd& matHeadCoil,
359 const Eigen::MatrixXd& matCoilsDev,
360 const Eigen::MatrixXd& matTrans);
361
362 SensorSet m_sensors;
363 SignalModel m_signalModel;
364
365};
366
367//=============================================================================================================
368// INLINE DEFINITIONS
369//=============================================================================================================
370
371} //NAMESPACE
372
373#ifndef metatype_HpiFitResult
374#define metatype_HpiFitResult
376#endif
377
378#endif // HPIFIT_H
FiffDigPointSet class declaration.
FiffChInfo class declaration.
FiffCoordTrans class declaration.
FiffDigPoint class declaration.
Q_DECLARE_METATYPE(Eigen::MatrixXf)
SensorSet class declaration.
SignalModel class declaration.
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
Coordinate transformation description.
Holds a set of digitizer points.
FIFF measurement file information.
Definition fiff_info.h:85
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.
Estimated dipole parameters (position, moment, goodness-of-fit) for a single HPI coil.
Definition hpifit.h:97
Eigen::VectorXd dpfiterror
Definition hpifit.h:100
Eigen::MatrixXd pos
Definition hpifit.h:98
Eigen::VectorXd dpfitnumitr
Definition hpifit.h:101
CoilParam(int iNumCoils)
Definition hpifit.h:103
Eigen::MatrixXd mom
Definition hpifit.h:99
Complete HPI fit output: per-coil dipole parameters, head-to-device transform, fit error,...
Definition hpifit.h:116
FIFFLIB::FiffCoordTrans devHeadTrans
Definition hpifit.h:119
Eigen::VectorXd GoF
Definition hpifit.h:121
FIFFLIB::FiffDigPointSet fittedCoils
Definition hpifit.h:118
QVector< double > errorDistances
Definition hpifit.h:120
QVector< int > hpiFreqs
Definition hpifit.h:117
void checkForUpdate(const SensorSet &sensorSet)
Definition hpifit.cpp:102
QSharedPointer< const HPIFit > ConstSPtr
Definition hpifit.h:143
void fit(const Eigen::MatrixXd &matProjectedData, const Eigen::MatrixXd &matProjectors, const HpiModelParameters &hpiModelParameters, const Eigen::MatrixXd &matCoilsHead, const bool bOrderFrequencies, HpiFitResult &hpiFitResult)
static void storeHeadPosition(float fTime, const Eigen::MatrixXf &matTransDevHead, Eigen::MatrixXd &matPosition, const Eigen::VectorXd &vecGoF, const QVector< double > &vecError)
Definition hpifit.cpp:527
void fit(const Eigen::MatrixXd &matProjectedData, const Eigen::MatrixXd &matProjectors, const HpiModelParameters &hpiModelParameters, const Eigen::MatrixXd &matCoilsHead, HpiFitResult &hpiFitResult)
QSharedPointer< HPIFit > SPtr
Definition hpifit.h:142
Configuration parameters for the HPI signal model (line frequency, coil frequencies,...
Stores MEG sensor geometry (positions, orientations, weights, coil count) for a single sensor type.
Definition sensorset.h:83
Generates the forward sinusoidal model matrix for HPI coil signals at known drive frequencies.
Definition signalmodel.h:74