MNE-CPP  0.1.9
A Framework for Electrophysiology
Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
INVERSELIB::RapMusic Class Reference

The RapMusic class provides the RAP MUSIC Algorithm CPU implementation. ToDo: Paper references. More...

#include <rapmusic.h>

Public Types

typedef QSharedPointer< RapMusicSPtr
 
typedef QSharedPointer< const RapMusicConstSPtr
 
typedef Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixXT
 
typedef Eigen::Matrix< double, Eigen::Dynamic, 6 > MatrixX6T
 
typedef Eigen::Matrix< double, 6, Eigen::Dynamic > Matrix6XT
 
typedef Eigen::Matrix< double, 6, 6 > Matrix6T
 
typedef Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorXT
 
typedef Eigen::Matrix< double, 6, 1 > Vector6T
 

Public Member Functions

 RapMusic ()
 
 RapMusic (MNELIB::MNEForwardSolution &p_pFwd, bool p_bSparsed, int p_iN=2, double p_dThr=0.5)
 
bool init (MNELIB::MNEForwardSolution &p_pFwd, bool p_bSparsed=false, int p_iN=2, double p_dThr=0.5)
 
virtual MNELIB::MNESourceEstimate calculateInverse (const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
 
virtual MNELIB::MNESourceEstimate calculateInverse (const Eigen::MatrixXd &data, float tmin, float tstep, bool pick_normal=false) const
 
virtual MNELIB::MNESourceEstimate calculateInverse (const Eigen::MatrixXd &p_matMeasurement, QList< DipolePair< double > > &p_RapDipoles) const
 
virtual const char * getName () const
 
virtual const MNELIB::MNESourceSpacegetSourceSpace () const
 
void setStcAttr (int p_iSampStcWin, float p_fStcOverlap)
 
- Public Member Functions inherited from INVERSELIB::IInverseAlgorithm
virtual ~IInverseAlgorithm ()
 

Protected Member Functions

int calcPhi_s (const MatrixXT &p_matMeasurement, MatrixXT *&p_pMatPhi_s) const
 
void calcOrthProj (const MatrixXT &p_matA_k_1, MatrixXT &p_matOrthProj) const
 
void calcPairCombinations (const int p_iNumPoints, const int p_iNumCombinations, Pair **p_ppPairIdxCombinations) const
 

Static Protected Member Functions

static double subcorr (MatrixX6T &p_matProj_G, const MatrixXT &p_pMatU_B)
 
static double subcorr (MatrixX6T &p_matProj_G, const MatrixXT &p_matU_B, Vector6T &p_vec_phi_k_1)
 
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)
 
static void getGainMatrixPair (const MatrixXT &p_matGainMarix, MatrixX6T &p_matGainMarix_Pair, int p_iIdx1, int p_iIdx2)
 
static void insertSource (int p_iDipoleIdx1, int p_iDipoleIdx2, const Vector6T &p_vec_phi_k_1, double p_valCor, QList< DipolePair< double > > &p_RapDipoles)
 
static int getRank (const MatrixXT &p_matSigma)
 
static int useFullRank (const MatrixXT &p_Mat, const MatrixXT &p_matSigma_src, MatrixXT &p_matFull_Rank, int type=NOT_TRANSPOSED)
 
static MatrixXT makeSquareMat (const MatrixXT &p_matF)
 

Protected Attributes

MNELIB::MNEForwardSolution m_ForwardSolution
 
int m_iN
 
double m_dThreshold
 
int m_iNumGridPoints
 
int m_iNumChannels
 
int m_iNumLeadFieldCombinations
 
Pair ** m_ppPairIdxCombinations
 
int m_iMaxNumThreads
 
bool m_bIsInit
 
int m_iSamplesStcWindow
 
float m_fStcOverlap
 

Detailed Description

The RapMusic class provides the RAP MUSIC Algorithm CPU implementation. ToDo: Paper references.

ToDo Detailed description

Definition at line 92 of file rapmusic.h.

Inheritance diagram for INVERSELIB::RapMusic:
Inheritance graph

Member Typedef Documentation

◆ ConstSPtr

