v2.0.0
Loading...
Searching...
No Matches
inv_hpi_fit.cpp
Go to the documentation of this file.
1//=============================================================================================================
37
38//=============================================================================================================
39// INCLUDES
40//=============================================================================================================
41
42#include "inv_hpi_fit.h"
43#include "inv_hpi_fit_data.h"
44#include "inv_sensor_set.h"
46#include "inv_signal_model.h"
48
49#include <utils/ioutils.h>
50
51#define _USE_MATH_DEFINES
52#include <cmath>
53#include <iostream>
54#include <vector>
55#include <numeric>
56#include <fiff/fiff_cov.h>
58#include <fstream>
59
60#include <fwd/fwd_coil_set.h>
61
62//=============================================================================================================
63// EIGEN INCLUDES
64//=============================================================================================================
65
66#include <Eigen/Dense>
67
68//=============================================================================================================
69// QT INCLUDES
70//=============================================================================================================
71
72#include <QVector>
73#include <QFuture>
74#include <QtConcurrent/QtConcurrent>
75
76//=============================================================================================================
77// USED NAMESPACES
78//=============================================================================================================
79
80using namespace Eigen;
81using namespace INVLIB;
82using namespace FIFFLIB;
83using namespace FWDLIB;
84
85//=============================================================================================================
86// DEFINE GLOBAL METHODS
87//=============================================================================================================
88
89//=============================================================================================================
90// DEFINE MEMBER METHODS
91//=============================================================================================================
92
93//=============================================================================================================
94
96 : m_sensors(sensorSet),
97 m_signalModel(InvSignalModel())
98{
99
100}
101
102//=============================================================================================================
103
105{
106 if(m_sensors != sensorSet) {
107 m_sensors = sensorSet;
108 }
109}
110
111//=============================================================================================================
112
113void InvHpiFit::fit(const MatrixXd& matProjectedData,
114 const MatrixXd& matProjectors,
115 const InvHpiModelParameters& hpiModelParameters,
116 const MatrixXd& matCoilsHead,
117 HpiFitResult& hpiFitResult)
118{
119 fit(matProjectedData,matProjectors,hpiModelParameters,matCoilsHead,false,hpiFitResult);
120}
121
122//=============================================================================================================
123
124void InvHpiFit::fit(const MatrixXd& matProjectedData,
125 const MatrixXd& matProjectors,
126 const InvHpiModelParameters& hpiModelParameters,
127 const MatrixXd& matCoilsHead,
128 const bool bOrderFrequencies,
129 HpiFitResult& hpiFitResult)
130{
131 if(matProjectedData.rows() != matProjectors.rows()) {
132 std::cout<< "InvHpiFit::fit - Projector and data dimensions do not match. Returning."<<std::endl;
133 return;
134 } else if(hpiModelParameters.iNHpiCoils()!= matCoilsHead.rows()) {
135 std::cout<< "InvHpiFit::fit - Number of coils and hpi digitizers do not match. Returning."<<std::endl;
136 return;
137 } else if(matProjectedData.rows()==0 || matProjectors.rows()==0) {
138 std::cout<< "InvHpiFit::fit - No data or Projectors passed. Returning."<<std::endl;
139 return;
140 } else if(m_sensors.ncoils() != matProjectedData.rows()) {
141 std::cout<< "InvHpiFit::fit - Number of channels in sensors and data do not match. Returning."<<std::endl;
142 return;
143 }
144
145 const MatrixXd matAmplitudes = computeAmplitudes(matProjectedData,
146 hpiModelParameters);
147
148 const MatrixXd matCoilsSeed = computeSeedPoints(matAmplitudes,
149 hpiFitResult.devHeadTrans,
150 hpiFitResult.errorDistances,
151 matCoilsHead);
152
153 CoilParam fittedCoilParams = dipfit(matCoilsSeed,
154 m_sensors,
155 matAmplitudes,
156 hpiModelParameters.iNHpiCoils(),
157 matProjectors,
158 500,
159 1e-9f);
160
161 if(bOrderFrequencies) {
162 const std::vector<int> vecOrder = findCoilOrder(fittedCoilParams.pos,
163 matCoilsHead);
164
165 fittedCoilParams.pos = order(vecOrder,fittedCoilParams.pos);
166 hpiFitResult.hpiFreqs = order(vecOrder,hpiModelParameters.vecHpiFreqs());
167 }
168
169 hpiFitResult.GoF = computeGoF(fittedCoilParams.dpfiterror);
170
171 hpiFitResult.fittedCoils = getFittedPointSet(fittedCoilParams.pos);
172
173 hpiFitResult.devHeadTrans = computeDeviceHeadTransformation(fittedCoilParams.pos,
174 matCoilsHead);
175
176 hpiFitResult.errorDistances = computeEstimationError(fittedCoilParams.pos,
177 matCoilsHead,
178 hpiFitResult.devHeadTrans);
179}
180
181//=============================================================================================================
182
183Eigen::MatrixXd InvHpiFit::computeAmplitudes(const Eigen::MatrixXd& matProjectedData,
184 const InvHpiModelParameters& hpiModelParameters)
185{
186 // fit model
187 MatrixXd matTopo = m_signalModel.fitData(hpiModelParameters,matProjectedData);
188 matTopo.transposeInPlace();
189
190 // split into sine and cosine amplitudes
191 const int iNumCoils = hpiModelParameters.iNHpiCoils();
192
193 MatrixXd matAmpSine(matProjectedData.cols(), iNumCoils);
194 MatrixXd matAmpCosine(matProjectedData.cols(), iNumCoils);
195
196 matAmpSine = matTopo.leftCols(iNumCoils);
197 matAmpCosine = matTopo.middleCols(iNumCoils,iNumCoils);
198
199 // Select sine or cosine component depending on their contributions to the amplitudes
200 for(int j = 0; j < iNumCoils; ++j) {
201 float fNS = 0.0;
202 float fNC = 0.0;
203 fNS = matAmpSine.col(j).array().square().sum();
204 fNC = matAmpCosine.col(j).array().square().sum();
205 if(fNC > fNS) {
206 matAmpSine.col(j) = matAmpCosine.col(j);
207 }
208 }
209
210 return matAmpSine;
211}
212
213//=============================================================================================================
214
215Eigen::MatrixXd InvHpiFit::computeSeedPoints(const Eigen::MatrixXd& matAmplitudes,
216 const FIFFLIB::FiffCoordTrans& transDevHead,
217 const QVector<double>& vecError,
218 const Eigen::MatrixXd& matCoilsHead)
219{
220 const int iNumCoils = matCoilsHead.rows();
221 MatrixXd matCoilsSeed = MatrixXd::Zero(iNumCoils,3);
222
223 const double dError = std::accumulate(vecError.begin(), vecError.end(), .0) / vecError.size();
224
225 if(transDevHead.trans != MatrixXd::Identity(4,4).cast<float>() && dError < 0.010) {
226 // if good last fit, use old trafo
227 matCoilsSeed = transDevHead.apply_inverse_trans(matCoilsHead.cast<float>()).cast<double>();
228 } else {
229 // if not, find max amplitudes in channels
230 VectorXi vecChIdcs(iNumCoils);
231
232 for (int j = 0; j < iNumCoils; j++) {
233 int iChIdx = 0;
234 VectorXd::Index indMax;
235 matAmplitudes.col(j).maxCoeff(&indMax);
236 if(indMax < m_sensors.ncoils()) {
237 iChIdx = indMax;
238 }
239 vecChIdcs(j) = iChIdx;
240 }
241 // and go 3 cm inwards from max channels
242 for (int j = 0; j < vecChIdcs.rows(); ++j) {
243 if(vecChIdcs(j) < m_sensors.ncoils()) {
244 Vector3d r0 = m_sensors.r0(vecChIdcs(j));
245 Vector3d ez = m_sensors.ez(vecChIdcs(j));
246 matCoilsSeed.row(j) = (-1 * ez * 0.03 + r0).transpose();
247 }
248 }
249 }
250 return matCoilsSeed;
251}
252
253//=============================================================================================================
254
255CoilParam InvHpiFit::dipfit(const MatrixXd matCoilsSeed,
256 const InvSensorSet& sensors,
257 const MatrixXd& matData,
258 const int iNumCoils,
259 const MatrixXd& matProjectors,
260 const int iMaxIterations,
261 const float fAbortError)
262{
263 //Do this in conncurrent mode
264 //Generate QList structure which can be handled by the QConcurrent framework
265 QList<InvHpiFitData> lCoilData;
266
267 for(qint32 i = 0; i < iNumCoils; ++i) {
268 InvHpiFitData coilData;
269 coilData.m_coilPos = matCoilsSeed.row(i);
270 coilData.m_sensorData = matData.col(i);
271 coilData.m_sensors = sensors;
272 coilData.m_matProjector = matProjectors;
273 coilData.m_iMaxIterations = iMaxIterations;
274 coilData.m_fAbortError = fAbortError;
275
276 lCoilData.append(coilData);
277 }
278 //Do the concurrent filtering
279 CoilParam coil(iNumCoils);
280
281 if(!lCoilData.isEmpty()) {
282 // //Do sequential
283 // for(int l = 0; l < lCoilData.size(); ++l) {
284 // doDipfitConcurrent(lCoilData[l]);
285 // }
286
287 //Do concurrent
288 QFuture<void> future = QtConcurrent::map(lCoilData,
290 future.waitForFinished();
291
292 //Transform results to final coil information
293 for(qint32 i = 0; i < lCoilData.size(); ++i) {
294 coil.pos.row(i) = lCoilData.at(i).m_coilPos;
295 coil.mom = lCoilData.at(i).m_errorInfo.moment.transpose();
296 coil.dpfiterror(i) = lCoilData.at(i).m_errorInfo.error;
297 coil.dpfitnumitr(i) = lCoilData.at(i).m_errorInfo.numIterations;
298
299 //std::cout<<std::endl<< "InvHpiFit::dipfit - Itr steps for coil " << i << " =" <<coil.dpfitnumitr(i);
300 }
301 }
302 return coil;
303}
304
305//=============================================================================================================
306
307std::vector<int> InvHpiFit::findCoilOrder(const MatrixXd& matCoilsDev,
308 const MatrixXd& matCoilsHead)
309{
310 // extract digitized and fitted coils
311 MatrixXd matCoilTemp = matCoilsDev;
312 const int iNumCoils = matCoilsDev.rows();
313
314 std::vector<int> vecOrder(iNumCoils);
315 std::iota(vecOrder.begin(), vecOrder.end(), 0);;
316
317 // maximum 10 mm mean error
318 const double dErrorMin = 0.010;
319 double dErrorActual = 0.0;
320 double dErrorBest = dErrorMin;
321
322 MatrixXd matTrans(4,4);
323 std::vector<int> vecOrderBest = vecOrder;
324
325 bool bSuccess = false;
326 // permutation
327 do {
328 for(int i = 0; i < iNumCoils; i++) {
329 matCoilTemp.row(i) = matCoilsDev.row(vecOrder[i]);
330 }
331 matTrans = computeTransformation(matCoilsHead,matCoilTemp);
332 dErrorActual = objectTrans(matCoilsHead,matCoilTemp,matTrans);
333 if(dErrorActual < dErrorMin && dErrorActual < dErrorBest) {
334 // exit
335 dErrorBest = dErrorActual;
336 vecOrderBest = vecOrder;
337 bSuccess = true;
338 }
339 } while (std::next_permutation(vecOrder.begin(), vecOrder.end()));
340 return vecOrderBest;
341}
342
343//=============================================================================================================
344
345double InvHpiFit::objectTrans(const MatrixXd& matHeadCoil,
346 const MatrixXd& matCoil,
347 const MatrixXd& matTrans)
348{
349 // Compute the fiducial registration error - the lower, the better.
350 const int iNumCoils = matHeadCoil.rows();
351 MatrixXd matTemp = matCoil;
352
353 // homogeneous coordinates
354 matTemp.conservativeResize(matCoil.rows(),matCoil.cols()+1);
355 matTemp.block(0,3,iNumCoils,1).setOnes();
356 matTemp.transposeInPlace();
357
358 // apply transformation
359 MatrixXd matTestPos = matTrans * matTemp;
360
361 // remove
362 MatrixXd matDiff = matTestPos.block(0,0,3,iNumCoils) - matHeadCoil.transpose();
363 VectorXd vecError = matDiff.colwise().norm();
364
365 // compute error
366 double dError = matDiff.colwise().norm().mean();;
367
368 return dError;
369}
370
371//=============================================================================================================
372
373Eigen::MatrixXd InvHpiFit::order(const std::vector<int>& vecOrder,
374 const Eigen::MatrixXd& matToOrder)
375{
376 const int iNumCoils = vecOrder.size();
377 MatrixXd matToOrderTemp = matToOrder;
378
379 for(int i = 0; i < iNumCoils; i++) {
380 matToOrderTemp.row(i) = matToOrder.row(vecOrder[i]);
381 }
382 return matToOrderTemp;
383}
384
385//=============================================================================================================
386
387QVector<int> InvHpiFit::order(const std::vector<int>& vecOrder,
388 const QVector<int>& vecToOrder)
389{
390 const int iNumCoils = vecOrder.size();
391 QVector<int> vecToOrderTemp = vecToOrder;
392
393 for(int i = 0; i < iNumCoils; i++) {
394 vecToOrderTemp[i] = vecToOrder[vecOrder[i]];
395 }
396 return vecToOrderTemp;
397}
398
399//=============================================================================================================
400
401Eigen::VectorXd InvHpiFit::computeGoF(const Eigen::VectorXd& vecDipFitError)
402{
403 VectorXd vecGoF(vecDipFitError.size());
404 for(int i = 0; i < vecDipFitError.size(); ++i) {
405 vecGoF(i) = 1 - vecDipFitError(i);
406 }
407 return vecGoF;
408}
409
410//=============================================================================================================
411
412FIFFLIB::FiffCoordTrans InvHpiFit::computeDeviceHeadTransformation(const Eigen::MatrixXd& matCoilsDev,
413 const Eigen::MatrixXd& matCoilsHead)
414{
415 const MatrixXd matTrans = computeTransformation(matCoilsHead,matCoilsDev);
416 return FiffCoordTrans(1,4,matTrans.cast<float>(),true);
417}
418
419//=============================================================================================================
420
421Eigen::Matrix4d InvHpiFit::computeTransformation(Eigen::MatrixXd matNH, MatrixXd matBT)
422{
423 MatrixXd matXdiff, matYdiff, matZdiff, matC, matQ;
424 Matrix4d matTransFinal = Matrix4d::Identity(4,4);
425 Matrix4d matRot = Matrix4d::Zero(4,4);
426 Matrix4d matTrans = Matrix4d::Identity(4,4);
427 double dMeanX,dMeanY,dMeanZ,dNormf;
428
429 for(int i = 0; i < 15; ++i) {
430 // Calculate mean translation for all points -> centroid of both data sets
431 matXdiff = matNH.col(0) - matBT.col(0);
432 matYdiff = matNH.col(1) - matBT.col(1);
433 matZdiff = matNH.col(2) - matBT.col(2);
434
435 dMeanX = matXdiff.mean();
436 dMeanY = matYdiff.mean();
437 dMeanZ = matZdiff.mean();
438
439 // Apply translation -> bring both data sets to the same center location
440 for (int j = 0; j < matBT.rows(); ++j) {
441 matBT(j,0) = matBT(j,0) + dMeanX;
442 matBT(j,1) = matBT(j,1) + dMeanY;
443 matBT(j,2) = matBT(j,2) + dMeanZ;
444 }
445
446 // Estimate rotation component
447 matC = matBT.transpose() * matNH;
448
449 JacobiSVD< MatrixXd > svd(matC ,Eigen::ComputeThinU | ComputeThinV);
450
451 matQ = svd.matrixU() * svd.matrixV().transpose();
452
453 //Handle special reflection case
454 if(matQ.determinant() < 0) {
455 matQ(0,2) = matQ(0,2) * -1;
456 matQ(1,2) = matQ(1,2) * -1;
457 matQ(2,2) = matQ(2,2) * -1;
458 }
459
460 // Apply rotation on translated points
461 matBT = matBT * matQ;
462
463 // Calculate GOF
464 dNormf = (matNH.transpose()-matBT.transpose()).norm();
465
466 // Store rotation part to transformation matrix
467 matRot(3,3) = 1;
468 for(int j = 0; j < 3; ++j) {
469 for(int k = 0; k < 3; ++k) {
470 matRot(j,k) = matQ(k,j);
471 }
472 }
473
474 // Store translation part to transformation matrix
475 matTrans(0,3) = dMeanX;
476 matTrans(1,3) = dMeanY;
477 matTrans(2,3) = dMeanZ;
478
479 // Safe rotation and translation to final matrix for next iteration step
480 // This step is safe to do since we change one of the input point sets (matBT)
481 // ToDo: Replace this for loop with a least square solution process
482 matTransFinal = matRot * matTrans * matTransFinal;
483 }
484 return matTransFinal;
485}
486
487//=============================================================================================================
488
489QVector<double> InvHpiFit::computeEstimationError(const Eigen::MatrixXd& matCoilsDev,
490 const Eigen::MatrixXd& matCoilsHead,
491 const FIFFLIB::FiffCoordTrans& transDevHead)
492{
493 //Calculate Error
494 MatrixXd matTemp = matCoilsDev;
495 MatrixXd matTestPos = transDevHead.apply_trans(matTemp.cast<float>()).cast<double>();
496 MatrixXd matDiffPos = matTestPos - matCoilsHead;
497
498 // compute error
499 int iNumCoils = matCoilsDev.rows();
500 QVector<double> vecError(iNumCoils);
501 for(int i = 0; i < matDiffPos.rows(); ++i) {
502 vecError[i] = matDiffPos.row(i).norm();
503 }
504 return vecError;
505}
506
507//=============================================================================================================
508
509FIFFLIB::FiffDigPointSet InvHpiFit::getFittedPointSet(const Eigen::MatrixXd& matCoilsDev)
510{
511 FiffDigPointSet fittedPointSet;
512 const int iNumCoils = matCoilsDev.rows();
513
514 for(int i = 0; i < iNumCoils; ++i) {
515 FiffDigPoint digPoint;
516 digPoint.kind = FIFFV_POINT_EEG; //Store as EEG so they have a different color
517 digPoint.ident = i;
518 digPoint.r[0] = matCoilsDev(i,0);
519 digPoint.r[1] = matCoilsDev(i,1);
520 digPoint.r[2] = matCoilsDev(i,2);
521
522 fittedPointSet << digPoint;
523 }
524 return fittedPointSet;
525}
526
527//=============================================================================================================
528
530 const Eigen::MatrixXf& transDevHead,
531 Eigen::MatrixXd& matPosition,
532 const Eigen::VectorXd& vecGoF,
533 const QVector<double>& vecError)
534
535{
536 Matrix3f matRot = transDevHead.block(0,0,3,3);
537
538 Eigen::Quaternionf quatHPI(matRot);
539 double dError = std::accumulate(vecError.begin(), vecError.end(), .0) / vecError.size(); // HPI estimation Error
540
541 matPosition.conservativeResize(matPosition.rows()+1, 10);
542 matPosition(matPosition.rows()-1,0) = fTime;
543 matPosition(matPosition.rows()-1,1) = quatHPI.x();
544 matPosition(matPosition.rows()-1,2) = quatHPI.y();
545 matPosition(matPosition.rows()-1,3) = quatHPI.z();
546 matPosition(matPosition.rows()-1,4) = transDevHead(0,3);
547 matPosition(matPosition.rows()-1,5) = transDevHead(1,3);
548 matPosition(matPosition.rows()-1,6) = transDevHead(2,3);
549 matPosition(matPosition.rows()-1,7) = vecGoF.mean();
550 matPosition(matPosition.rows()-1,8) = dError;
551 matPosition(matPosition.rows()-1,9) = 0;
552}
553
554//=============================================================================================================
555
556bool InvHpiFit::compareTransformation(const MatrixX4f& mDevHeadT,
557 const MatrixX4f& mDevHeadDest,
558 const float& fThreshRot,
559 const float& fThreshTrans)
560{
561 bool bState = false;
562
563 Matrix3f mRot = mDevHeadT.block(0,0,3,3);
564 Matrix3f mRotDest = mDevHeadDest.block(0,0,3,3);
565
566 VectorXf vTrans = mDevHeadT.block(0,3,3,1);
567 VectorXf vTransDest = mDevHeadDest.block(0,3,3,1);
568
569 Quaternionf quat(mRot);
570 Quaternionf quatNew(mRotDest);
571
572 // Compare Rotation
573 float fAngle = quat.angularDistance(quatNew);
574 fAngle = fAngle * 180 / M_PI;
575
576 // Compare translation
577 float fMove = (vTrans-vTransDest).norm();
578
579 // compare to thresholds and update
580 if(fMove > fThreshTrans) {
581 qInfo() << "Large movement: " << fMove*1000 << "mm";
582 bState = true;
583 } else if (fAngle > fThreshRot) {
584 qInfo() << "Large rotation: " << fAngle << "degree";
585 bState = true;
586 } else {
587 bState = false;
588 }
589
590 return bState;
591}
592
#define M_PI
IOUtils class declaration.
InvHpiModelParameters class declaration.
InvHpiFit class declaration.
InvHpiFitData class declaration.
InvHpiDataUpdater class declaration.
InvSignalModel class declaration.
InvSensorSet class declaration.
FwdCoilSet class declaration.
FiffDigPointSet class declaration.
FiffCov class declaration.
#define FIFFV_POINT_EEG
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:91
Eigen::MatrixX3f apply_inverse_trans(const Eigen::MatrixX3f &rr, bool do_move=true) const
Eigen::MatrixX3f apply_trans(const Eigen::MatrixX3f &rr, bool do_move=true) const
Eigen::Matrix< float, 4, 4, Eigen::DontAlign > trans
Estimated dipole parameters (position, moment, goodness-of-fit) for a single HPI coil.
Definition inv_hpi_fit.h:99
Eigen::MatrixXd pos
Eigen::VectorXd dpfiterror
Complete HPI fit output: per-coil dipole parameters, head-to-device transform, fit error,...
FIFFLIB::FiffCoordTrans devHeadTrans
QVector< int > hpiFreqs
QVector< double > errorDistances
Eigen::VectorXd GoF
FIFFLIB::FiffDigPointSet fittedCoils
static bool compareTransformation(const Eigen::MatrixX4f &mDevHeadT, const Eigen::MatrixX4f &mDevHeadDest, const float &fThreshRot, const float &fThreshTrans)
void fit(const Eigen::MatrixXd &matProjectedData, const Eigen::MatrixXd &matProjectors, const InvHpiModelParameters &hpiModelParameters, const Eigen::MatrixXd &matCoilsHead, HpiFitResult &hpiFitResult)
void checkForUpdate(const InvSensorSet &sensorSet)
static void storeHeadPosition(float fTime, const Eigen::MatrixXf &matTransDevHead, Eigen::MatrixXd &matPosition, const Eigen::VectorXd &vecGoF, const QVector< double > &vecError)
Eigen::RowVectorXd m_sensorData
Eigen::MatrixXd m_matProjector
Eigen::MatrixXd m_coilPos
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.
Generates the forward sinusoidal model matrix for HPI coil signals at known drive frequencies.