MNE-CPP  0.1.9
A Framework for Electrophysiology
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Friends | List of all members
MNELIB::MNEInverseOperator Class Reference

Inverse operator. More...

#include <mne_inverse_operator.h>

Public Types

typedef QSharedPointer< MNEInverseOperatorSPtr
 
typedef QSharedPointer< const MNEInverseOperatorConstSPtr
 

Public Member Functions

 MNEInverseOperator ()
 
 MNEInverseOperator (QIODevice &p_IODevice)
 
 MNEInverseOperator (const FIFFLIB::FiffInfo &info, const MNEForwardSolution &forward, const FIFFLIB::FiffCov &p_noise_cov, float loose=0.2f, float depth=0.8f, bool fixed=false, bool limit_depth_chs=true)
 
 MNEInverseOperator (const MNEInverseOperator &p_MNEInverseOperator)
 
 ~MNEInverseOperator ()
 
bool assemble_kernel (const FSLIB::Label &label, QString method, bool pick_normal, Eigen::MatrixXd &K, Eigen::SparseMatrix< double > &noise_norm, QList< Eigen::VectorXi > &vertno)
 
bool check_ch_names (const FIFFLIB::FiffInfo &info) const
 
Eigen::MatrixXd cluster_kernel (const FSLIB::AnnotationSet &p_AnnotationSet, qint32 p_iClusterSize, Eigen::MatrixXd &p_D, QString p_sMethod="cityblock") const
 
Eigen::MatrixXd & getKernel ()
 
Eigen::MatrixXd getKernel () const
 
bool isFixedOrient () const
 
MNEInverseOperator prepare_inverse_operator (qint32 nave, float lambda2, bool dSPM, bool sLORETA=false) const
 
void write (QIODevice &p_IODevice)
 
void writeToStream (FIFFLIB::FiffStream *p_pStream)
 

Static Public Member Functions

static MNEInverseOperator make_inverse_operator (const FIFFLIB::FiffInfo &info, MNEForwardSolution forward, const FIFFLIB::FiffCov &p_noise_cov, float loose=0.2f, float depth=0.8f, bool fixed=false, bool limit_depth_chs=true)
 
static bool read_inverse_operator (QIODevice &p_IODevice, MNEInverseOperator &inv)
 

Public Attributes

FIFFLIB::FiffInfoBase info
 
FIFFLIB::fiff_int_t methods
 
FIFFLIB::fiff_int_t source_ori
 
FIFFLIB::fiff_int_t nsource
 
FIFFLIB::fiff_int_t nchan
 
FIFFLIB::fiff_int_t coord_frame
 
Eigen::MatrixXf source_nn
 
Eigen::VectorXd sing
 
bool eigen_leads_weighted
 
FIFFLIB::FiffNamedMatrix::SDPtr eigen_leads
 
FIFFLIB::FiffNamedMatrix::SDPtr eigen_fields
 
FIFFLIB::FiffCov::SDPtr noise_cov
 
FIFFLIB::FiffCov::SDPtr source_cov
 
FIFFLIB::FiffCov::SDPtr orient_prior
 
FIFFLIB::FiffCov::SDPtr depth_prior
 
FIFFLIB::FiffCov::SDPtr fmri_prior
 
MNESourceSpace src
 
FIFFLIB::FiffCoordTrans mri_head_t
 
FIFFLIB::fiff_int_t nave
 
QList< FIFFLIB::FiffProjprojs
 
Eigen::MatrixXd proj
 
Eigen::MatrixXd whitener
 
Eigen::VectorXd reginv
 
Eigen::SparseMatrix< double > noisenorm
 

Friends

std::ostream & operator<< (std::ostream &out, const MNELIB::MNEInverseOperator &p_MNEInverseOperator)
 

Detailed Description

Inverse operator.

Inverse operator

Definition at line 141 of file mne_inverse_operator.h.

Member Typedef Documentation

◆ ConstSPtr

Const shared pointer type for MNEInverseOperator.

Definition at line 145 of file mne_inverse_operator.h.

◆ SPtr

Shared pointer type for MNEInverseOperator.