typedef QSharedPointer<const RapMusic> INVERSELIB::RapMusic::ConstSPtr

Const shared pointer type for RapMusic.

Definition at line 96 of file rapmusic.h.

◆ Matrix6T

typedef Eigen::Matrix<double, 6, 6> INVERSELIB::RapMusic::Matrix6T

Defines Eigen::Matrix<T, 6, 6> as Matrix6T type.

Definition at line 109 of file rapmusic.h.

◆ Matrix6XT

typedef Eigen::Matrix<double, 6, Eigen::Dynamic> INVERSELIB::RapMusic::Matrix6XT

Defines Eigen::Matrix<T, 6, Eigen::Dynamic> as Matrix6XT type.

Definition at line 107 of file rapmusic.h.

◆ MatrixX6T

typedef Eigen::Matrix<double, Eigen::Dynamic, 6> INVERSELIB::RapMusic::MatrixX6T

Defines Eigen::Matrix<T, Eigen::Dynamic, 6> as MatrixX6T type.

Definition at line 105 of file rapmusic.h.

◆ MatrixXT

typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> INVERSELIB::RapMusic::MatrixXT

Defines Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic> as MatrixXT type.

Definition at line 103 of file rapmusic.h.

◆ SPtr

typedef QSharedPointer<RapMusic> INVERSELIB::RapMusic::SPtr

Shared pointer type for RapMusic.

Definition at line 95 of file rapmusic.h.

◆ Vector6T

typedef Eigen::Matrix<double, 6, 1> INVERSELIB::RapMusic::Vector6T

Defines Eigen::Matrix<T, 6, 1> as Vector6T type.

Definition at line 113 of file rapmusic.h.

◆ VectorXT

typedef Eigen::Matrix<double, Eigen::Dynamic, 1> INVERSELIB::RapMusic::VectorXT

Defines Eigen::Matrix<T, Eigen::Dynamic, 1> as VectorXT type.

Definition at line 111 of file rapmusic.h.

Constructor & Destructor Documentation

◆ RapMusic() [1/2]

RapMusic::RapMusic ( )

Default constructor creates an empty RapMusic algorithm which still needs to be initialized.

Definition at line 61 of file rapmusic.cpp.

◆ RapMusic() [2/2]

RapMusic::RapMusic ( MNELIB::MNEForwardSolution p_pFwd,
bool  p_bSparsed,
int  p_iN = 2,
double  p_dThr = 0.5 
)

Constructor which initializes the RapMusic algorithm with the given model.

Parameters
[in]p_FwdThe model which contains the gain matrix and its corresponding grid matrix.
[in]p_bSparsedTrue when sparse matrices should be used.
[in]p_iNThe number (default 2) of uncorrelated sources, which should be found. Starting with. the strongest.
[in]p_dThrThe correlation threshold (default 0.5) at which the search for sources stops.

Definition at line 77 of file rapmusic.cpp.

Member Function Documentation

◆ calcA_k_1()

void RapMusic::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 
)
staticprotected

Calculates the accumulated manifold vectors A_{k1}

Parameters
[in]p_matG_k_1The Lead Field combination for the currently best correlated pair.
[in]p_matPhi_k_1Is equal to u_k_1 in the paper and it is the direction of the currently best. correlated pair.
[in]p_iIdxk_1The current position in the manifold vector array A_k_1.
[out]p_matA_k_1The array of the manifold vectors.

Definition at line 740 of file rapmusic.cpp.

◆ calcOrthProj()

void RapMusic::calcOrthProj ( const MatrixXT p_matA_k_1,
MatrixXT p_matOrthProj 
) const
protected

Calculates the orthogonal projector Phi_A_k_1 like in the paper Mosher 1999 (13)

Parameters
[in]p_matA_k_1The array of the manifold vectors.
[out]p_matOrthProjThe orthogonal projector.

Definition at line 755 of file rapmusic.cpp.

◆ calcPairCombinations()

void RapMusic::calcPairCombinations ( const int  p_iNumPoints,
const int  p_iNumCombinations,
Pair **  p_ppPairIdxCombinations 
) const
protected

