v2.0.0
Loading...
Searching...
No Matches
rtsensorinterpolationmatworker.cpp
Go to the documentation of this file.
1//=============================================================================================================
34
35//=============================================================================================================
36// INCLUDES
37//=============================================================================================================
38
40#include <fwd/fwd_field_map.h>
41
42#include <fiff/fiff_ch_info.h>
43#include <fiff/fiff_constants.h>
45#include <fwd/fwd_coil_set.h>
46#include <fs/surface.h>
47
48#include <QCoreApplication>
49#include <QDebug>
50#include <QVector3D>
51#include <QMatrix4x4>
52#include <cmath>
53
54using namespace FIFFLIB;
55using namespace DISP3DRHILIB;
56
57//=============================================================================================================
58// ANONYMOUS HELPERS
59//=============================================================================================================
60
61namespace
62{
63
67Eigen::Vector3f applyTransform(const Eigen::Vector3f &point,
68 const FiffCoordTrans &trans)
69{
70 if (trans.isEmpty()) return point;
71 float r[3] = {point.x(), point.y(), point.z()};
73 return Eigen::Vector3f(r[0], r[1], r[2]);
74}
75
79QMatrix4x4 toQMatrix4x4(const Eigen::Matrix4f &m)
80{
81 QMatrix4x4 q;
82 for (int r = 0; r < 4; ++r)
83 for (int c = 0; c < 4; ++c)
84 q(r, c) = m(r, c);
85 return q;
86}
87
88} // anonymous namespace
89
90//=============================================================================================================
91// MEMBER METHODS
92//=============================================================================================================
93
95 : QObject(parent)
96{
97}
98
99//=============================================================================================================
100
102{
103 QMutexLocker locker(&m_mutex);
104 m_evoked = evoked;
105 m_hasEvoked = (evoked.data.rows() > 0 && evoked.data.cols() > 0);
106}
107
108//=============================================================================================================
109
111 bool applySensorTrans)
112{
113 QMutexLocker locker(&m_mutex);
114 m_headToMriTrans = trans;
115 m_applySensorTrans = applySensorTrans;
116}
117
118//=============================================================================================================
119
121{
122 QMutexLocker locker(&m_mutex);
123 m_megOnHead = onHead;
124}
125
126//=============================================================================================================
127
129 const Eigen::MatrixX3f &vertices,
130 const Eigen::MatrixX3f &normals,
131 const Eigen::MatrixX3i &triangles)
132{
133 QMutexLocker locker(&m_mutex);
134 m_megSurfaceKey = surfaceKey;
135 m_megVertices = vertices;
136 m_megNormals = normals;
137 m_megTriangles = triangles;
138 m_hasMegSurface = (vertices.rows() > 0);
139}
140
141//=============================================================================================================
142
144 const Eigen::MatrixX3f &vertices)
145{
146 QMutexLocker locker(&m_mutex);
147 m_eegSurfaceKey = surfaceKey;
148 m_eegVertices = vertices;
149 m_hasEegSurface = (vertices.rows() > 0);
150}
151
152//=============================================================================================================
153
155{
156 QMutexLocker locker(&m_mutex);
157 m_bads = bads;
158}
159
160//=============================================================================================================
161
163{
164 // ── Snapshot parameters under lock ─────────────────────────────────
165 FiffEvoked evoked;
166 FiffCoordTrans headToMriTrans;
167 bool applySensorTrans;
168 bool megOnHead;
169 QString megSurfaceKey, eegSurfaceKey;
170 Eigen::MatrixX3f megVerts, megNorms, eegVerts;
171 Eigen::MatrixX3i megTris;
172 bool hasMegSurface, hasEegSurface, hasEvoked;
173 QStringList bads;
174
175 {
176 QMutexLocker locker(&m_mutex);
177 evoked = m_evoked;
178 hasEvoked = m_hasEvoked;
179 headToMriTrans = m_headToMriTrans;
180 applySensorTrans = m_applySensorTrans;
181 megOnHead = m_megOnHead;
182 megSurfaceKey = m_megSurfaceKey;
183 eegSurfaceKey = m_eegSurfaceKey;
184 megVerts = m_megVertices;
185 megNorms = m_megNormals;
186 megTris = m_megTriangles;
187 hasMegSurface = m_hasMegSurface;
188 eegVerts = m_eegVertices;
189 hasEegSurface = m_hasEegSurface;
190 bads = m_bads;
191 }
192
193 if (!hasEvoked) {
194 qDebug() << "RtSensorInterpolationMatWorker: No evoked data set, skipping.";
195 return;
196 }
197
198 qDebug() << "RtSensorInterpolationMatWorker: Computing mapping matrices...";
199
200 // ── Build coordinate transforms ────────────────────────────────────
201 bool hasDevHead = false;
202 QMatrix4x4 devHeadQt;
203 if (!evoked.info.dev_head_t.isEmpty() &&
206 !evoked.info.dev_head_t.trans.isIdentity()) {
207 hasDevHead = true;
208 devHeadQt = toQMatrix4x4(evoked.info.dev_head_t.trans);
209 }
210
211 QMatrix4x4 headToMri;
212 if (applySensorTrans && !headToMriTrans.isEmpty()) {
213 headToMri = toQMatrix4x4(headToMriTrans.trans);
214 }
215
216 // ── Classify channels ──────────────────────────────────────────────
217 QList<FiffChInfo> megChs, eegChs;
218 QVector<int> megPick, eegPick;
219
220 for (int k = 0; k < evoked.info.chs.size(); ++k) {
221 const auto &ch = evoked.info.chs[k];
222 if (bads.contains(ch.ch_name)) continue;
223
224 QVector3D pos(ch.chpos.r0(0), ch.chpos.r0(1), ch.chpos.r0(2));
225
226 if (ch.kind == FIFFV_MEG_CH) {
227 if (hasDevHead) pos = devHeadQt.map(pos);
228 if (applySensorTrans && !headToMriTrans.isEmpty()) pos = headToMri.map(pos);
229 megPick.append(k);
230 megChs.append(ch);
231 } else if (ch.kind == FIFFV_EEG_CH) {
232 if (applySensorTrans && !headToMriTrans.isEmpty()) pos = headToMri.map(pos);
233 eegPick.append(k);
234 eegChs.append(ch);
235 }
236 }
237
238 // ── Constants (matching MNE-Python) ────────────────────────────────
239 constexpr float kIntrad = 0.06f;
240 constexpr float kMegMiss = 1e-4f;
241 constexpr float kEegMiss = 1e-3f;
242 const Eigen::Vector3f defaultOrigin(0.0f, 0.0f, 0.04f);
243
244 FiffCoordTrans headMri = (applySensorTrans && !headToMriTrans.isEmpty())
245 ? headToMriTrans : FiffCoordTrans();
246 FiffCoordTrans devHead = (!evoked.info.dev_head_t.isEmpty() &&
249 ? evoked.info.dev_head_t : FiffCoordTrans();
250
251 // ── MEG mapping ────────────────────────────────────────────────────
252 if (hasMegSurface && !megChs.isEmpty()) {
253 Eigen::MatrixX3f norms = megNorms;
254
255 // Recompute normals if missing
256 if (norms.rows() != megVerts.rows()) {
257 if (megTris.rows() > 0) {
258 norms = FSLIB::Surface::compute_normals(megVerts, megTris);
259 }
260 }
261
262 if (megVerts.rows() > 0 && norms.rows() == megVerts.rows()) {
263 const QString coilPath = QCoreApplication::applicationDirPath()
264 + "/../resources/general/coilDefinitions/coil_def.dat";
265 std::unique_ptr<FWDLIB::FwdCoilSet> templates(
267
268 if (templates) {
269 FiffCoordTrans devToTarget;
270 if (megOnHead && !headMri.isEmpty()) {
271 if (!devHead.isEmpty()) {
272 devToTarget = FiffCoordTrans::combine(
274 devHead, headMri);
275 }
276 } else if (!devHead.isEmpty()) {
277 devToTarget = devHead;
278 }
279
280 Eigen::Vector3f origin = defaultOrigin;
281 if (megOnHead && !headMri.isEmpty()) {
282 origin = applyTransform(origin, headMri);
283 }
284
285 std::unique_ptr<FWDLIB::FwdCoilSet> coils(templates->create_meg_coils(
286 megChs, megChs.size(), FWD_COIL_ACCURACY_NORMAL, devToTarget));
287
288 if (coils && coils->ncoil > 0) {
290 *coils, megVerts, norms, origin, kIntrad, kMegMiss);
291
292 if (mat && mat->rows() > 0) {
293 qDebug() << "RtSensorInterpolationMatWorker: MEG mapping computed:"
294 << mat->rows() << "x" << mat->cols();
295 emit newMegMappingAvailable(megSurfaceKey, mat, megPick);
296 }
297 }
298 }
299 }
300 }
301
302 // ── EEG mapping ────────────────────────────────────────────────────
303 if (hasEegSurface && !eegChs.isEmpty()) {
304 if (eegVerts.rows() > 0) {
305 Eigen::Vector3f origin = defaultOrigin;
306 if (!headMri.isEmpty()) origin = applyTransform(origin, headMri);
307
308 std::unique_ptr<FWDLIB::FwdCoilSet> eegCoils(
310 eegChs, eegChs.size(), headMri));
311
312 if (eegCoils && eegCoils->ncoil > 0) {
314 *eegCoils, eegVerts, origin, kIntrad, kEegMiss);
315
316 if (mat && mat->rows() > 0) {
317 qDebug() << "RtSensorInterpolationMatWorker: EEG mapping computed:"
318 << mat->rows() << "x" << mat->cols();
319 emit newEegMappingAvailable(eegSurfaceKey, mat, eegPick);
320 }
321 }
322 }
323 }
324
325 qDebug() << "RtSensorInterpolationMatWorker: Mapping computation complete.";
326}
Fiff constants.
#define FIFFV_EEG_CH
#define FIFFV_COORD_DEVICE
#define FIFFV_MEG_CH
#define FIFFV_COORD_HEAD
#define FIFFV_COORD_MRI
#define FIFFV_MOVE
FiffChInfo class declaration.
FiffCoordTrans class declaration.
RtSensorInterpolationMatWorker class declaration.
Surface class declaration.
Sphere-model field mapping.
#define FWD_COIL_ACCURACY_NORMAL
Definition fwd_coil.h:70
FwdCoilSet class declaration.
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
QMatrix4x4 toQMatrix4x4(const Eigen::Matrix4f &m)
3-D brain visualisation using the Qt RHI rendering backend.
void newMegMappingAvailable(const QString &surfaceKey, QSharedPointer< Eigen::MatrixXf > mappingMat, const QVector< int > &pick)
void newEegMappingAvailable(const QString &surfaceKey, QSharedPointer< Eigen::MatrixXf > mappingMat, const QVector< int > &pick)
void setMegSurface(const QString &surfaceKey, const Eigen::MatrixX3f &vertices, const Eigen::MatrixX3f &normals, const Eigen::MatrixX3i &triangles)
void setTransform(const FIFFLIB::FiffCoordTrans &trans, bool applySensorTrans)
void setEegSurface(const QString &surfaceKey, const Eigen::MatrixX3f &vertices)
Coordinate transformation description.
Eigen::Matrix< float, 4, 4, Eigen::DontAlign > trans
Eigen::MatrixXd data
QList< FiffChInfo > chs
FiffCoordTrans dev_head_t
static Eigen::MatrixX3f compute_normals(const Eigen::MatrixX3f &rr, const Eigen::MatrixX3i &tris)
Definition surface.cpp:133
static FwdCoilSet * create_eeg_els(const QList< FIFFLIB::FiffChInfo > &chs, int nch, const FIFFLIB::FiffCoordTrans &t=FIFFLIB::FiffCoordTrans())
static FwdCoilSet * read_coil_defs(const QString &name)
static QSharedPointer< Eigen::MatrixXf > computeEegMapping(const FwdCoilSet &coils, const Eigen::MatrixX3f &vertices, const Eigen::Vector3f &origin, float intrad=0.06f, float miss=1e-3f)
static QSharedPointer< Eigen::MatrixXf > computeMegMapping(const FwdCoilSet &coils, const Eigen::MatrixX3f &vertices, const Eigen::MatrixX3f &normals, const Eigen::Vector3f &origin, float intrad=0.06f, float miss=1e-4f)
static FiffCoordTrans combine(int from, int to, const FiffCoordTrans &t1, const FiffCoordTrans &t2)
Eigen::MatrixX3f apply_trans(const Eigen::MatrixX3f &rr, bool do_move=true) const