Definition at line 144 of file mne_inverse_operator.h.

Constructor & Destructor Documentation

◆ MNEInverseOperator() [1/4]

MNEInverseOperator::MNEInverseOperator ( )

Default constructor

Definition at line 74 of file mne_inverse_operator.cpp.

◆ MNEInverseOperator() [2/4]

MNEInverseOperator::MNEInverseOperator ( QIODevice &  p_IODevice)

Constructs an inverse operator, by reading from a IO device.

Parameters
[in]p_IODeviceIO device to read from the evoked data set.

Definition at line 96 of file mne_inverse_operator.cpp.

◆ MNEInverseOperator() [3/4]

MNEInverseOperator::MNEInverseOperator ( const FIFFLIB::FiffInfo info,
const MNEForwardSolution forward,
const FIFFLIB::FiffCov p_noise_cov,
float  loose = 0.2f,
float  depth = 0.8f,
bool  fixed = false,
bool  limit_depth_chs = true 
)

Constructs an inverse operator by making one.

Parameters
[in]infoThe measurement info to specify the channels to include. Bad channels in info['bads'] are not used.
[in]forwardForward operator.
[in]p_noise_covThe noise covariance matrix.
[in]loosefloat in [0, 1]. Value that weights the source variances of the dipole components defining the tangent space of the cortical surfaces.
[in]depthfloat in [0, 1]. Depth weighting coefficients. If None, no depth weighting is performed.
[in]fixedUse fixed source orientations normal to the cortical mantle. If True, the loose parameter is ignored.
[in]limit_depth_chsIf True, use only grad channels in depth weighting (equivalent to MNE C code). If grad chanels aren't present, only mag channels will be used (if no mag, then eeg). If False, use all channels.

Definition at line 105 of file mne_inverse_operator.cpp.

◆ MNEInverseOperator() [4/4]

MNEInverseOperator::MNEInverseOperator ( const MNEInverseOperator p_MNEInverseOperator)

Copy constructor.

Parameters
[in]p_MNEInverseOperatorMNE forward solution.

Definition at line 120 of file mne_inverse_operator.cpp.

◆ ~MNEInverseOperator()

MNEInverseOperator::~MNEInverseOperator ( )

Destroys the MNEInverseOperator.

Definition at line 152 of file mne_inverse_operator.cpp.

Member Function Documentation

◆ assemble_kernel()

bool MNEInverseOperator::assemble_kernel ( const FSLIB::Label label,
QString  method,
bool  pick_normal,
Eigen::MatrixXd &  K,
Eigen::SparseMatrix< double > &  noise_norm,
QList< Eigen::VectorXi > &  vertno 
)

Simple matrix multiplication followed by combination of the current components

This does all the data transformations to compute the weights for the eigenleads

Parameters
[in]labellabels.
[in]methodThe applied normals. ("MNE" | "dSPM" | "sLORETA").
[in]pick_normalPick normals.
[out]KKernel.
[out]noise_normNoise normals.
[out]vertnoVertices of the hemispheres.
Returns
the assembled kernel.

Definition at line 158 of file mne_inverse_operator.cpp.

◆ check_ch_names()

bool MNEInverseOperator::check_ch_names ( const FIFFLIB::FiffInfo info) const

Check that channels in inverse operator are measurements.

Parameters
[in]infoThe measurement info.
Returns
true when successful, false otherwise.

Definition at line 307 of file mne_inverse_operator.cpp.

◆ cluster_kernel()

MatrixXd MNEInverseOperator::cluster_kernel ( const FSLIB::AnnotationSet p_AnnotationSet,
qint32  p_iClusterSize,
Eigen::MatrixXd &  p_D,
QString  p_sMethod = "cityblock" 
) const

Clusters the current kernel

Parameters
[in]p_AnnotationSetAnnotation set containing the annotation of left & right hemisphere.
[in]p_iClusterSizeMaximal cluster size per roi.
[out]p_DThe cluster operator.
[in]p_sMethod"cityblock" or "sqeuclidean".
Returns
the clustered kernel.

Definition at line 352 of file mne_inverse_operator.cpp.