Pre-Calculates the gain matrix index combinations to search for a two dipole independent topography (IT = source).

Parameters
[in]p_iNumPointsThe number of Lead Field points -> for dimension check.
[in]p_iNumCombinationsThe number of pair index combinations.
[out]p_ppPairIdxCombinationsThe destination which contains pointer to pointer of index. combinations of Lead Field indices -> Number of pointers = Combination (number of grid points over 2 = Num + 1 C 2)

Definition at line 793 of file rapmusic.cpp.

◆ calcPhi_s()

int RapMusic::calcPhi_s ( const MatrixXT p_matMeasurement,
MatrixXT *&  p_pMatPhi_s 
) const
protected

Computes the signal subspace Phi_s out of the measurement F.

Parameters
[in]p_pMatMeasurementThe current measured data to process (for best performance it should have. the dimension channels x samples with samples = number of channels)
[out]p_pMatPhi_sThe calculated signal subspace.
Returns
The rank of the measurement F (named r lt. Mosher 1998, 1999).

Definition at line 581 of file rapmusic.cpp.

◆ calculateInverse() [1/2]

virtual MNELIB::MNESourceEstimate INVERSELIB::RapMusic::calculateInverse ( const Eigen::MatrixXd &  data,
float  tmin,
float  tstep,
bool  pick_normal = false 
) const
virtual

Applies the inverse algorithm to input data and returns a source estimate.

Parameters
[in]p_fiffEvokedEvoked data.
[in]tminMinimal time point.
[in]tminTime between two samples.
[in]pick_normalIf True, rather than pooling the orientations by taking the norm, only the. radial component is kept. This is only applied when working with loose orientations.
Returns
the calculated source estimation.

Implements INVERSELIB::IInverseAlgorithm.

◆ calculateInverse() [2/2]

MNESourceEstimate RapMusic::calculateInverse ( const FIFFLIB::FiffEvoked p_fiffEvoked,
bool  pick_normal = false 
)
virtual

Applies the inverse algorithm to input data and returns a source estimate.

Parameters
[in]p_fiffEvokedEvoked data.
[in]pick_normalIf True, rather than pooling the orientations by taking the norm, only the. radial component is kept. This is only applied when working with loose orientations.
Returns
the calculated source estimation.

Implements INVERSELIB::IInverseAlgorithm.

Reimplemented in INVERSELIB::PwlRapMusic.

Definition at line 200 of file rapmusic.cpp.

◆ getGainMatrixPair()

void RapMusic::getGainMatrixPair ( const MatrixXT p_matGainMarix,
MatrixX6T p_matGainMarix_Pair,
int  p_iIdx1,
int  p_iIdx2 
)
staticprotected

Returns a gain matrix pair for the given indices

Parameters
[in]p_matGainMarixThe Lead Field matrix.
[out]p_matGainMarix_PairLead Field combination (dimension: m x 6).
[in]p_iIdx1first Lead Field index point.
[in]p_iIdx2second Lead Field index point.

Definition at line 835 of file rapmusic.cpp.

◆ getName()

const char * RapMusic::getName ( ) const
virtual

Returns the algorithm name

Returns
the algorithm name.

Implements INVERSELIB::IInverseAlgorithm.

Reimplemented in INVERSELIB::PwlRapMusic.

Definition at line 186 of file rapmusic.cpp.

◆ getPointPair()

void RapMusic::getPointPair ( const int  p_iPoints,
const int  p_iCurIdx,
int &  p_iIdx1,
int &  p_iIdx2 
)
staticprotected

Calculates the combination indices Idx1 and Idx2 of n points.
[ (0,0) (0,1) (0,2) ... (0,n-1) (0,n)
(1,1) (1,2) ... (1,n-1) (1,n)
(2,2) ... (2,n-1) (2,n)

(n-1,n-1) (n-1,n)
(n,n)]

Parameters
[in]p_iPointsThe number of points n which are combined with each other.
[in]p_iCurIdxThe current combination index (between 0 and nchoosek(n+1,2)).
[out]p_iIdx1The resulting index 1.
[out]p_iIdx2The resulting index 2.

