v2.0.0
Loading...
Searching...
No Matches
inv_rap_music.cpp
Go to the documentation of this file.
1//=============================================================================================================
35
36//=============================================================================================================
37// INCLUDES
38//=============================================================================================================
39
40#include "inv_rap_music.h"
41
42#include <math/numerics.h>
43
44#include <QDebug>
45
46#ifdef _OPENMP
47#include <omp.h>
48#endif
49
50#include <stdexcept>
51//=============================================================================================================
52// USED NAMESPACES
53//=============================================================================================================
54
55using namespace INVLIB;
56using namespace MNELIB;
57using namespace FIFFLIB;
58using namespace UTILSLIB;
59
60//=============================================================================================================
61// DEFINE MEMBER METHODS
62//=============================================================================================================
63
76
77//=============================================================================================================
78
79InvRapMusic::InvRapMusic(MNEForwardSolution& p_pFwd, bool p_bSparsed, int p_iN, double p_dThr)
80: m_iN(0)
81, m_dThreshold(0)
86, m_bIsInit(false)
88, m_fStcOverlap(-1)
89{
90 //Init
91 init(p_pFwd, p_bSparsed, p_iN, p_dThr);
92}
93
94//=============================================================================================================
95
99
100//=============================================================================================================
101
102bool InvRapMusic::init(MNEForwardSolution& p_pFwd, bool p_bSparsed, int p_iN, double p_dThr)
103{
104 //Get available thread number
105 #ifdef _OPENMP
106 qDebug() << "OpenMP enabled";
107 m_iMaxNumThreads = omp_get_max_threads();
108 #else
109 qDebug() << "OpenMP disabled (to enable it: VS2010->Project Properties->C/C++->Language, then modify OpenMP Support)";
111 #endif
112 qDebug() << "Available Threats: " << m_iMaxNumThreads;
113
114 //Initialize RAP MUSIC
115 qDebug() << "##### Initialization RAP MUSIC started ######\n\n";
116
117 m_iN = p_iN;
118 m_dThreshold = p_dThr;
119
120// //Grid check
121// if(p_pMatGrid != nullptr)
122// {
123// if ( p_pMatGrid->rows() != p_pMatLeadField->cols() / 3 )
124// {
125// qDebug() << "Grid does not fit to given Lead Field!\n";
126// return false;
127// }
128// }
129
130// m_pMatGrid = p_pMatGrid;
131
132 //Lead Field check
133 if ( p_pFwd.sol->data.cols() % 3 != 0 )
134 {
135 qDebug() << "Gain matrix is not associated with a 3D grid!\n";
136 return false;
137 }
138
139 m_iNumGridPoints = p_pFwd.sol->data.cols()/3;
140
141 m_iNumChannels = p_pFwd.sol->data.rows();
142
143// m_pMappedMatLeadField = new Eigen::Map<MatrixXT>
144// ( p_pMatLeadField->data(),
145// p_pMatLeadField->rows(),
146// p_pMatLeadField->cols() );
147
148 m_ForwardSolution = p_pFwd;
149
150 //##### Calc lead field combination #####
151
152 qDebug() << "Calculate gain matrix combinations. \n";
153
155
157
159
160 qDebug() << "Gain matrix combinations calculated. \n\n";
161
162 //##### Calc lead field combination end #####
163
164 qDebug() << "Number of grid points: " << m_iNumGridPoints << "\n\n";
165
166 qDebug() << "Number of combinated points: " << m_iNumLeadFieldCombinations << "\n\n";
167
168 qDebug() << "Number of sources to find: " << m_iN << "\n\n";
169
170 qDebug() << "Threshold: " << m_dThreshold << "\n\n";
171
172 //Init end
173
174 qDebug() << "##### Initialization RAP MUSIC completed ######\n\n\n";
175
176 Q_UNUSED(p_bSparsed);
177
178 m_bIsInit = true;
179
180 return m_bIsInit;
181}
182
183//=============================================================================================================
184
185const char* InvRapMusic::getName() const
186{
187 return "RAP MUSIC";
188}
189
190//=============================================================================================================
191
193{
194 return m_ForwardSolution.src;
195}
196
197//=============================================================================================================
198
199InvSourceEstimate InvRapMusic::calculateInverse(const FiffEvoked &p_fiffEvoked, bool pick_normal)
200{
201 Q_UNUSED(pick_normal);
202
203 InvSourceEstimate p_sourceEstimate;
204
205 if(p_fiffEvoked.data.rows() != m_iNumChannels)
206 {
207 qDebug() << "Number of FiffEvoked channels (" << p_fiffEvoked.data.rows() << ") doesn't match the number of channels (" << m_iNumChannels << ") of the forward solution.";
208 return p_sourceEstimate;
209 }
210// else
211// qDebug() << "Number of FiffEvoked channels (" << p_fiffEvoked.data.rows() << ") matchs the number of channels (" << m_iNumChannels << ") of the forward solution.";
212
213 //
214 // Rap MUSIC Source estimate
215 //
216 p_sourceEstimate.data = MatrixXd::Zero(m_ForwardSolution.nsource, p_fiffEvoked.data.cols());
217
218 //Results
219 p_sourceEstimate.vertices = VectorXi(m_ForwardSolution.src[0].vertno.size() + m_ForwardSolution.src[1].vertno.size());
220 p_sourceEstimate.vertices << m_ForwardSolution.src[0].vertno, m_ForwardSolution.src[1].vertno;
221
222 p_sourceEstimate.times = p_fiffEvoked.times;
223 p_sourceEstimate.tmin = p_fiffEvoked.times[0];
224 p_sourceEstimate.tstep = p_fiffEvoked.times[1] - p_fiffEvoked.times[0];
225
226 if(m_iSamplesStcWindow <= 3) //if samples per stc aren't set -> use full window
227 {
228 QList< InvDipolePair<double> > t_RapDipoles;
229 calculateInverse(p_fiffEvoked.data, t_RapDipoles);
230
231 for(qint32 i = 0; i < t_RapDipoles.size(); ++i)
232 {
233 double dip1 = sqrt( pow(t_RapDipoles[i].m_Dipole1.phi_x(),2) +
234 pow(t_RapDipoles[i].m_Dipole1.phi_y(),2) +
235 pow(t_RapDipoles[i].m_Dipole1.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
236
237 double dip2 = sqrt( pow(t_RapDipoles[i].m_Dipole2.phi_x(),2) +
238 pow(t_RapDipoles[i].m_Dipole2.phi_y(),2) +
239 pow(t_RapDipoles[i].m_Dipole2.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
240
241 RowVectorXd dip1Time = RowVectorXd::Constant(p_fiffEvoked.data.cols(), dip1);
242 RowVectorXd dip2Time = RowVectorXd::Constant(p_fiffEvoked.data.cols(), dip2);
243
244 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx1, 0, 1, p_fiffEvoked.data.cols()) = dip1Time;
245 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx2, 0, 1, p_fiffEvoked.data.cols()) = dip2Time;
246 }
247 }
248 else
249 {
250 bool first = true;
251 bool last = false;
252
253 qint32 t_iNumSensors = p_fiffEvoked.data.rows();
254 qint32 t_iNumSteps = p_fiffEvoked.data.cols();
255
256 qint32 t_iSamplesOverlap = static_cast<qint32>(floor(static_cast<float>(m_iSamplesStcWindow)*m_fStcOverlap));
257 qint32 t_iSamplesDiscard = t_iSamplesOverlap/2;
258
259 MatrixXd data = MatrixXd::Zero(t_iNumSensors, m_iSamplesStcWindow);
260
261 qint32 curSample = 0;
262 qint32 curResultSample = 0;
263 qint32 stcWindowSize = m_iSamplesStcWindow - 2*t_iSamplesDiscard;
264
265 while(!last)
266 {
267 QList< InvDipolePair<double> > t_RapDipoles;
268
269 //Data
270 if(curSample + m_iSamplesStcWindow >= t_iNumSteps) //last
271 {
272 last = true;
273 data = p_fiffEvoked.data.block(0, p_fiffEvoked.data.cols()-m_iSamplesStcWindow, t_iNumSensors, m_iSamplesStcWindow);
274 }
275 else
276 data = p_fiffEvoked.data.block(0, curSample, t_iNumSensors, m_iSamplesStcWindow);
277
278 curSample += (m_iSamplesStcWindow - t_iSamplesOverlap);
279 if(first)
280 curSample -= t_iSamplesDiscard; //shift on start t_iSamplesDiscard backwards
281
282 //Calculate
283 calculateInverse(data, t_RapDipoles);
284
285 //Assign Result
286 if(last)
287 stcWindowSize = p_sourceEstimate.data.cols() - curResultSample;
288
289 for(qint32 i = 0; i < t_RapDipoles.size(); ++i)
290 {
291 double dip1 = sqrt( pow(t_RapDipoles[i].m_Dipole1.phi_x(),2) +
292 pow(t_RapDipoles[i].m_Dipole1.phi_y(),2) +
293 pow(t_RapDipoles[i].m_Dipole1.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
294
295 double dip2 = sqrt( pow(t_RapDipoles[i].m_Dipole2.phi_x(),2) +
296 pow(t_RapDipoles[i].m_Dipole2.phi_y(),2) +
297 pow(t_RapDipoles[i].m_Dipole2.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
298
299 RowVectorXd dip1Time = RowVectorXd::Constant(stcWindowSize, dip1);
300 RowVectorXd dip2Time = RowVectorXd::Constant(stcWindowSize, dip2);
301
302 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx1, curResultSample, 1, stcWindowSize) = dip1Time;
303 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx2, curResultSample, 1, stcWindowSize) = dip2Time;
304
305 }
306
307 curResultSample += stcWindowSize;
308
309 if(first)
310 first = false;
311 }
312 }
313
314 return p_sourceEstimate;
315}
316
317//=============================================================================================================
318
319InvSourceEstimate InvRapMusic::calculateInverse(const MatrixXd &data, float tmin, float tstep, bool pick_normal) const
320{
321 Q_UNUSED(pick_normal);
322
323 InvSourceEstimate p_sourceEstimate;
324
325 if(data.rows() != m_iNumChannels)
326 {
327 qDebug() << "Number of FiffEvoked channels (" << data.rows() << ") doesn't match the number of channels (" << m_iNumChannels << ") of the forward solution.";
328 return p_sourceEstimate;
329 }
330// else
331// qDebug() << "Number of FiffEvoked channels (" << data.rows() << ") matchs the number of channels (" << m_iNumChannels << ") of the forward solution.";
332
333 //
334 // Rap MUSIC Source estimate
335 //
336 p_sourceEstimate.data = MatrixXd::Zero(m_ForwardSolution.nsource, data.cols());
337
338 //Results
339 p_sourceEstimate.vertices = VectorXi(m_ForwardSolution.src[0].vertno.size() + m_ForwardSolution.src[1].vertno.size());
340 p_sourceEstimate.vertices << m_ForwardSolution.src[0].vertno, m_ForwardSolution.src[1].vertno;
341
342 p_sourceEstimate.times = RowVectorXf::Zero(data.cols());
343 p_sourceEstimate.times[0] = tmin;
344 for(qint32 i = 1; i < p_sourceEstimate.times.size(); ++i)
345 p_sourceEstimate.times[i] = p_sourceEstimate.times[i-1] + tstep;
346 p_sourceEstimate.tmin = tmin;
347 p_sourceEstimate.tstep = tstep;
348
349 QList< InvDipolePair<double> > t_RapDipoles;
350 calculateInverse(data, t_RapDipoles);
351
352 for(qint32 i = 0; i < t_RapDipoles.size(); ++i)
353 {
354 double dip1 = sqrt( pow(t_RapDipoles[i].m_Dipole1.phi_x(),2) +
355 pow(t_RapDipoles[i].m_Dipole1.phi_y(),2) +
356 pow(t_RapDipoles[i].m_Dipole1.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
357
358 double dip2 = sqrt( pow(t_RapDipoles[i].m_Dipole2.phi_x(),2) +
359 pow(t_RapDipoles[i].m_Dipole2.phi_y(),2) +
360 pow(t_RapDipoles[i].m_Dipole2.phi_z(),2) ) * t_RapDipoles[i].m_vCorrelation;
361
362 RowVectorXd dip1Time = RowVectorXd::Constant(data.cols(), dip1);
363 RowVectorXd dip2Time = RowVectorXd::Constant(data.cols(), dip2);
364
365 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx1, 0, 1, data.cols()) = dip1Time;
366 p_sourceEstimate.data.block(t_RapDipoles[i].m_iIdx2, 0, 1, data.cols()) = dip2Time;
367 }
368
369 return p_sourceEstimate;
370}
371
372//=============================================================================================================
373
374InvSourceEstimate InvRapMusic::calculateInverse(const MatrixXd& p_matMeasurement, QList< InvDipolePair<double> > &p_RapDipoles) const
375{
376 InvSourceEstimate p_SourceEstimate;
377
378 //if not initialized -> break
379 if(!m_bIsInit)
380 {
381 throw std::logic_error("RAP MUSIC was not initialized");
382 }
383
384 //Test if data are correct
385 if(p_matMeasurement.rows() != m_iNumChannels)
386 {
387 throw std::invalid_argument("Lead field channels do not match number of measurement channels");
388 }
389
390 //Inits
391 //Stop the time for benchmark purpose
392 clock_t start, end;
393 start = clock();
394
395// //Map HPCMatrix to Eigen Matrix
396// Eigen::Map<MatrixXT>
397// t_MappedMatMeasurement( p_pMatMeasurement->data(),
398// p_pMatMeasurement->rows(),
399// p_pMatMeasurement->cols() );
400
401 //Calculate the signal subspace (t_pMatPhi_s)
402 MatrixXT* t_pMatPhi_s = nullptr;//(m_iNumChannels, m_iN < t_r ? m_iN : t_r);
403 int t_r = calcPhi_s(/*(MatrixXT)*/p_matMeasurement, t_pMatPhi_s);
404
405 int t_iMaxSearch = m_iN < t_r ? m_iN : t_r; //The smallest of Rank and Iterations
406
407 if (t_r < m_iN)
408 {
409 qDebug() << "Warning: Rank " << t_r << " of the measurement data is smaller than the " << m_iN;
410 qDebug() << " sources to find.";
411 qDebug() << " Searching now for " << t_iMaxSearch << " correlated sources.";
412 qDebug();
413 }
414
415 //Create Orthogonal Projector
416 //OrthProj
418 t_matOrthProj.setIdentity();
419
420 //A_k_1
421 MatrixXT t_matA_k_1(m_iNumChannels, t_iMaxSearch);
422 t_matA_k_1.setZero();
423
424// if (m_pMatGrid != nullptr)
425// {
426// if(p_pRapDipoles != nullptr)
427// p_pRapDipoles->initRapDipoles(m_pMatGrid);
428// else
429// p_pRapDipoles = new RapDipoles<T>(m_pMatGrid);
430// }
431// else
432// {
433// if(p_pRapDipoles != nullptr)
434// delete p_pRapDipoles;
435
436// p_pRapDipoles = new RapDipoles<T>();
437// }
438 p_RapDipoles.clear();
439
440 qDebug() << "##### Calculation of RAP MUSIC started ######\n\n";
441
442 MatrixXT t_matProj_Phi_s(t_matOrthProj.rows(), t_pMatPhi_s->cols());
443 //new Version: Calculate projection before
444 MatrixXT t_matProj_LeadField(m_ForwardSolution.sol->data.rows(), m_ForwardSolution.sol->data.cols());
445
446 for(int r = 0; r < t_iMaxSearch ; ++r)
447 {
448 t_matProj_Phi_s = t_matOrthProj*(*t_pMatPhi_s);
449
450 //new Version: Calculating Projection before
451 t_matProj_LeadField = t_matOrthProj * m_ForwardSolution.sol->data;//Subtract the found sources from the current found source
452
453 //###First Option###
454 //Step 1: lt. Mosher 1998 -> Maybe tmp_Proj_Phi_S is already orthogonal -> so no SVD needed -> U_B = tmp_Proj_Phi_S;
455 Eigen::JacobiSVD< MatrixXT > t_svdProj_Phi_S(t_matProj_Phi_s, Eigen::ComputeThinU);
456 MatrixXT t_matU_B;
457 useFullRank(t_svdProj_Phi_S.matrixU(), t_svdProj_Phi_S.singularValues().asDiagonal(), t_matU_B);
458
459 //Inits
461 t_vecRoh.setZero();
462
463 //subcorr benchmark
464 //Stop the time
465 clock_t start_subcorr, end_subcorr;
466 start_subcorr = clock();
467
468 //Multithreading correlation calculation
469 #ifdef _OPENMP
470 #pragma omp parallel num_threads(m_iMaxNumThreads)
471 #endif
472 {
473 #ifdef _OPENMP
474 #pragma omp for
475 #endif
476 for(int i = 0; i < m_iNumLeadFieldCombinations; i++)
477 {
478 //new Version: calculate matrix multiplication before
479 //Create Lead Field combinations -> It would be better to use a pointer construction, to increase performance
480 MatrixX6T t_matProj_G(t_matProj_LeadField.rows(),6);
481
482 int idx1 = m_ppPairIdxCombinations[i].x1;
483 int idx2 = m_ppPairIdxCombinations[i].x2;
484
485 InvRapMusic::getGainMatrixPair(t_matProj_LeadField, t_matProj_G, idx1, idx2);
486
487 t_vecRoh(i) = InvRapMusic::subcorr(t_matProj_G, t_matU_B);//t_vecRoh holds the correlations roh_k
488 }
489 }
490
491// if(r==0)
492// {
493// std::fstream filestr;
494// std::stringstream filename;
495// filename << "Roh_gold.txt";
496//
497// filestr.open ( filename.str().c_str(), std::fstream::out);
498// for(int i = 0; i < m_iNumLeadFieldCombinations; ++i)
499// {
500// filestr << t_vecRoh(i) << "\n";
501// }
502// filestr.close();
503//
504// //exit(0);
505// }
506
507 //subcorr benchmark
508 end_subcorr = clock();
509
510 float t_fSubcorrElapsedTime = ( static_cast<float>(end_subcorr-start_subcorr) / static_cast<float>(CLOCKS_PER_SEC) ) * 1000.0f;
511 qDebug() << "Time Elapsed: " << t_fSubcorrElapsedTime << " ms";
512
513 //Find the maximum of correlation - can't put this in the for loop because it's running in different threads.
514 double t_val_roh_k;
515
516 VectorXT::Index t_iMaxIdx;
517
518 t_val_roh_k = t_vecRoh.maxCoeff(&t_iMaxIdx);//p_vecCor = ^roh_k
519
520 //get positions in sparsed leadfield from index combinations;
521 int t_iIdx1 = m_ppPairIdxCombinations[t_iMaxIdx].x1;
522 int t_iIdx2 = m_ppPairIdxCombinations[t_iMaxIdx].x2;
523
524 // (Idx+1) because of MATLAB positions -> starting with 1 not with 0
525 qDebug() << "Iteration: " << r+1 << " of " << t_iMaxSearch
526 << "; Correlation: " << t_val_roh_k<< "; Position (Idx+1): " << t_iIdx1+1 << " - " << t_iIdx2+1 <<"\n\n";
527
528 //Calculations with the max correlated dipole pair G_k_1 -> ToDo Obsolet when taking direkt Projected Lead Field
529 MatrixX6T t_matG_k_1(m_ForwardSolution.sol->data.rows(),6);
530 InvRapMusic::getGainMatrixPair(m_ForwardSolution.sol->data, t_matG_k_1, t_iIdx1, t_iIdx2);
531
532 MatrixX6T t_matProj_G_k_1(t_matOrthProj.rows(), t_matG_k_1.cols());
533 t_matProj_G_k_1 = t_matOrthProj * t_matG_k_1;//Subtract the found sources from the current found source
534// MatrixX6T t_matProj_G_k_1(t_matProj_LeadField.rows(), 6);
535// getLeadFieldPair(t_matProj_LeadField, t_matProj_G_k_1, t_iIdx1, t_iIdx2);
536
537 //Calculate source direction
538 //source direction (p_pMatPhi) for current source r (phi_k_1)
539 Vector6T t_vec_phi_k_1(6);
540 InvRapMusic::subcorr(t_matProj_G_k_1, t_matU_B, t_vec_phi_k_1);//Correlate the current source to calculate the direction
541
542 //Set return values
543 InvRapMusic::insertSource(t_iIdx1, t_iIdx2, t_vec_phi_k_1, t_val_roh_k, p_RapDipoles);
544
545 //Stop Searching when Correlation is smaller then the Threshold
546 if (t_val_roh_k < m_dThreshold)
547 {
548 qDebug() << "Searching stopped, last correlation " << t_val_roh_k;
549 qDebug() << " is smaller then the given threshold " << m_dThreshold;
550 break;
551 }
552
553 //Calculate A_k_1 = [a_theta_1..a_theta_k_1] matrix for subtraction of found source
554 InvRapMusic::calcA_k_1(t_matG_k_1, t_vec_phi_k_1, r, t_matA_k_1);
555
556 //Calculate new orthogonal Projector (Pi_k_1)
557 calcOrthProj(t_matA_k_1, t_matOrthProj);
558
559 //garbage collecting
560 //ToDo
561 }
562
563 qDebug() << "##### Calculation of RAP MUSIC completed ######";
564
565 end = clock();
566
567 float t_fElapsedTime = ( static_cast<float>(end-start) / static_cast<float>(CLOCKS_PER_SEC) ) * 1000.0f;
568 qDebug() << "Total Time Elapsed: " << t_fElapsedTime << " ms";
569
570 //garbage collecting
571 delete t_pMatPhi_s;
572
573 return p_SourceEstimate;
574}
575
576//=============================================================================================================
577
578int InvRapMusic::calcPhi_s(const MatrixXT& p_matMeasurement, MatrixXT* &p_pMatPhi_s) const
579{
580 //Calculate p_pMatPhi_s
581 MatrixXT t_matF;//t_matF = makeSquareMat(p_pMatMeasurement); //FF^T -> ToDo Check this
582 if (p_matMeasurement.cols() > p_matMeasurement.rows())
583 t_matF = makeSquareMat(p_matMeasurement); //FF^T
584 else
585 t_matF = MatrixXT(p_matMeasurement);
586
587 Eigen::JacobiSVD<MatrixXT> t_svdF(t_matF, Eigen::ComputeThinU);
588
589 int t_r = getRank(t_svdF.singularValues().asDiagonal());
590
591 int t_iCols = t_r;//t_r < m_iN ? m_iN : t_r;
592
593 if (p_pMatPhi_s != nullptr)
594 delete p_pMatPhi_s;
595
596 //m_iNumChannels has to be equal to t_svdF.matrixU().rows()
597 p_pMatPhi_s = new MatrixXT(m_iNumChannels, t_iCols);
598
599 //assign the signal subspace
600 memcpy(p_pMatPhi_s->data(), t_svdF.matrixU().data(), sizeof(double) * m_iNumChannels * t_iCols);
601
602 return t_r;
603}
604
605//=============================================================================================================
606
607double InvRapMusic::subcorr(MatrixX6T& p_matProj_G, const MatrixXT& p_matU_B)
608{
609 //Orthogonalisierungstest wegen performance weggelassen -> ohne is es viel schneller
610
611 Matrix6T t_matSigma_A(6, 6);
612 Matrix6XT t_matU_A_T(6, p_matProj_G.rows()); //rows and cols are changed, because of CV_SVD_U_T
613
614 Eigen::JacobiSVD<MatrixXT> t_svdProj_G(p_matProj_G, Eigen::ComputeThinU);
615
616 t_matSigma_A = t_svdProj_G.singularValues().asDiagonal();
617 t_matU_A_T = t_svdProj_G.matrixU().transpose();
618
619 //lt. Mosher 1998 ToDo: Only Retain those Components of U_A and U_B that correspond to nonzero singular values
620 //for U_A and U_B the number of columns corresponds to their ranks
621 MatrixXT t_matU_A_T_full;
622 //reduce to rank only when directions aren't calculated, otherwise use the full t_matU_A_T
623
624 useFullRank(t_matU_A_T, t_matSigma_A, t_matU_A_T_full, IS_TRANSPOSED);//lt. Mosher 1998: Only Retain those Components of U_A that correspond to nonzero singular values -> for U_A the number of columns corresponds to their ranks
625
626 MatrixXT t_matCor(t_matU_A_T_full.rows(), p_matU_B.cols());
627
628 //Step 2: compute the subspace correlation
629 t_matCor = t_matU_A_T_full*p_matU_B;//lt. Mosher 1998: C = U_A^T * U_B
630
631 VectorXT t_vecSigma_C;
632
633 if (t_matCor.cols() > t_matCor.rows())
634 {
635 MatrixXT t_matCor_H = t_matCor.adjoint(); //for complex it has to be adjoint
636
637 Eigen::JacobiSVD<MatrixXT> t_svdCor_H(t_matCor_H);
638
639 t_vecSigma_C = t_svdCor_H.singularValues();
640 }
641 else
642 {
643 Eigen::JacobiSVD<MatrixXT> t_svdCor(t_matCor);
644
645 t_vecSigma_C = t_svdCor.singularValues();
646 }
647
648 //Step 3
649 double t_dRetSigma_C;
650 t_dRetSigma_C = t_vecSigma_C(0); //Take only the correlation of the first principal components
651
652 //garbage collecting
653 //ToDo
654
655 return t_dRetSigma_C;
656}
657
658//=============================================================================================================
659
660double InvRapMusic::subcorr(MatrixX6T& p_matProj_G, const MatrixXT& p_matU_B, Vector6T& p_vec_phi_k_1)
661{
662 //Orthogonalisierungstest wegen performance weggelassen -> ohne is es viel schneller
663
664 Matrix6T sigma_A(6, 6);
665 Matrix6XT U_A_T(6, p_matProj_G.rows()); //rows and cols are changed, because of CV_SVD_U_T
666 Matrix6T V_A(6, 6);
667
668 Eigen::JacobiSVD<MatrixXT> svdOfProj_G(p_matProj_G, Eigen::ComputeThinU | Eigen::ComputeThinV);
669
670 sigma_A = svdOfProj_G.singularValues().asDiagonal();
671 U_A_T = svdOfProj_G.matrixU().transpose();
672 V_A = svdOfProj_G.matrixV();
673
674 //lt. Mosher 1998 ToDo: Only Retain those Components of U_A and U_B that correspond to nonzero singular values
675 //for U_A and U_B the number of columns corresponds to their ranks
676 //-> reduce to rank only when directions aren't calculated, otherwise use the full U_A_T
677
678 Matrix6XT t_matCor(6, p_matU_B.cols());
679
680 //Step 2: compute the subspace correlation
681 t_matCor = U_A_T*p_matU_B;//lt. Mosher 1998: C = U_A^T * U_B
682
683 VectorXT sigma_C;
684
685 //Step 4
686 Matrix6XT U_C;
687
688 if (t_matCor.cols() > t_matCor.rows())
689 {
690 MatrixX6T Cor_H(t_matCor.cols(), 6);
691 Cor_H = t_matCor.adjoint(); //for complex it has to be adjunct
692
693 Eigen::JacobiSVD<MatrixXT> svdOfCor_H(Cor_H, Eigen::ComputeThinV);
694
695 U_C = svdOfCor_H.matrixV(); //because t_matCor Hermitesch U and V are exchanged
696 sigma_C = svdOfCor_H.singularValues();
697 }
698 else
699 {
700 Eigen::JacobiSVD<MatrixXT> svdOfCor(t_matCor, Eigen::ComputeThinU);
701
702 U_C = svdOfCor.matrixU();
703 sigma_C = svdOfCor.singularValues();
704 }
705
706 Matrix6T sigma_a_inv;
707 sigma_a_inv = sigma_A.inverse();
708
709 Matrix6XT X;
710 X = (V_A*sigma_a_inv)*U_C;//X = V_A*Sigma_A^-1*U_C
711
712 Vector6T X_max;//only for the maximum c - so instead of X->cols use 1
713 X_max = X.col(0);
714
715 double norm_X = 1/(X_max.norm());
716
717 //Multiply a scalar with an Array -> linear transform
718 p_vec_phi_k_1 = X_max*norm_X;//u1 = x1/||x1|| this is the orientation
719
720 //garbage collecting
721 //ToDo
722
723 //Step 3
724 double ret_sigma_C;
725 ret_sigma_C = sigma_C(0); //Take only the correlation of the first principal components
726
727 //garbage collecting
728 //ToDo
729
730 return ret_sigma_C;
731}
732
733//=============================================================================================================
734
735void InvRapMusic::calcA_k_1( const MatrixX6T& p_matG_k_1,
736 const Vector6T& p_matPhi_k_1,
737 const int p_iIdxk_1,
738 MatrixXT& p_matA_k_1)
739{
740 //Calculate A_k_1 = [a_theta_1..a_theta_k_1] matrix for subtraction of found source
741 VectorXT t_vec_a_theta_k_1(p_matG_k_1.rows(),1);
742
743 t_vec_a_theta_k_1 = p_matG_k_1*p_matPhi_k_1; // a_theta_k_1 = G_k_1*phi_k_1 this corresponds to the normalized signal component in subspace r
744
745 p_matA_k_1.block(0,p_iIdxk_1,p_matA_k_1.rows(),1) = t_vec_a_theta_k_1;
746}
747
748//=============================================================================================================
749
750void InvRapMusic::calcOrthProj(const MatrixXT& p_matA_k_1, MatrixXT& p_matOrthProj) const
751{
752 //Calculate OrthProj=I-A_k_1*(A_k_1'*A_k_1)^-1*A_k_1' //Wetterling -> A_k_1 = Gain
753
754 MatrixXT t_matA_k_1_tmp(p_matA_k_1.cols(), p_matA_k_1.cols());
755 t_matA_k_1_tmp = p_matA_k_1.adjoint()*p_matA_k_1;//A_k_1'*A_k_1 = A_k_1_tmp -> A_k_1' has to be adjoint for complex
756
757 int t_size = t_matA_k_1_tmp.cols();
758
759 while (!t_matA_k_1_tmp.block(0,0,t_size,t_size).fullPivLu().isInvertible())
760 {
761 --t_size;
762 }
763
764 MatrixXT t_matA_k_1_tmp_inv(t_matA_k_1_tmp.rows(), t_matA_k_1_tmp.cols());
765 t_matA_k_1_tmp_inv.setZero();
766
767 t_matA_k_1_tmp_inv.block(0,0,t_size,t_size) = t_matA_k_1_tmp.block(0,0,t_size,t_size).inverse();//(A_k_1_tmp)^-1 = A_k_1_tmp_inv
768
769 t_matA_k_1_tmp = MatrixXT::Zero(p_matA_k_1.rows(), p_matA_k_1.cols());
770
771 t_matA_k_1_tmp = p_matA_k_1*t_matA_k_1_tmp_inv;//(A_k_1*A_k_1_tmp_inv) = A_k_1_tmp
772
773 MatrixXT t_matA_k_1_tmp2(p_matA_k_1.rows(), p_matA_k_1.rows());
774
775 t_matA_k_1_tmp2 = t_matA_k_1_tmp*p_matA_k_1.adjoint();//(A_k_1_tmp)*A_k_1' -> here A_k_1' is only transposed - it has to be adjoint
776
778 I.setIdentity();
779
780 p_matOrthProj = I-t_matA_k_1_tmp2; //OrthProj=I-A_k_1*(A_k_1'*A_k_1)^-1*A_k_1';
781
782 //garbage collecting
783 //ToDo
784}
785
786//=============================================================================================================
787
788void InvRapMusic::calcPairCombinations( const int p_iNumPoints,
789 const int p_iNumCombinations,
790 std::vector<Pair>& p_pairIdxCombinations) const
791{
792 int idx1 = 0;
793 int idx2 = 0;
794
795 //Process Code in {m_max_num_threads} threads -> When compile with Intel Compiler -> probably obsolete
796 #ifdef _OPENMP
797 #pragma omp parallel num_threads(m_iMaxNumThreads) private(idx1, idx2)
798 #endif
799 {
800 #ifdef _OPENMP
801 #pragma omp for
802 #endif
803 for (int i = 0; i < p_iNumCombinations; ++i)
804 {
805 InvRapMusic::getPointPair(p_iNumPoints, i, idx1, idx2);
806
807 Pair t_pairCombination;
808 t_pairCombination.x1 = idx1;
809 t_pairCombination.x2 = idx2;
810
811 p_pairIdxCombinations[i] = t_pairCombination;
812 }
813 }
814}
815
816//=============================================================================================================
817
818void InvRapMusic::getPointPair(const int p_iPoints, const int p_iCurIdx, int &p_iIdx1, int &p_iIdx2)
819{
820 int ii = p_iPoints*(p_iPoints+1)/2-1-p_iCurIdx;
821 int K = static_cast<int>(floor((sqrt(static_cast<double>(8*ii+1))-1)/2));
822
823 p_iIdx1 = p_iPoints-1-K;
824
825 p_iIdx2 = (p_iCurIdx-p_iPoints*(p_iPoints+1)/2 + (K+1)*(K+2)/2)+p_iIdx1;
826}
827
828//=============================================================================================================
829//ToDo don't make a real copy
830void InvRapMusic::getGainMatrixPair( const MatrixXT& p_matGainMarix,
831 MatrixX6T& p_matGainMarix_Pair,
832 int p_iIdx1, int p_iIdx2)
833{
834 p_matGainMarix_Pair.block(0,0,p_matGainMarix.rows(),3) = p_matGainMarix.block(0, p_iIdx1*3, p_matGainMarix.rows(), 3);
835
836 p_matGainMarix_Pair.block(0,3,p_matGainMarix.rows(),3) = p_matGainMarix.block(0, p_iIdx2*3, p_matGainMarix.rows(), 3);
837}
838
839//=============================================================================================================
840
841void InvRapMusic::insertSource( int p_iDipoleIdx1, int p_iDipoleIdx2,
842 const Vector6T &p_vec_phi_k_1,
843 double p_valCor,
844 QList< InvDipolePair<double> > &p_RapDipoles)
845{
846 InvDipolePair<double> t_pRapDipolePair;
847
848 t_pRapDipolePair.m_iIdx1 = p_iDipoleIdx1; //p_iDipoleIdx1+1 because of MATLAB index
849 t_pRapDipolePair.m_iIdx2 = p_iDipoleIdx2;
850
851 t_pRapDipolePair.m_Dipole1.x() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx1, 0) : 0;
852 t_pRapDipolePair.m_Dipole1.y() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx1, 1) : 0;
853 t_pRapDipolePair.m_Dipole1.z() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx1, 2) : 0;
854
855 t_pRapDipolePair.m_Dipole2.x() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx2, 0) : 0;
856 t_pRapDipolePair.m_Dipole2.y() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx2, 1) : 0;
857 t_pRapDipolePair.m_Dipole2.z() = 0; //m_bGridInitialized ? (*m_pMatGrid)(p_iDipoleIdx2, 2) : 0;
858
859 t_pRapDipolePair.m_Dipole1.phi_x() = p_vec_phi_k_1[0];
860 t_pRapDipolePair.m_Dipole1.phi_y() = p_vec_phi_k_1[1];
861 t_pRapDipolePair.m_Dipole1.phi_z() = p_vec_phi_k_1[2];
862
863 t_pRapDipolePair.m_Dipole2.phi_x() = p_vec_phi_k_1[3];
864 t_pRapDipolePair.m_Dipole2.phi_y() = p_vec_phi_k_1[4];
865 t_pRapDipolePair.m_Dipole2.phi_z() = p_vec_phi_k_1[5];
866
867 t_pRapDipolePair.m_vCorrelation = p_valCor;
868
869 p_RapDipoles.append(t_pRapDipolePair);
870}
871
872//=============================================================================================================
873
874void InvRapMusic::setStcAttr(int p_iSampStcWin, float p_fStcOverlap)
875{
876 m_iSamplesStcWindow = p_iSampStcWin;
877 m_fStcOverlap = p_fStcOverlap;
878}
Numerics class declaration.
InvRapMusic algorithm class declaration.
#define IS_TRANSPOSED
constexpr int X
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Eigen::RowVectorXf times
Eigen::MatrixXd data
Pair of correlated dipole indices and orientations found by the RAP MUSIC scanning step.
Definition inv_dipole.h:59
InvDipole< T > m_Dipole1
Definition inv_dipole.h:61
InvDipole< T > m_Dipole2
Definition inv_dipole.h:64
Index pair representing two grid points evaluated together in the RAP MUSIC subspace scan.
Eigen::Matrix< double, 6, 1 > Vector6T
void calcPairCombinations(const int p_iNumPoints, const int p_iNumCombinations, std::vector< Pair > &p_pairIdxCombinations) const
static void calcA_k_1(const MatrixX6T &p_matG_k_1, const Vector6T &p_matPhi_k_1, const int p_iIdxk_1, MatrixXT &p_matA_k_1)
static void getPointPair(const int p_iPoints, const int p_iCurIdx, int &p_iIdx1, int &p_iIdx2)
virtual InvSourceEstimate calculateInverse(const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
static int getRank(const MatrixXT &p_matSigma)
virtual const char * getName() const
static MatrixXT makeSquareMat(const MatrixXT &p_matF)
static void insertSource(int p_iDipoleIdx1, int p_iDipoleIdx2, const Vector6T &p_vec_phi_k_1, double p_valCor, QList< InvDipolePair< double > > &p_RapDipoles)
void setStcAttr(int p_iSampStcWin, float p_fStcOverlap)
bool init(MNELIB::MNEForwardSolution &p_pFwd, bool p_bSparsed=false, int p_iN=2, double p_dThr=0.5)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixXT
Eigen::Matrix< double, 6, 6 > Matrix6T
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorXT
void calcOrthProj(const MatrixXT &p_matA_k_1, MatrixXT &p_matOrthProj) const
int calcPhi_s(const MatrixXT &p_matMeasurement, MatrixXT *&p_pMatPhi_s) const
static int useFullRank(const MatrixXT &p_Mat, const MatrixXT &p_matSigma_src, MatrixXT &p_matFull_Rank, int type=NOT_TRANSPOSED)
Eigen::Matrix< double, Eigen::Dynamic, 6 > MatrixX6T
MNELIB::MNEForwardSolution m_ForwardSolution
virtual const MNELIB::MNESourceSpaces & getSourceSpace() const
std::vector< Pair > m_ppPairIdxCombinations
Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6XT
static void getGainMatrixPair(const MatrixXT &p_matGainMarix, MatrixX6T &p_matGainMarix_Pair, int p_iIdx1, int p_iIdx2)
static double subcorr(MatrixX6T &p_matProj_G, const MatrixXT &p_pMatU_B)
static int nchoose2(int n)
Definition numerics.cpp:92
MNELIB::MNESourceSpaces src
FIFFLIB::FiffNamedMatrix::SDPtr sol
Source Space descritpion.