v2.0.0
Loading...
Searching...
No Matches
hpifitdata.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef HPIFITDATA_H
38#define HPIFITDATA_H
39
40//=============================================================================================================
41// INCLUDES
42//=============================================================================================================
43
44#include "../inverse_global.h"
45#include "hpifit.h"
46#include "sensorset.h"
47
48//=============================================================================================================
49// EIGEN INCLUDES
50//=============================================================================================================
51
52#include <Eigen/Core>
53
54//=============================================================================================================
55// QT INCLUDES
56//=============================================================================================================
57
58#include <QSharedPointer>
59
60//=============================================================================================================
61// FORWARD DECLARATIONS
62//=============================================================================================================
63
64namespace FIFFLIB{
65 class FiffInfo;
66 class FiffCoordTrans;
67 class FiffDigPointSet;
68}
69
70//=============================================================================================================
71// DEFINE NAMESPACE INVERSELIB
72//=============================================================================================================
73
74namespace INVERSELIB
75{
76
77//=============================================================================================================
78// Declare all structures to be used
79//=============================================================================================================
86 double error;
87 Eigen::MatrixXd moment;
89};
90
91//=========================================================================================================
98 double base_arr;
99 int idx;
100};
101
102//=============================================================================================================
103// FORWARD DECLARATIONS
104//=============================================================================================================
105
106//=============================================================================================================
113{
114
115public:
116 typedef QSharedPointer<HPIFitData> SPtr;
117 typedef QSharedPointer<const HPIFitData> ConstSPtr;
118
119 //=========================================================================================================
123 explicit HPIFitData();
124
125 //=========================================================================================================
129 void doDipfitConcurrent();
130
131 Eigen::MatrixXd m_coilPos;
132 Eigen::RowVectorXd m_sensorData;
135 Eigen::MatrixXd m_matProjector;
136
139
140protected:
141 //=========================================================================================================
146 Eigen::MatrixXd magnetic_dipole(Eigen::MatrixXd matPos,
147 Eigen::MatrixXd matPnt,
148 Eigen::MatrixXd matOri);
149
150 //=========================================================================================================
159 Eigen::MatrixXd compute_leadfield(const Eigen::MatrixXd& matPos,
160 const SensorSet& sensors);
161
162 //=========================================================================================================
169 DipFitError dipfitError(const Eigen::MatrixXd& matPos,
170 const Eigen::MatrixXd& matData,
171 const SensorSet& sensors,
172 const Eigen::MatrixXd& matProjectors);
173
174 //=========================================================================================================
178 static bool compare(HPISortStruct a, HPISortStruct b);
179
180 //=========================================================================================================
186 Eigen::MatrixXd fminsearch(const Eigen::MatrixXd& matPos,
187 int iMaxiter,
188 int iMaxfun,
189 int iDisplay,
190 const Eigen::MatrixXd& matData,
191 const Eigen::MatrixXd& matProjectors,
192 const SensorSet& sensors,
193 int &iSimplexNumitr);
194};
195
196//=============================================================================================================
197// INLINE DEFINITIONS
198//=============================================================================================================
199} //NAMESPACE
200
201#endif // HPIFITDATA_H
SensorSet class declaration.
HPIFit 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).
Coordinate transformation description.
Holds a set of digitizer points.
FIFF measurement file information.
Definition fiff_info.h:85
Residual error and moment vector from a single magnetic dipole fit iteration.
Definition hpifitdata.h:85
Eigen::MatrixXd moment
Definition hpifitdata.h:87
Helper for sorting HPI coil dipole fits by matching each fit to the nearest expected coil position.
Definition hpifitdata.h:97
Eigen::MatrixXd m_coilPos
Definition hpifitdata.h:131
Eigen::MatrixXd fminsearch(const Eigen::MatrixXd &matPos, int iMaxiter, int iMaxfun, int iDisplay, const Eigen::MatrixXd &matData, const Eigen::MatrixXd &matProjectors, const SensorSet &sensors, int &iSimplexNumitr)
Eigen::MatrixXd m_matProjector
Definition hpifitdata.h:135
QSharedPointer< const HPIFitData > ConstSPtr
Definition hpifitdata.h:117
static bool compare(HPISortStruct a, HPISortStruct b)
Eigen::MatrixXd magnetic_dipole(Eigen::MatrixXd matPos, Eigen::MatrixXd matPnt, Eigen::MatrixXd matOri)
Eigen::RowVectorXd m_sensorData
Definition hpifitdata.h:132
QSharedPointer< HPIFitData > SPtr
Definition hpifitdata.h:116
Eigen::MatrixXd compute_leadfield(const Eigen::MatrixXd &matPos, const SensorSet &sensors)
DipFitError dipfitError(const Eigen::MatrixXd &matPos, const Eigen::MatrixXd &matData, const SensorSet &sensors, const Eigen::MatrixXd &matProjectors)
DipFitError m_errorInfo
Definition hpifitdata.h:133
Stores MEG sensor geometry (positions, orientations, weights, coil count) for a single sensor type.
Definition sensorset.h:83