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