Definition at line 823 of file rapmusic.cpp.

◆ getRank()

int INVERSELIB::RapMusic::getRank ( const MatrixXT p_matSigma)
inlinestaticprotected

Returns the rank r of a singular value matrix based on non-zero singular values (singular value > epsilon = 10^-5)

Parameters
[in]p_matSigmadiagonal matrix which contains the Singular values (Dimension n x n).
Returns
The rank r.

Definition at line 358 of file rapmusic.h.

◆ getSourceSpace()

const MNESourceSpace & RapMusic::getSourceSpace ( ) const
virtual

Returns the current mne source space on which the inverse algorithm is performing on. Either from inverse operator (minimum norm estimate), or from forward solution (beamformers)

Returns
the mne source space information.

Implements INVERSELIB::IInverseAlgorithm.

Definition at line 193 of file rapmusic.cpp.

◆ init()

bool RapMusic::init ( MNELIB::MNEForwardSolution p_pFwd,
bool  p_bSparsed = false,
int  p_iN = 2,
double  p_dThr = 0.5 
)

Initializes the RAP MUSIC algorithm with the given model.

Parameters
[in]p_FwdThe model which contains the gain matrix and its corresponding Grid matrix.
[in]p_bSparsedTrue when sparse matrices should be used.
[in]p_iNThe number (default 2) of uncorrelated sources, which should be found. Starting with. the strongest.
[in]p_dThrThe correlation threshold (default 0.5) at which the search for sources stops.
Returns
true if successful initialized, false otherwise.

Definition at line 103 of file rapmusic.cpp.

◆ insertSource()

void RapMusic::insertSource ( int  p_iDipoleIdx1,
int  p_iDipoleIdx2,
const Vector6T p_vec_phi_k_1,
double  p_valCor,
QList< DipolePair< double > > &  p_RapDipoles 
)
staticprotected

Adds a new correlated dipole pair to th RapDipoles. This function is called by the RAP MUSIC Algorithm.

Parameters
[in]p_iDipoleIdx1Index (Lead Field grid index) of the first dipole.
[in]p_iDipoleIdx2Index (Lead Field grid index) of the second dipole.
[in]p_vec_phi_k_1Array of the dipole directories (phi_x1, phi_y1, phi_z1, phi_x2, phi_y2, phi_z2).
[in]p_valCorCorrelation value of the dipole pair.
[out]p_RapDipolesthe list of dipole pairs.

Definition at line 846 of file rapmusic.cpp.

◆ makeSquareMat()

RapMusic::MatrixXT INVERSELIB::RapMusic::makeSquareMat ( const MatrixXT p_matF)
inlinestaticprotected

Performs F * F^Transposed, is used when n > m

Parameters
[in]p_matFThe matrix which should be transformed.
Returns
F * F^Transposed (we call it FFT ;)).

Definition at line 391 of file rapmusic.h.

◆ setStcAttr()

void RapMusic::setStcAttr ( int  p_iSampStcWin,
float  p_fStcOverlap 
)

Sets the source estimate attributes.

Parameters
[in]p_iSampStcWinSamples per source localization window (default - 1 = not set).
[in]p_fStcOverlapPercentage of localization window overlap.

Definition at line 879 of file rapmusic.cpp.

◆ subcorr() [1/2]

double RapMusic::subcorr ( MatrixX6T p_matProj_G,
const MatrixXT p_matU_B,
Vector6T p_vec_phi_k_1 
)
staticprotected

Computes the subspace correlation between the projected G_rho and the projected signal subspace Phi_s, as well as the resulting direction. For speed-up: we calculate the decomposition of the projected Phi_s before this function. So the argument for this function is U_B instead of Phi_s.

Parameters
[in]p_matProj_GThe projected Lead Field combination. This is a m x 6 matrix composed of the. Lead Field combination of two Points for all m channels and 3 orthogonal components (x y z).
[in]p_matU_BThe matrix U is the subspace projection of the orthogonal projected Phi_s.
[out]p_vec_phi_k_1Returns the orientation for a correlated dipole pair. (phi_x1, phi_y1, phi_z1, phi_x2, phi_y2, phi_z2)
Returns
The maximal correlation c_1 of the subspace correlation of the current projected Lead Field. combination and the projected measurement.

