MNE-CPP  beta 1.0
Public Types | Public Member Functions | Static Public Member Functions | List of all members
UTILSLIB::MNEMath Class Reference

Math methods. More...

#include <mnemath.h>

Public Types

typedef std::pair< int, int > IdxIntValue
 

Public Member Functions

virtual ~MNEMath ()
 

Static Public Member Functions

static VectorXd * combine_xyz (const VectorXd &vec)
 
static double getConditionNumber (const MatrixXd &A, VectorXd &s)
 
static double getConditionSlope (const MatrixXd &A, VectorXd &s)
 
static void get_whitener (MatrixXd &A, bool pca, QString ch_type, VectorXd &eig, MatrixXd &eigvec)
 
static VectorXi intersect (const VectorXi &v1, const VectorXi &v2, VectorXi &idx_sel)
 
static bool issparse (VectorXd &v)
 
static MatrixXd legendre (qint32 n, const VectorXd &X, QString normalize=QString("unnorm"))
 
static SparseMatrix< double > * make_block_diag (const MatrixXd &A, qint32 n)
 
static int nchoose2 (int n)
 
static qint32 rank (const MatrixXd &A, double tol=1e-8)
 
static MatrixXd rescale (const MatrixXd &data, const RowVectorXf &times, QPair< QVariant, QVariant > baseline, QString mode)
 
template<typename T >
static VectorXi sort (Matrix< T, Dynamic, 1 > &v, bool desc=true)
 
template<typename T >
static VectorXi sort (Matrix< T, Dynamic, 1 > &v_prime, Matrix< T, Dynamic, Dynamic > &mat, bool desc=true)
 
template<typename T >
static std::vector< Triplet< T > > sortrows (const std::vector< Triplet< T > > &A, qint32 column=0)
 
template<typename T >
static bool compareIdxValuePairBiggerThan (const std::pair< int, T > &lhs, const std::pair< int, T > &rhs)
 
template<typename T >
static bool compareIdxValuePairSmallerThan (const std::pair< int, T > &lhs, const std::pair< int, T > &rhs)
 
template<typename T >
static bool compareTripletFirstEntry (const Triplet< T > &lhs, const Triplet< T > &rhs)
 
template<typename T >
static bool compareTripletSecondEntry (const Triplet< T > &lhs, const Triplet< T > &rhs)
 
template<typename T >
static double log2 (const T d)
 

Detailed Description

Math methods.

ToDo make this a template class Generalized math methods used by mne methods

Definition at line 101 of file mnemath.h.

Member Typedef Documentation

typedef std::pair<int,int> UTILSLIB::MNEMath::IdxIntValue

Typedef of a pair of ints.

Definition at line 104 of file mnemath.h.

Constructor & Destructor Documentation

virtual UTILSLIB::MNEMath::~MNEMath ( )
inlinevirtual

Destroys the MNEMath object

Definition at line 110 of file mnemath.h.

Member Function Documentation

VectorXd * MNEMath::combine_xyz ( const VectorXd &  vec)
static

ToDo make this a template function

mne_combine_xyz

MNE toolbox root function

Compute the three Cartesian components of a vector together

Parameters
[in]vecInput row vector [ x1 y1 z1 ... x_n y_n z_n ]
Returns
Output vector [x1^2+y1^2+z1^2 ... x_n^2+y_n^2+z_n^2 ]

Definition at line 79 of file mnemath.cpp.

template<typename T >
bool UTILSLIB::MNEMath::compareIdxValuePairBiggerThan ( const std::pair< int, T > &  lhs,
const std::pair< int, T > &  rhs 
)
inlinestatic

Compares two index-value-pairs.

Parameters
[in]lhsleft hand side of the comparison
[in]rhsright hand side of the comparison
Returns
true if value of lhs is bigger than value of rhs

Definition at line 439 of file mnemath.h.

template<typename T >
bool UTILSLIB::MNEMath::compareIdxValuePairSmallerThan ( const std::pair< int, T > &  lhs,
const std::pair< int, T > &  rhs 
)
inlinestatic

Compares two index-value-pairs.

Parameters
[in]lhsleft hand side of the comparison
[in]rhsright hand side of the comparison
Returns
true if value of lhs is smaller than value of rhs

Definition at line 448 of file mnemath.h.

template<typename T >
bool UTILSLIB::MNEMath::compareTripletFirstEntry ( const Triplet< T > &  lhs,
const Triplet< T > &  rhs 
)
inlinestatic

Compares triplet first entry

Parameters
[in]lhsleft hand side of the comparison
[in]rhsright hand side of the comparison
Returns
true if value of lhs is smaller than value of rhs

Definition at line 457 of file mnemath.h.

template<typename T >
bool UTILSLIB::MNEMath::compareTripletSecondEntry ( const Triplet< T > &  lhs,
const Triplet< T > &  rhs 
)
inlinestatic

Compares triplet second entry

Parameters
[in]lhsleft hand side of the comparison
[in]rhsright hand side of the comparison
Returns
true if value of lhs is smaller than value of rhs

Definition at line 466 of file mnemath.h.

void MNEMath::get_whitener ( MatrixXd &  A,
bool  pca,
QString  ch_type,
VectorXd &  eig,
MatrixXd &  eigvec 
)
static

Returns the whitener of a given matrix.

Parameters
[in]AMatrix to compute the whitener from
[in]pcaperform a pca
Returns
rank of matrix A

Definition at line 129 of file mnemath.cpp.

double MNEMath::getConditionNumber ( const MatrixXd &  A,
VectorXd &  s 
)
static

MNE toolbox root function ###: Implementation of the mne_block_diag function - decoding part

Returns the condition number of a given matrix.

Parameters
[in]AMatrix to compute the condition number from
Returns
the condition number

Definition at line 103 of file mnemath.cpp.

double MNEMath::getConditionSlope ( const MatrixXd &  A,
VectorXd &  s 
)
static

Returns the condition slope of a given matrix.

Parameters
[in]AMatrix to compute the condition number from
Returns
the condition slope

Definition at line 116 of file mnemath.cpp.

VectorXi MNEMath::intersect ( const VectorXi &  v1,
const VectorXi &  v2,
VectorXi &  idx_sel 
)
static

Find the intersection of two vectors

Parameters
[in]v1Input vector 1
[in]v2Input vector 2
[out]idx_selIndex of intersection based on v1
Returns
the sorted, unique values that are in both of the input arrays.

Definition at line 158 of file mnemath.cpp.

bool MNEMath::issparse ( VectorXd &  v)
static

Determines if a given data (stored as vector v) are representing a sparse matrix. ToDo: status is experimental -> needs to be increased in speed.

Parameters
[in]vdata to be tested
Returns
true if sparse false otherwise;

Definition at line 255 of file mnemath.cpp.

MatrixXd MNEMath::legendre ( qint32  n,
const VectorXd &  X,
QString  normalize = QString("unnorm") 
)
static

LEGENDRE Associated Legendre function.

P = LEGENDRE(N,X) computes the associated Legendre functions of degree N and order M = 0, 1, ..., N, evaluated for each element of X. N must be a scalar integer and X must contain real values between -1 <= X <= 1.

Returns
associated Legendre functions

Definition at line 277 of file mnemath.cpp.

template<typename T >
double UTILSLIB::MNEMath::log2 ( const T  d)
inlinestatic

Compute log2 of given number

Parameters
[in]dinput value
Returns
double result of log2 operation

Definition at line 475 of file mnemath.h.

SparseMatrix< double > * MNEMath::make_block_diag ( const MatrixXd &  A,
qint32  n 
)
static

ToDo make this a template function

MNE toolbox root function ###: Implementation of the mne_block_diag function - encoding part

Make a sparse block diagonal matrix

Returns a sparse block diagonal, diagonalized from the elements in "A". "A" is ma x na, comprising bdn=(na/"n") blocks of submatrices. Each submatrix is ma x "n", and these submatrices are placed down the diagonal of the matrix.

Parameters
[in,out]AMatrix which should be diagonlized
[in,out]nColumns of the submatrices
Returns
A sparse block diagonal, diagonalized from the elements in "A".

Definition at line 295 of file mnemath.cpp.

int MNEMath::nchoose2 ( int  n)
static

Calculates the combination of n over 2 (nchoosek(n,2))

Parameters
[in]nThe number of elements which should be combined with each other (n over 2)
Returns
The number of combinations

Definition at line 335 of file mnemath.cpp.

qint32 MNEMath::rank ( const MatrixXd &  A,
double  tol = 1e-8 
)
static

ToDo make this a template function

Returns the rank of a matrix A.

Parameters
[in]AMatrix to get the rank from
[in]tolrealtive threshold: biggest singualr value multiplied with tol is smallest singular value considered non-zero
Returns
rank of matrix A

Definition at line 348 of file mnemath.cpp.

MatrixXd MNEMath::rescale ( const MatrixXd &  data,
const RowVectorXf &  times,
QPair< QVariant, QVariant >  baseline,
QString  mode 
)
static

ToDo: Maybe new processing class

Rescale aka baseline correct data

Parameters
[in]dataData Matrix (m x n_time)
[in]timesTime instants is seconds.
[in]baselineIf baseline is (a, b) the interval is between "a (s)" and "b (s)". If a is invalid the beginning of the data is used and if b is invalid then b is set to the end of the interval. If baseline is equal to (invalid, invalid) all the time interval is used.
[in]baseline_usageSee description of parameter baseline.
[in]modeDo baseline correction with ratio (power is divided by mean power during baseline) or zscore (power is divided by standard deviatio of power during baseline after substracting the mean, power = [power - mean(power_baseline)] / std(power_baseline)). ("logratio" | "ratio" | "zscore" | "mean" | "percent")
Returns
rescaled data matrix rescaling.

Definition at line 363 of file mnemath.cpp.

template<typename T >
VectorXi UTILSLIB::MNEMath::sort ( Matrix< T, Dynamic, 1 > &  v,
bool  desc = true 
)
static

Sorts a vector (ascending order) in place and returns the track of the original indeces

Parameters
[in,out]vvector to sort; it's sorted in place
[in]descif true its sorted in a descending order, otherwise ascending (optional, default = true)
Returns
Vector of the original indeces in the new order

Definition at line 368 of file mnemath.h.

template<typename T >
VectorXi UTILSLIB::MNEMath::sort ( Matrix< T, Dynamic, 1 > &  v_prime,
Matrix< T, Dynamic, Dynamic > &  mat,
bool  desc = true 
)
static

Sorts a vector (ascending order) and a corresponding matrix in place and returns the track of the original indeces The matrix is sorted along the columns using the vector values for comparison.

Parameters
[in,out]v_primevector to sort (sorted in place)
[in,out]matmatrix to sort (sorted in place)
[in]descif true its sorted in a descending order, otherwise ascending (optional, default = true)
Returns
Vector of the original indeces in the new order

Definition at line 400 of file mnemath.h.

template<typename T >
std::vector< Triplet< T > > UTILSLIB::MNEMath::sortrows ( const std::vector< Triplet< T > > &  A,
qint32  column = 0 
)
static

Sort rows in ascending order

Parameters
[in]Atriplet vector to sort (sorted in place)
[in]columnsorts the triplet vector based on the column specified
Returns
Vector of the original indeces in the new order

Definition at line 420 of file mnemath.h.


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