v2.0.0
Loading...
Searching...
No Matches
inv_rap_music.h
Go to the documentation of this file.
1//=============================================================================================================
35
36#ifndef INV_RAP_MUSIC_H
37#define INV_RAP_MUSIC_H
38
39//=============================================================================================================
40// INCLUDES
41//=============================================================================================================
42
43#include "../inv_global.h"
44
45#include "inv_dipole.h"
46
49#include <time.h>
50
51#include <QVector>
52
53#include <vector>
54
55//=============================================================================================================
56// EIGEN INCLUDES
57//=============================================================================================================
58
59#include <Eigen/Core>
60#include <Eigen/SVD>
61#include <Eigen/LU>
62
63//=============================================================================================================
64// DEFINE NAMESPACE INVLIB
65//=============================================================================================================
66
67namespace INVLIB
68{
69
70//=============================================================================================================
71// SOME DEFINES
72//=============================================================================================================
73
74#define NOT_TRANSPOSED 0
75#define IS_TRANSPOSED 1
76
77//=============================================================================================================
83struct Pair
84{
85 int x1;
86 int x2;
87};
88
89//=============================================================================================================
102{
103public:
104 typedef QSharedPointer<InvRapMusic> SPtr;
105 typedef QSharedPointer<const InvRapMusic> ConstSPtr;
106
107 //*********************************************************************************************************
108 //=========================================================================================================
109 // TYPEDEFS
110 //=========================================================================================================
111
112 typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> MatrixXT;
114 typedef Eigen::Matrix<double, Eigen::Dynamic, 6> MatrixX6T;
116 typedef Eigen::Matrix<double, 6, Eigen::Dynamic> Matrix6XT;
118 typedef Eigen::Matrix<double, 6, 6> Matrix6T;
120 typedef Eigen::Matrix<double, Eigen::Dynamic, 1> VectorXT;
122 typedef Eigen::Matrix<double, 6, 1> Vector6T;
124
125 //=========================================================================================================
129 InvRapMusic();
130
131 //=========================================================================================================
141 InvRapMusic(MNELIB::MNEForwardSolution& p_pFwd, bool p_bSparsed, int p_iN = 2, double p_dThr = 0.5);
142
143 virtual ~InvRapMusic();
144
145 //=========================================================================================================
156 bool init(MNELIB::MNEForwardSolution& p_pFwd, bool p_bSparsed = false, int p_iN = 2, double p_dThr = 0.5);
157
158 virtual InvSourceEstimate calculateInverse(const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal = false);
159
160 virtual InvSourceEstimate calculateInverse(const Eigen::MatrixXd &data, float tmin, float tstep, bool pick_normal = false) const;
161
162 virtual InvSourceEstimate calculateInverse(const Eigen::MatrixXd& p_matMeasurement, QList< InvDipolePair<double> > &p_RapDipoles) const;
163
164 virtual const char* getName() const;
165
166 virtual const MNELIB::MNESourceSpaces& getSourceSpace() const;
167
168 //=========================================================================================================
175 void setStcAttr(int p_iSampStcWin, float p_fStcOverlap);
176
177protected:
178 //=========================================================================================================
187 int calcPhi_s(const MatrixXT& p_matMeasurement, MatrixXT* &p_pMatPhi_s) const;
188
189 //=========================================================================================================
202 static double subcorr(MatrixX6T& p_matProj_G, const MatrixXT& p_pMatU_B);
203
204 //=========================================================================================================
220 static double subcorr(MatrixX6T& p_matProj_G, const MatrixXT& p_matU_B, Vector6T& p_vec_phi_k_1);
221
222 //=========================================================================================================
232 static void calcA_k_1( const MatrixX6T& p_matG_k_1,
233 const Vector6T& p_matPhi_k_1,
234 const int p_iIdxk_1,
235 MatrixXT& p_matA_k_1);
236
237 //=========================================================================================================
244 void calcOrthProj(const MatrixXT& p_matA_k_1, MatrixXT& p_matOrthProj) const;
245
246 //=========================================================================================================
257 void calcPairCombinations( const int p_iNumPoints,
258 const int p_iNumCombinations,
259 std::vector<Pair>& p_pairIdxCombinations) const;
260
261 //=========================================================================================================
277 static void getPointPair(const int p_iPoints, const int p_iCurIdx, int &p_iIdx1, int &p_iIdx2);
278
279 //=========================================================================================================
288 static void getGainMatrixPair( const MatrixXT& p_matGainMarix,
289 MatrixX6T& p_matGainMarix_Pair,
290 int p_iIdx1, int p_iIdx2);
291
292 //=========================================================================================================
302 static void insertSource( int p_iDipoleIdx1, int p_iDipoleIdx2,
303 const Vector6T &p_vec_phi_k_1,
304 double p_valCor,
305 QList< InvDipolePair<double> > &p_RapDipoles);
306
308
309 int m_iN;
313
317
318 std::vector<Pair> m_ppPairIdxCombinations;
319
321
323
324 //Stc stuff
327
328 //=========================================================================================================
336 static inline int getRank(const MatrixXT& p_matSigma);
337
338 //=========================================================================================================
348 static inline int useFullRank( const MatrixXT& p_Mat,
349 const MatrixXT& p_matSigma_src,
350 MatrixXT& p_matFull_Rank,
351 int type = NOT_TRANSPOSED);
352
353 //=========================================================================================================
360 static inline MatrixXT makeSquareMat(const MatrixXT& p_matF);
361};
362
363//=============================================================================================================
364// INLINE DEFINITIONS
365//=============================================================================================================
366
367inline int InvRapMusic::getRank(const MatrixXT& p_matSigma)
368{
369 int t_iRank;
370 //if once a singularvalue is smaller than epsilon = 10^-5 the following values are also smaller
371 // -> because Singular values are ordered
372 for(t_iRank = p_matSigma.rows()-1; t_iRank > 0; t_iRank--)
373 if (p_matSigma(t_iRank, t_iRank) > 0.00001)
374 break;
375
376 t_iRank++;//rank corresponding to epsilon
377
378 return t_iRank;
379}
380
381//=============================================================================================================
382
383inline int InvRapMusic::useFullRank( const MatrixXT& p_Mat,
384 const MatrixXT& p_matSigma_src,
385 MatrixXT& p_matFull_Rank,
386 int type)
387{
388 int rank = getRank(p_matSigma_src);
389
390 if (type == NOT_TRANSPOSED)
391 p_matFull_Rank = p_Mat.block(0,0,p_Mat.rows(),rank);
392 else
393 p_matFull_Rank = p_Mat.block(0,0,rank,p_Mat.cols());
394
395 return rank;
396}
397
398//=============================================================================================================
399
401{
402 //Make rectangular - p_matF*p_matF^T
403 //MatrixXT FFT = p_matF*p_matF.transpose();
404
405 MatrixXT mat = p_matF.transpose();
406
407 return p_matF*mat;
408}
409} //NAMESPACE
410
411#endif // INV_RAP_MUSIC_H
InvSourceEstimate class declaration.
inverse library export/import macros.
#define INVSHARED_EXPORT
Definition inv_global.h:52
InvDipole class declaration for RAP MUSIC dipole results.
#define NOT_TRANSPOSED
MNEForwardSolution class declaration.
Inverse source estimation (MNE, dSPM, sLORETA, dipole fitting).
Pair of correlated dipole indices and orientations found by the RAP MUSIC scanning step.
Definition inv_dipole.h:59
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
QSharedPointer< InvRapMusic > SPtr
virtual InvSourceEstimate calculateInverse(const Eigen::MatrixXd &p_matMeasurement, QList< InvDipolePair< double > > &p_RapDipoles) const
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
virtual InvSourceEstimate calculateInverse(const Eigen::MatrixXd &data, float tmin, float tstep, bool pick_normal=false) const
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)
QSharedPointer< const InvRapMusic > ConstSPtr
Source Space descritpion.