◆ getKernel() [1/2]

Eigen::MatrixXd MNELIB::MNEInverseOperator::getKernel ( )
inline

Returns the current kernel

Returns
the current kernel.

Definition at line 389 of file mne_inverse_operator.h.

◆ getKernel() [2/2]

Eigen::MatrixXd MNELIB::MNEInverseOperator::getKernel ( ) const
inline

Returns the current kernel

Returns
the current kernel.

◆ isFixedOrient()

bool MNELIB::MNEInverseOperator::isFixedOrient ( ) const
inline

Has inverse operator fixed orientation?

Returns
true if inverse operator has fixed orientation, false otherwise.

Definition at line 403 of file mne_inverse_operator.h.

◆ make_inverse_operator()

MNEInverseOperator MNEInverseOperator::make_inverse_operator ( const FIFFLIB::FiffInfo info,
MNEForwardSolution  forward,
const FIFFLIB::FiffCov p_noise_cov,
float  loose = 0.2f,
float  depth = 0.8f,
bool  fixed = false,
bool  limit_depth_chs = true 
)
static

Assembles the inverse operator.

Parameters
[in]infoThe measurement info to specify the channels to include. Bad channels in info['bads'] are not used.
[in]forwardForward operator.
[in]p_noise_covThe noise covariance matrix.
[in]loosefloat in [0, 1]. Value that weights the source variances of the dipole components defining the tangent space of the cortical surfaces.
[in]depthfloat in [0, 1]. Depth weighting coefficients. If None, no depth weighting is performed.
[in]fixedUse fixed source orientations normal to the cortical mantle. If True, the loose parameter is ignored.
[in]limit_depth_chsIf True, use only grad channels in depth weighting (equivalent to MNE C code). If grad chanels aren't present, only mag channels will be used (if no mag, then eeg). If False, use all channels.
Returns
the assembled inverse operator.

Definition at line 682 of file mne_inverse_operator.cpp.

◆ prepare_inverse_operator()

MNEInverseOperator MNEInverseOperator::prepare_inverse_operator ( qint32  nave,
float  lambda2,
bool  dSPM,
bool  sLORETA = false 
) const

mne_prepare_inverse_operator

MNE toolbox root function

Prepare for actually computing the inverse

Parameters
[in]naveNumber of averages (scales the noise covariance).
[in]lambda2The 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 926 of file mne_inverse_operator.cpp.

◆ read_inverse_operator()

bool MNEInverseOperator::read_inverse_operator ( QIODevice &  p_IODevice,
MNEInverseOperator inv 
)
static

mne_read_inverse_operator

MNE toolbox root function

Reads the inverse operator decomposition from a fif file

Parameters
[in]p_IODeviceA fiff IO device like a fiff QFile or QTCPSocket.
[in,out]invThe read inverse operator.
Returns
true if succeeded, false otherwise.

Definition at line 1100 of file mne_inverse_operator.cpp.

◆ write()

void MNEInverseOperator::write ( QIODevice &  p_IODevice)

write_inverse_operator

Writes an inverse operator to a fif file

Parameters
[in]p_IODeviceIO device to write the inverse operator to.

Definition at line 1345 of file mne_inverse_operator.cpp.

◆ writeToStream()

void MNEInverseOperator::writeToStream ( FIFFLIB::FiffStream p_pStream)

Writes the inverse operator to a FIFF stream

Parameters
[in]p_pStreamThe stream to write to.

Definition at line 1360 of file mne_inverse_operator.cpp.

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const MNELIB::MNEInverseOperator p_MNEInverseOperator 
)
friend

overloading the stream out operator<<

Parameters
[in]outThe stream to which the fiff covariance should be assigned to.
[in]p_MNEInverseOperatorMNEInverseOperator which should be assigned to the stream.
Returns
the stream with the attached fiff covariance matrix.

Definition at line 410 of file mne_inverse_operator.h.

Member Data Documentation

◆ coord_frame

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::coord_frame

Coordinate system definition.

Definition at line 361 of file mne_inverse_operator.h.

◆ depth_prior

FIFFLIB::FiffCov::SDPtr MNELIB::MNEInverseOperator::depth_prior