Definition at line 665 of file rapmusic.cpp.

◆ subcorr() [2/2]

double RapMusic::subcorr ( MatrixX6T p_matProj_G,
const MatrixXT p_pMatU_B 
)
staticprotected

Computes the subspace correlation between the projected G_rho and the projected signal subspace Phi_s. For speed-up: we calculate the decomposition of the projected Phi_s before this function. So the argument for this function is U_B instead of Phi_s.

Parameters
[in]p_matProj_GThe projected Lead Field combination. This is a m x 6 matrix composed of the. Lead Field combination of two Points for all m channels and 3 orthogonal components (x y z).
[in]p_matU_BThe matrix U is the subspace projection of the orthogonal projected Phi_s.
Returns
The maximal correlation c_1 of the subspace correlation of the current projected Lead Field. combination and the projected measurement.

Definition at line 610 of file rapmusic.cpp.

◆ useFullRank()

int INVERSELIB::RapMusic::useFullRank ( const MatrixXT p_Mat,
const MatrixXT p_matSigma_src,
MatrixXT p_matFull_Rank,
int  type = NOT_TRANSPOSED 
)
inlinestaticprotected

lt. Mosher 1998 -> Only Retain those Components of U_A and U_B that correspond to nonzero singular values for U_A and U_B the number of columns corresponds to their ranks

Parameters
[in]p_MatThe Matrix which should be reduced to its rank.
[in]p_matSigma_srcThe singular values of the matrix.
[out]p_matFull_RankThe corresponding full rank matrix.
[in]typeWhether p_Mat is transposed, than rows and columns are changed.

Definition at line 374 of file rapmusic.h.

Member Data Documentation

◆ m_bIsInit

bool INVERSELIB::RapMusic::m_bIsInit
protected

Whether the algorithm is initialized.

Definition at line 313 of file rapmusic.h.

◆ m_dThreshold

double INVERSELIB::RapMusic::m_dThreshold
protected

Threshold which defines the minimal correlation. Is the correlation of the found dipole pair smaller as this threshold than the RAP MUSIC calculation is stopped.

Definition at line 301 of file rapmusic.h.

◆ m_ForwardSolution

MNELIB::MNEForwardSolution INVERSELIB::RapMusic::m_ForwardSolution
protected

The Forward operator which should be scanned through

Definition at line 298 of file rapmusic.h.

◆ m_fStcOverlap

float INVERSELIB::RapMusic::m_fStcOverlap
protected

Percentage of localization window overlap.

Definition at line 317 of file rapmusic.h.

◆ m_iMaxNumThreads

int INVERSELIB::RapMusic::m_iMaxNumThreads
protected

Number of available CPU threads.

Definition at line 311 of file rapmusic.h.

◆ m_iN

int INVERSELIB::RapMusic::m_iN
protected

Number of Sources to find

Definition at line 300 of file rapmusic.h.

◆ m_iNumChannels

int INVERSELIB::RapMusic::m_iNumChannels
protected

Number of channels.

Definition at line 306 of file rapmusic.h.

◆ m_iNumGridPoints

int INVERSELIB::RapMusic::m_iNumGridPoints
protected

Number of Grid points.

Definition at line 305 of file rapmusic.h.

◆ m_iNumLeadFieldCombinations

int INVERSELIB::RapMusic::m_iNumLeadFieldCombinations
protected

Number of Lead Filed combinations (grid points + 1 over 2)

Definition at line 307 of file rapmusic.h.

◆ m_iSamplesStcWindow

int INVERSELIB::RapMusic::m_iSamplesStcWindow
protected

Number of samples per localization window.

Definition at line 316 of file rapmusic.h.

◆ m_ppPairIdxCombinations

Pair** INVERSELIB::RapMusic::m_ppPairIdxCombinations
protected

Index combination vector with grid pair indices.

Definition at line 309 of file rapmusic.h.


The documentation for this class was generated from the following files: