v2.0.0
Loading...
Searching...
No Matches
INVLIB::InvMinimumNorm Class Reference

Minimum norm estimation. More...

#include <inv_minimum_norm.h>

Public Types

typedef QSharedPointer< InvMinimumNormSPtr
typedef QSharedPointer< const InvMinimumNormConstSPtr

Public Member Functions

 InvMinimumNorm (const MNELIB::MNEInverseOperator &p_inverseOperator, float lambda, const QString method)
 InvMinimumNorm (const MNELIB::MNEInverseOperator &p_inverseOperator, float lambda, bool dSPM, bool sLORETA)
virtual ~InvMinimumNorm ()
virtual InvSourceEstimate calculateInverse (const FIFFLIB::FiffEvoked &p_fiffEvoked, bool pick_normal=false)
virtual InvSourceEstimate calculateInverse (const Eigen::MatrixXd &data, float tmin, float tstep, bool pick_normal=false) const
virtual void doInverseSetup (qint32 nave, bool pick_normal=false)
virtual const char * getName () const
virtual const MNELIB::MNESourceSpacesgetSourceSpace () const
MNELIB::MNEInverseOperatorgetPreparedInverseOperator ()
void setMethod (QString method)
void setMethod (bool dSPM, bool sLORETA, bool eLoreta=false)
void setRegularization (float lambda)
void setELoretaOptions (int maxIter=20, double eps=1e-6, bool forceEqual=false)
Eigen::MatrixXd & getKernel ()

Detailed Description

Minimum norm estimation.

Minimum norm estimation algorithm.

Computes L2 minimum-norm, dSPM, or sLORETA source estimates from MEG/EEG data using a pre-computed inverse operator.

References:

  • Hamalainen & Ilmoniemi, Med. & Biol. Eng. & Comput. 32, 35-42, 1994.
  • Dale et al., Neuron 26, 55-67, 2000 (dSPM).
  • Pascual-Marqui, Methods Find. Exp. Clin. Pharmacol. 24D, 5-12, 2002 (sLORETA).

Definition at line 77 of file inv_minimum_norm.h.

Member Typedef Documentation

◆ ConstSPtr

typedef QSharedPointer<const InvMinimumNorm> INVLIB::InvMinimumNorm::ConstSPtr

Const shared pointer type for InvMinimumNorm.

Definition at line 81 of file inv_minimum_norm.h.

◆ SPtr

Shared pointer type for InvMinimumNorm.

Definition at line 80 of file inv_minimum_norm.h.

Constructor & Destructor Documentation

◆ InvMinimumNorm() [1/2]

InvMinimumNorm::InvMinimumNorm ( const MNELIB::MNEInverseOperator & p_inverseOperator,
float lambda,
const QString method )
explicit

Constructs minimum norm inverse algorithm

Parameters
[in]p_inverseOperatorThe inverse operator.
[in]lambdaThe regularization factor.
[in]methodUse mininum norm, dSPM or sLORETA. ("MNE" | "dSPM" | "sLORETA").
Returns
the prepared inverse operator.

Definition at line 74 of file inv_minimum_norm.cpp.

◆ InvMinimumNorm() [2/2]

InvMinimumNorm::InvMinimumNorm ( const MNELIB::MNEInverseOperator & p_inverseOperator,
float lambda,
bool dSPM,
bool sLORETA )
explicit

Constructs minimum norm inverse algorithm

Parameters
[in]p_inverseOperatorThe inverse operator.
[in]lambdaThe regularization factor.
[in]dSPMCompute the noise-normalization factors for dSPM?.
[in]sLORETACompute the noise-normalization factors for sLORETA?.
Returns
the prepared inverse operator.

Definition at line 88 of file inv_minimum_norm.cpp.

◆ ~InvMinimumNorm()

virtual INVLIB::InvMinimumNorm::~InvMinimumNorm ( )
inlinevirtual

Definition at line 108 of file inv_minimum_norm.h.

Member Function Documentation

◆ calculateInverse() [1/2]

virtual InvSourceEstimate INVLIB::InvMinimumNorm::calculateInverse ( const Eigen::MatrixXd & data,
float tmin,
float tstep,
bool pick_normal = false ) const
virtual

◆ calculateInverse() [2/2]

InvSourceEstimate InvMinimumNorm::calculateInverse ( const FIFFLIB::FiffEvoked & p_fiffEvoked,
bool pick_normal = false )
virtual

Computes a L2-norm inverse solution Actual code using these principles might be different because the inverse operator is often reused across data sets.

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.

Definition at line 102 of file inv_minimum_norm.cpp.

◆ doInverseSetup()

void InvMinimumNorm::doInverseSetup ( qint32 nave,
bool pick_normal = false )
virtual

Perform the inverse setup: Prepares this inverse operator and assembles the kernel.

Parameters
[in]naveNumber of averages to use.
[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.

Definition at line 256 of file inv_minimum_norm.cpp.

◆ getKernel()

Eigen::MatrixXd & INVLIB::InvMinimumNorm::getKernel ( )
inline

Get the assembled kernel

Returns
the assembled kernel.

Definition at line 240 of file inv_minimum_norm.h.

◆ getName()

const char * InvMinimumNorm::getName ( ) const
virtual

Get the name of the inverse operator.

Returns
the name of the inverse operator.

Definition at line 278 of file inv_minimum_norm.cpp.

◆ getPreparedInverseOperator()

MNELIB::MNEInverseOperator & INVLIB::InvMinimumNorm::getPreparedInverseOperator ( )
inline

Get the prepared inverse operator.

Returns
the prepared inverse operator.

Definition at line 247 of file inv_minimum_norm.h.

◆ getSourceSpace()

const MNESourceSpaces & InvMinimumNorm::getSourceSpace ( ) const
virtual

Get the source space corresponding to this inverse operator.

Returns
the source space corresponding to this inverse operator.

Definition at line 285 of file inv_minimum_norm.cpp.

◆ setELoretaOptions()

void InvMinimumNorm::setELoretaOptions ( int maxIter = 20,
double eps = 1e-6,
bool forceEqual = false )

Set eLORETA options.

Parameters
[in]maxIterMaximum number of eLORETA weight iterations (default: 20).
[in]epsConvergence threshold for weight fitting (default: 1e-6).
[in]forceEqualIf true, use uniform orientation weights (default: false).

Definition at line 348 of file inv_minimum_norm.cpp.

◆ setMethod() [1/2]

void InvMinimumNorm::setMethod ( bool dSPM,
bool sLORETA,
bool eLoreta = false )

Set minimum norm algorithm method ("MNE" | "dSPM" | "sLORETA" | "eLORETA")

Parameters
[in]dSPMCompute the noise-normalization factors for dSPM?.
[in]sLORETACompute the noise-normalization factors for sLORETA?.
[in]eLoretaCompute the eLORETA source weights?.

Definition at line 314 of file inv_minimum_norm.cpp.

◆ setMethod() [2/2]

void InvMinimumNorm::setMethod ( QString method)

Set minimum norm algorithm method ("MNE" | "dSPM" | "sLORETA" | "eLORETA")

Parameters
[in]methodUse minimum norm, dSPM, sLORETA, or eLORETA.

Definition at line 292 of file inv_minimum_norm.cpp.

◆ setRegularization()

void InvMinimumNorm::setRegularization ( float lambda)

Set regularization factor

Parameters
[in]lambdaThe regularization factor.

Definition at line 341 of file inv_minimum_norm.cpp.


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