Depth priors.

Definition at line 370 of file mne_inverse_operator.h.

◆ eigen_fields

FIFFLIB::FiffNamedMatrix::SDPtr MNELIB::MNEInverseOperator::eigen_fields

Eigen fields.

Definition at line 366 of file mne_inverse_operator.h.

◆ eigen_leads

FIFFLIB::FiffNamedMatrix::SDPtr MNELIB::MNEInverseOperator::eigen_leads

Eigen leads.

Definition at line 365 of file mne_inverse_operator.h.

◆ eigen_leads_weighted

bool MNELIB::MNEInverseOperator::eigen_leads_weighted

If eigen lead are weighted.

Definition at line 364 of file mne_inverse_operator.h.

◆ fmri_prior

FIFFLIB::FiffCov::SDPtr MNELIB::MNEInverseOperator::fmri_prior

fMRI priors.

Definition at line 371 of file mne_inverse_operator.h.

◆ info

FIFFLIB::FiffInfoBase MNELIB::MNEInverseOperator::info

light weighted measurement info.

Definition at line 356 of file mne_inverse_operator.h.

◆ methods

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::methods

MEG, EEG or both.

Definition at line 357 of file mne_inverse_operator.h.

◆ mri_head_t

FIFFLIB::FiffCoordTrans MNELIB::MNEInverseOperator::mri_head_t

MRI head coordinate transformation.

Definition at line 373 of file mne_inverse_operator.h.

◆ nave

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::nave

Number of averages used to regularize the solution. Set to 1 on single Epoch by default.

Definition at line 374 of file mne_inverse_operator.h.

◆ nchan

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::nchan

Number of channels.

Definition at line 360 of file mne_inverse_operator.h.

◆ noise_cov

FIFFLIB::FiffCov::SDPtr MNELIB::MNEInverseOperator::noise_cov

Noise covariance matrix.

Definition at line 367 of file mne_inverse_operator.h.

◆ noisenorm

Eigen::SparseMatrix<double> MNELIB::MNEInverseOperator::noisenorm

These are the noise-normalization factors.

Definition at line 379 of file mne_inverse_operator.h.

◆ nsource

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::nsource

Number of source points.

Definition at line 359 of file mne_inverse_operator.h.

◆ orient_prior

FIFFLIB::FiffCov::SDPtr MNELIB::MNEInverseOperator::orient_prior

Orientation priors.

Definition at line 369 of file mne_inverse_operator.h.

◆ proj

Eigen::MatrixXd MNELIB::MNEInverseOperator::proj

The projector to apply to the data.

Definition at line 376 of file mne_inverse_operator.h.

◆ projs

QList<FIFFLIB::FiffProj> MNELIB::MNEInverseOperator::projs

SSP operator.

Definition at line 375 of file mne_inverse_operator.h.

◆ reginv

Eigen::VectorXd MNELIB::MNEInverseOperator::reginv

The diagonal matrix implementing. regularization and the inverse.

Definition at line 378 of file mne_inverse_operator.h.

◆ sing

Eigen::VectorXd MNELIB::MNEInverseOperator::sing

Singular values.

Definition at line 363 of file mne_inverse_operator.h.

◆ source_cov

FIFFLIB::FiffCov::SDPtr MNELIB::MNEInverseOperator::source_cov

Source covariance matrix.

Definition at line 368 of file mne_inverse_operator.h.

◆ source_nn

Eigen::MatrixXf MNELIB::MNEInverseOperator::source_nn

Source normals.

Definition at line 362 of file mne_inverse_operator.h.

◆ source_ori

FIFFLIB::fiff_int_t MNELIB::MNEInverseOperator::source_ori

Source orientation: f.

Definition at line 358 of file mne_inverse_operator.h.

◆ src

MNESourceSpace MNELIB::MNEInverseOperator::src

Source Space.

Definition at line 372 of file mne_inverse_operator.h.

◆ whitener

Eigen::MatrixXd MNELIB::MNEInverseOperator::whitener

Whitens the data.

Definition at line 377 of file mne_inverse_operator.h.


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