v2.0.0
Loading...
Searching...
No Matches
pwlrapmusic.cpp
Go to the documentation of this file.
1//=============================================================================================================
35
36//=============================================================================================================
37// INCLUDES
38//=============================================================================================================
39
40#include "pwlrapmusic.h"
41
42#ifdef _OPENMP
43#include <omp.h>
44#endif
45
46//=============================================================================================================
47// STL INCLUDES
48//=============================================================================================================
49
50#include <iostream>
51
52//=============================================================================================================
53// USED NAMESPACES
54//=============================================================================================================
55
56using namespace INVERSELIB;
57using namespace MNELIB;
58using namespace FIFFLIB;
59
60//=============================================================================================================
61// DEFINE MEMBER METHODS
62//=============================================================================================================
63
68
69//=============================================================================================================
70
71PwlRapMusic::PwlRapMusic(MNEForwardSolution& p_pFwd, bool p_bSparsed, int p_iN, double p_dThr)
72: RapMusic(p_pFwd, p_bSparsed, p_iN, p_dThr)
73{
74 //Init
75 init(p_pFwd, p_bSparsed, p_iN, p_dThr);
76}
77
78//=============================================================================================================
79
83
84//=============================================================================================================
85
86const char* PwlRapMusic::getName() const
87{
88 return "Powell RAP MUSIC";
89}
90
91//=============================================================================================================
92
93MNESourceEstimate PwlRapMusic::calculateInverse(const FiffEvoked &p_fiffEvoked, bool pick_normal)
94{
95 return RapMusic::calculateInverse(p_fiffEvoked, pick_normal);
96}
97
98//=============================================================================================================
99
100MNESourceEstimate PwlRapMusic::calculateInverse(const MatrixXd &data, float tmin, float tstep) const
101{
102 return RapMusic::calculateInverse(data, tmin, tstep);
103}
104
105//=============================================================================================================
106
107MNESourceEstimate PwlRapMusic::calculateInverse(const MatrixXd& p_matMeasurement, QList< DipolePair<double> > &p_RapDipoles) const
108{
109 MNESourceEstimate p_SourceEstimate;
110
111 //if not initialized -> break
112 if(!m_bIsInit)
113 {
114 std::cout << "RAP MUSIC wasn't initialized!"; //ToDo: catch this earlier
115 return p_SourceEstimate; //false
116 }
117
118 //Test if data are correct
119 if(p_matMeasurement.rows() != m_iNumChannels)
120 {
121 std::cout << "Lead Field channels do not fit to number of measurement channels!"; //ToDo: catch this earlier
122 return p_SourceEstimate;
123 }
124
125 //Inits
126 //Stop the time for benchmark purpose
127 clock_t start, end;
128 start = clock();
129
130// //Map HPCMatrix to Eigen Matrix
131// Eigen::Map<MatrixXT>
132// t_MappedMatMeasurement( p_pMatMeasurement->data(),
133// p_pMatMeasurement->rows(),
134// p_pMatMeasurement->cols() );
135
136 //Calculate the signal subspace (t_pMatPhi_s)
137 MatrixXT* t_pMatPhi_s = NULL;//(m_iNumChannels, m_iN < t_r ? m_iN : t_r);
138 int t_r = calcPhi_s(/*(MatrixXT)*/p_matMeasurement, t_pMatPhi_s);
139
140 int t_iMaxSearch = m_iN < t_r ? m_iN : t_r; //The smallest of Rank and Iterations
141
142 if (t_r < m_iN)
143 {
144 std::cout << "Warning: Rank " << t_r << " of the measurement data is smaller than the " << m_iN;
145 std::cout << " sources to find." << std::endl;
146 std::cout << " Searching now for " << t_iMaxSearch << " correlated sources.";
147 std::cout << std::endl << std::endl;
148 }
149
150 //Create Orthogonal Projector
151 //OrthProj
153 t_matOrthProj.setIdentity();
154
155 //A_k_1
156 MatrixXT t_matA_k_1(m_iNumChannels, t_iMaxSearch);
157 t_matA_k_1.setZero();
158
159// if (m_pMatGrid != NULL)
160// {
161// if(p_pRapDipoles != NULL)
162// p_pRapDipoles->initRapDipoles(m_pMatGrid);
163// else
164// p_pRapDipoles = new RapDipoles<T>(m_pMatGrid);
165// }
166// else
167// {
168// if(p_pRapDipoles != NULL)
169// delete p_pRapDipoles;
170
171// p_pRapDipoles = new RapDipoles<T>();
172// }
173 p_RapDipoles.clear();
174
175 std::cout << "##### Calculation of PWL RAP MUSIC started ######\n\n";
176
177 MatrixXT t_matProj_Phi_s(t_matOrthProj.rows(), t_pMatPhi_s->cols());
178 //new Version: Calculate projection before
179 MatrixXT t_matProj_LeadField(m_ForwardSolution.sol->data.rows(), m_ForwardSolution.sol->data.cols());
180
181 for(int r = 0; r < t_iMaxSearch ; ++r)
182 {
183 t_matProj_Phi_s = t_matOrthProj*(*t_pMatPhi_s);
184
185 //new Version: Calculating Projection before
186 t_matProj_LeadField = t_matOrthProj * m_ForwardSolution.sol->data;//Subtract the found sources from the current found source
187
188 //###First Option###
189 //Step 1: lt. Mosher 1998 -> Maybe tmp_Proj_Phi_S is already orthogonal -> so no SVD needed -> U_B = tmp_Proj_Phi_S;
190 Eigen::JacobiSVD< MatrixXT > t_svdProj_Phi_S(t_matProj_Phi_s, Eigen::ComputeThinU);
191 MatrixXT t_matU_B;
192 useFullRank(t_svdProj_Phi_S.matrixU(), t_svdProj_Phi_S.singularValues().asDiagonal(), t_matU_B);
193
194 //Inits
196 t_vecRoh.setZero();
197
198 //subcorr benchmark
199 //Stop the time
200 clock_t start_subcorr, end_subcorr;
201 start_subcorr = clock();
202
203 double t_val_roh_k;
204
205 //Powell
206 int t_iCurrentRow = 2;
207
208 int t_iIdx1 = -1;
209 int t_iIdx2 = -1;
210
211 int t_iMaxIdx_old = -1;
212
213 int t_iMaxFound = 0;
214
215 Eigen::VectorXi t_pVecIdxElements(m_iNumGridPoints);
216
217 PowellIdxVec(t_iCurrentRow, m_iNumGridPoints, t_pVecIdxElements);
218
219 int t_iNumVecElements = m_iNumGridPoints;
220
221 while(t_iMaxFound == 0)
222 {
223
224 //Multithreading correlation calculation
225 #ifdef _OPENMP
226 #pragma omp parallel num_threads(m_iMaxNumThreads)
227 #endif
228 {
229 #ifdef _OPENMP
230 #pragma omp for
231 #endif
232 for(int i = 0; i < t_iNumVecElements; i++)
233 {
234 int k = t_pVecIdxElements(i);
235 //new Version: calculate matrix multiplication before
236 //Create Lead Field combinations -> It would be better to use a pointer construction, to increase performance
237 MatrixX6T t_matProj_G(t_matProj_LeadField.rows(),6);
238
239 int idx1 = m_ppPairIdxCombinations[k]->x1;
240 int idx2 = m_ppPairIdxCombinations[k]->x2;
241
242 RapMusic::getGainMatrixPair(t_matProj_LeadField, t_matProj_G, idx1, idx2);
243
244 t_vecRoh(k) = RapMusic::subcorr(t_matProj_G, t_matU_B);//t_vecRoh holds the correlations roh_k
245 }
246 }
247
248 // if(r==0)
249 // {
250 // std::fstream filestr;
251 // std::stringstream filename;
252 // filename << "Roh_gold.txt";
253 //
254 // filestr.open ( filename.str().c_str(), std::fstream::out);
255 // for(int i = 0; i < m_iNumLeadFieldCombinations; ++i)
256 // {
257 // filestr << t_vecRoh(i) << "\n";
258 // }
259 // filestr.close();
260 //
261 // //exit(0);
262 // }
263
264 //Find the maximum of correlation - can't put this in the for loop because it's running in different threads.
265
266 VectorXT::Index t_iMaxIdx;
267
268 t_val_roh_k = t_vecRoh.maxCoeff(&t_iMaxIdx);//p_vecCor = ^roh_k
269
270 if((int)t_iMaxIdx == t_iMaxIdx_old)
271 {
272 t_iMaxFound = 1;
273 break;
274 }
275 else
276 {
277 t_iMaxIdx_old = t_iMaxIdx;
278 //get positions in sparsed leadfield from index combinations;
279 t_iIdx1 = m_ppPairIdxCombinations[t_iMaxIdx]->x1;
280 t_iIdx2 = m_ppPairIdxCombinations[t_iMaxIdx]->x2;
281 }
282
283 //set new index
284 if(t_iIdx1 == t_iCurrentRow)
285 t_iCurrentRow = t_iIdx2;
286 else
287 t_iCurrentRow = t_iIdx1;
288
289 PowellIdxVec(t_iCurrentRow, m_iNumGridPoints, t_pVecIdxElements);
290 }
291
292 //subcorr benchmark
293 end_subcorr = clock();
294
295 float t_fSubcorrElapsedTime = ( (float)(end_subcorr-start_subcorr) / (float)CLOCKS_PER_SEC ) * 1000.0f;
296 std::cout << "Time Elapsed: " << t_fSubcorrElapsedTime << " ms" << std::endl;
297
298 // (Idx+1) because of MATLAB positions -> starting with 1 not with 0
299 std::cout << "Iteration: " << r+1 << " of " << t_iMaxSearch
300 << "; Correlation: " << t_val_roh_k<< "; Position (Idx+1): " << t_iIdx1+1 << " - " << t_iIdx2+1 <<"\n\n";
301
302 //Calculations with the max correlated dipole pair G_k_1
303 MatrixX6T t_matG_k_1(m_ForwardSolution.sol->data.rows(),6);
304 RapMusic::getGainMatrixPair(m_ForwardSolution.sol->data, t_matG_k_1, t_iIdx1, t_iIdx2);
305
306 MatrixX6T t_matProj_G_k_1(t_matOrthProj.rows(), t_matG_k_1.cols());
307 t_matProj_G_k_1 = t_matOrthProj * t_matG_k_1;//Subtract the found sources from the current found source
308// MatrixX6T t_matProj_G_k_1(t_matProj_LeadField.rows(), 6);
309// getLeadFieldPair(t_matProj_LeadField, t_matProj_G_k_1, t_iIdx1, t_iIdx2);
310
311 //Calculate source direction
312 //source direction (p_pMatPhi) for current source r (phi_k_1)
313 Vector6T t_vec_phi_k_1(6, 1);
314 RapMusic::subcorr(t_matProj_G_k_1, t_matU_B, t_vec_phi_k_1);//Correlate the current source to calculate the direction
315
316 //Set return values
317 RapMusic::insertSource(t_iIdx1, t_iIdx2, t_vec_phi_k_1, t_val_roh_k, p_RapDipoles);
318
319 //Stop Searching when Correlation is smaller then the Threshold
320 if (t_val_roh_k < m_dThreshold)
321 {
322 std::cout << "Searching stopped, last correlation " << t_val_roh_k;
323 std::cout << " is smaller then the given threshold " << m_dThreshold << std::endl << std::endl;
324 break;
325 }
326
327 //Calculate A_k_1 = [a_theta_1..a_theta_k_1] matrix for subtraction of found source
328 RapMusic::calcA_k_1(t_matG_k_1, t_vec_phi_k_1, r, t_matA_k_1);
329
330 //Calculate new orthogonal Projector (Pi_k_1)
331 calcOrthProj(t_matA_k_1, t_matOrthProj);
332
333 //garbage collecting
334 //ToDo
335 }
336
337 std::cout << "##### Calculation of PWL RAP MUSIC completed ######"<< std::endl << std::endl << std::endl;
338
339 end = clock();
340
341 float t_fElapsedTime = ( (float)(end-start) / (float)CLOCKS_PER_SEC ) * 1000.0f;
342 std::cout << "Total Time Elapsed: " << t_fElapsedTime << " ms" << std::endl << std::endl;
343
344 //garbage collecting
345 delete t_pMatPhi_s;
346
347 return p_SourceEstimate;
348}
349
350//=============================================================================================================
351
352int PwlRapMusic::PowellOffset(int p_iRow, int p_iNumPoints)
353{
354
355 return p_iRow*p_iNumPoints - (( (p_iRow-1)*p_iRow) / 2); //triangular series 1 3 6 10 ... = (num_pairs*(num_pairs+1))/2
356}
357
358//=============================================================================================================
359
360void PwlRapMusic::PowellIdxVec(int p_iRow, int p_iNumPoints, Eigen::VectorXi& p_pVecElements)
361{
362
363 // if(p_pVecElements != NULL)
364 // delete[] p_pVecElements;
365 //
366 // p_pVecElements = new int(p_iNumPoints);
367
368 //col combination index
369 for(int i = 0; i <= p_iRow; ++i)//=p_iNumPoints-1
370 p_pVecElements(i) = PwlRapMusic::PowellOffset(i+1,p_iNumPoints)-(p_iNumPoints-p_iRow);
371
372 //row combination index
373 int off = PwlRapMusic::PowellOffset(p_iRow,p_iNumPoints);
374 int length = p_iNumPoints - p_iRow;
375 int k=0;
376 for(int i = p_iRow; i < p_iRow+length; ++i)//=p_iNumPoints-1
377 {
378 p_pVecElements(i) = off+k;
379 k = k + 1;
380 }
381}
PwlRapMusic algorithm class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Pair of correlated dipole indices and orientations found by the RAP MUSIC scanning step.
Definition dipole.h:59
virtual const char * getName() const
static void PowellIdxVec(int p_iRow, int p_iNumPoints, Eigen::VectorXi &p_pVecElements)
virtual MNELIB::MNESourceEstimate calculateInverse(const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
static int PowellOffset(int p_iRow, int p_iNumPoints)
static double subcorr(MatrixX6T &p_matProj_G, const MatrixXT &p_pMatU_B)
Definition rapmusic.cpp:616
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorXT
Definition rapmusic.h:113
Eigen::Matrix< double, Eigen::Dynamic, 6 > MatrixX6T
Definition rapmusic.h:107
bool init(MNELIB::MNEForwardSolution &p_pFwd, bool p_bSparsed=false, int p_iN=2, double p_dThr=0.5)
Definition rapmusic.cpp:109
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixXT
Definition rapmusic.h:105
int m_iNumLeadFieldCombinations
Definition rapmusic.h:309
static void insertSource(int p_iDipoleIdx1, int p_iDipoleIdx2, const Vector6T &p_vec_phi_k_1, double p_valCor, QList< DipolePair< double > > &p_RapDipoles)
Definition rapmusic.cpp:852
void calcOrthProj(const MatrixXT &p_matA_k_1, MatrixXT &p_matOrthProj) const
Definition rapmusic.cpp:761
Pair ** m_ppPairIdxCombinations
Definition rapmusic.h:311
virtual MNELIB::MNESourceEstimate calculateInverse(const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
Definition rapmusic.cpp:206
static int useFullRank(const MatrixXT &p_Mat, const MatrixXT &p_matSigma_src, MatrixXT &p_matFull_Rank, int type=NOT_TRANSPOSED)
Definition rapmusic.h:376
MNELIB::MNEForwardSolution m_ForwardSolution
Definition rapmusic.h:300
static void getGainMatrixPair(const MatrixXT &p_matGainMarix, MatrixX6T &p_matGainMarix_Pair, int p_iIdx1, int p_iIdx2)
Definition rapmusic.cpp:841
int calcPhi_s(const MatrixXT &p_matMeasurement, MatrixXT *&p_pMatPhi_s) const
Definition rapmusic.cpp:587
Eigen::Matrix< double, 6, 1 > Vector6T
Definition rapmusic.h:115
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)
Definition rapmusic.cpp:746