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::MNEForwardSolution Class Reference

Forward operator. More...

#include <mne_forwardsolution.h>

Public Types

typedef QSharedPointer< MNEForwardSolutionSPtr
 
typedef QSharedPointer< const MNEForwardSolutionConstSPtr
 

Public Member Functions

 MNEForwardSolution ()
 
 MNEForwardSolution (QIODevice &p_IODevice, bool force_fixed=false, bool surf_ori=false, const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList, bool bExcludeBads=false)
 
 MNEForwardSolution (const MNEForwardSolution &p_MNEForwardSolution)
 
 ~MNEForwardSolution ()
 
void clear ()
 
MNEForwardSolution cluster_forward_solution (const FSLIB::AnnotationSet &p_AnnotationSet, qint32 p_iClusterSize, Eigen::MatrixXd &p_D=defaultD, const FIFFLIB::FiffCov &p_pNoise_cov=defaultCov, const FIFFLIB::FiffInfo &p_pInfo=defaultInfo, QString p_sMethod="cityblock") const
 
FIFFLIB::FiffCov compute_orient_prior (float loose=0.2)
 
bool isClustered () const
 
bool isEmpty () const
 
bool isFixedOrient () const
 
MNEForwardSolution pick_channels (const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList) const
 
MNEForwardSolution pick_regions (const QList< FSLIB::Label > &p_qListLabels) const
 
MNEForwardSolution pick_types (bool meg, bool eeg, const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList) const
 
void prepare_forward (const FIFFLIB::FiffInfo &p_info, const FIFFLIB::FiffCov &p_noise_cov, bool p_pca, FIFFLIB::FiffInfo &p_outFwdInfo, Eigen::MatrixXd &gain, FIFFLIB::FiffCov &p_outNoiseCov, Eigen::MatrixXd &p_outWhitener, qint32 &p_outNumNonZero) const
 
Eigen::VectorXi tripletSelection (const Eigen::VectorXi &p_vecIdxSelection) const
 
MNEForwardSolution reduce_forward_solution (qint32 p_iNumDipoles, Eigen::MatrixXd &p_D) const
 
void to_fixed_ori ()
 
Eigen::MatrixX3f getSourcePositionsByLabel (const QList< FSLIB::Label > &lPickedLabels, const FSLIB::SurfaceSet &tSurfSetInflated)
 

Static Public Member Functions

static FIFFLIB::FiffCov compute_depth_prior (const Eigen::MatrixXd &Gain, const FIFFLIB::FiffInfo &gain_info, bool is_fixed_ori, double exp=0.8, double limit=10.0, const Eigen::MatrixXd &patch_areas=FIFFLIB::defaultConstMatrixXd, bool limit_depth_chs=false)
 
static bool read (QIODevice &p_IODevice, MNEForwardSolution &fwd, bool force_fixed=false, bool surf_ori=false, const QStringList &include=FIFFLIB::defaultQStringList, const QStringList &exclude=FIFFLIB::defaultQStringList, bool bExcludeBads=true)
 
static void restrict_gain_matrix (Eigen::MatrixXd &G, const FIFFLIB::FiffInfo &info)
 

Public Attributes

FIFFLIB::FiffInfoBase info
 
FIFFLIB::fiff_int_t source_ori
 
bool surf_ori
 
FIFFLIB::fiff_int_t coord_frame
 
FIFFLIB::fiff_int_t nsource
 
FIFFLIB::fiff_int_t nchan
 
FIFFLIB::FiffNamedMatrix::SDPtr sol
 
FIFFLIB::FiffNamedMatrix::SDPtr sol_grad
 
FIFFLIB::FiffCoordTrans mri_head_t
 
MNESourceSpace src
 
Eigen::MatrixX3f source_rr
 
Eigen::MatrixX3f source_nn
 

Friends

std::ostream & operator<< (std::ostream &out, const MNELIB::MNEForwardSolution &p_MNEForwardSolution)
 
bool operator== (const MNEForwardSolution &a, const MNEForwardSolution &b)
 

Detailed Description

Forward operator.

Forward operator

Definition at line 170 of file mne_forwardsolution.h.

Member Typedef Documentation

◆ ConstSPtr

Const shared pointer type for MNEForwardSolution.

Definition at line 174 of file mne_forwardsolution.h.

◆ SPtr

Shared pointer type for MNEForwardSolution.

Definition at line 173 of file mne_forwardsolution.h.

Constructor & Destructor Documentation

◆ MNEForwardSolution() [1/3]

MNEForwardSolution::MNEForwardSolution ( )

Default constructor.

Definition at line 82 of file mne_forwardsolution.cpp.

◆ MNEForwardSolution() [2/3]

MNEForwardSolution::MNEForwardSolution ( QIODevice &  p_IODevice,
bool  force_fixed = false,
bool  surf_ori = false,
const QStringList &  include = FIFFLIB::defaultQStringList,
const QStringList &  exclude = FIFFLIB::defaultQStringList,
bool  bExcludeBads = false 
)

Constructs a forward operator, by reading from a IO device.

Parameters
[in]p_IODeviceIO device to read from the forward operator.
[in]force_fixedForce fixed source orientation mode? (optional).
[in]surf_oriUse surface based source coordinate system? (optional).
[in]includeInclude these channels (optional).
[in]excludeExclude these channels (optional).
[in]bExcludeBadsIf true bads are also read; default = false (optional).

Definition at line 99 of file mne_forwardsolution.cpp.

◆ MNEForwardSolution() [3/3]

MNEForwardSolution::MNEForwardSolution ( const MNEForwardSolution p_MNEForwardSolution)

Copy constructor.

Parameters
[in]p_MNEForwardSolutionMNE forward solution.

Definition at line 121 of file mne_forwardsolution.cpp.

◆ ~MNEForwardSolution()

MNEForwardSolution::~MNEForwardSolution ( )

Destroys the MNEForwardSolution.

Definition at line 139 of file mne_forwardsolution.cpp.

Member Function Documentation

◆ clear()

void MNEForwardSolution::clear ( )

Initializes the MNE forward solution.

Definition at line 145 of file mne_forwardsolution.cpp.

◆ cluster_forward_solution()

MNEForwardSolution MNEForwardSolution::cluster_forward_solution ( const FSLIB::AnnotationSet p_AnnotationSet,
qint32  p_iClusterSize,
Eigen::MatrixXd &  p_D = defaultD,
const FIFFLIB::FiffCov p_pNoise_cov = defaultCov,
const FIFFLIB::FiffInfo p_pInfo = defaultInfo,
QString  p_sMethod = "cityblock" 
) const

Cluster the forward solution and stores the result to p_fwdOut. The clustering is done by using the provided annotations

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_pNoise_cov.
[in]p_pInfo.
[in]p_sMethod"cityblock" or "sqeuclidean".
Returns
clustered MNE forward solution.

Definition at line 163 of file mne_forwardsolution.cpp.

◆ compute_depth_prior()

FiffCov MNEForwardSolution::compute_depth_prior ( const Eigen::MatrixXd &  Gain,
const FIFFLIB::FiffInfo gain_info,
bool  is_fixed_ori,
double  exp = 0.8,
double  limit = 10.0,
const Eigen::MatrixXd &  patch_areas = FIFFLIB::defaultConstMatrixXd,
bool  limit_depth_chs = false 
)
static

Compute weighting for depth prior. ToDo move this to FiffCov

Parameters
[in]Gaingain matrix.
[in]gain_infoThe measurement info to specify the channels to include.
[in]is_fixed_oriFixed orientation?.
[in]expfloat in [0, 1]. Depth weighting coefficients. If None, no depth weighting is performed. (optional; default = 0.8).
[in]limit(optional; default = 10.0).
[in]patch_areas(optional).
[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. (optional).
Returns
the depth prior.

Definition at line 807 of file mne_forwardsolution.cpp.

◆ compute_orient_prior()

FiffCov MNEForwardSolution::compute_orient_prior ( float  loose = 0.2)

Compute orientation prior

Parameters
[in]looseThe loose orientation parameter.
Returns
Orientation priors.

Definition at line 916 of file mne_forwardsolution.cpp.

◆ getSourcePositionsByLabel()

MatrixX3f MNEForwardSolution::getSourcePositionsByLabel ( const QList< FSLIB::Label > &  lPickedLabels,
const FSLIB::SurfaceSet tSurfSetInflated 
)

Returns the positions of the specified sources based on their beloning labels

Parameters
[in]lPickedLabelsThe stream to which the MNE forward solution should be assigned to.
[in]tSurfSetInflatedThe surface used to pick the source from, based on their index specified bzy this forward solution.
Returns
the source position in 3D space.

Definition at line 1835 of file mne_forwardsolution.cpp.

◆ isClustered()

bool MNELIB::MNEForwardSolution::isClustered ( ) const
inline

Indicates whether fwd conatins a clustered forward solution.

Returns
true if forward solution is clustered, false otherwise.

Definition at line 538 of file mne_forwardsolution.h.

◆ isEmpty()

bool MNELIB::MNEForwardSolution::isEmpty ( ) const
inline

True if FIFF measurement file information is empty.

Returns
true if FIFF measurement file information is empty.

Definition at line 545 of file mne_forwardsolution.h.

◆ isFixedOrient()

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

Has forward operator fixed orientation?

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

Definition at line 552 of file mne_forwardsolution.h.

◆ pick_channels()

MNEForwardSolution MNEForwardSolution::pick_channels ( const QStringList &  include = FIFFLIB::defaultQStringList,
const QStringList &  exclude = FIFFLIB::defaultQStringList 
) const

mne.fiff.pick_channels_forward

Pick channels from forward operator

Parameters
[in]includeList of channels to include. (if None, include all available).
[in]excludeChannels to exclude (if None, do not exclude any).
Returns
Forward solution restricted to selected channel types.

Definition at line 965 of file mne_forwardsolution.cpp.

◆ pick_regions()

MNEForwardSolution MNEForwardSolution::pick_regions ( const QList< FSLIB::Label > &  p_qListLabels) const

Reduces a forward solution to selected regions

Parameters
[in]p_qListLabelsROIs.
Returns
the reduced forward solution.

Definition at line 1029 of file mne_forwardsolution.cpp.

◆ pick_types()

MNEForwardSolution MNEForwardSolution::pick_types ( bool  meg,
bool  eeg,
const QStringList &  include = FIFFLIB::defaultQStringList,
const QStringList &  exclude = FIFFLIB::defaultQStringList 
) const

mne.fiff.pick_types_forward

Pick by channel type and names from a forward operator

Parameters
[in]megInclude MEG channels.
[in]eegInclude EEG channels.
[in]includeAdditional channels to include (if empty, do not add any).
[in]excludeChannels to exclude (if empty, do not exclude any).
Returns
Forward solution restricted to selected channel types.

Definition at line 1080 of file mne_forwardsolution.cpp.

◆ prepare_forward()

void MNEForwardSolution::prepare_forward ( const FIFFLIB::FiffInfo p_info,
const FIFFLIB::FiffCov p_noise_cov,
bool  p_pca,
FIFFLIB::FiffInfo p_outFwdInfo,
Eigen::MatrixXd &  gain,
FIFFLIB::FiffCov p_outNoiseCov,
Eigen::MatrixXd &  p_outWhitener,
qint32 &  p_outNumNonZero 
) const

Prepare forward for assembling the inverse operator

Parameters
[in]p_infoThe measurement info to specify the channels to include. Bad channels in info['bads'] are not used.
[in]p_noise_covThe noise covariance matrix.
[in]p_pcaCalculate pca or not.
[out]ch_namesSelected channel names.
[out]gainGain matrix.
[out]p_outNoiseCovnoise covariance matrix.
[out]p_outWhitenerWhitener.
[out]p_outNumNonZerothe rank (non zeros).

Definition at line 1093 of file mne_forwardsolution.cpp.

◆ read()

bool MNEForwardSolution::read ( QIODevice &  p_IODevice,
MNEForwardSolution fwd,
bool  force_fixed = false,
bool  surf_ori = false,
const QStringList &  include = FIFFLIB::defaultQStringList,
const QStringList &  exclude = FIFFLIB::defaultQStringList,
bool  bExcludeBads = true 
)
static

MNE toolbox root function ###: Definition of the mne_read_forward_solution function

Reads a forward solution from a fif file

Parameters
[in]p_IODeviceA fiff IO device like a fiff QFile or QTCPSocket.
[out]fwdA forward solution from a fif file.
[in]force_fixedForce fixed source orientation mode? (optional).
[in]surf_oriUse surface based source coordinate system? (optional).
[in]includeInclude these channels (optional).
[in]excludeExclude these channels (optional).
[in]bExcludeBadsIf true bads are also read; default = false (optional).
Returns
true if succeeded, false otherwise.

Definition at line 1192 of file mne_forwardsolution.cpp.

◆ reduce_forward_solution()

MNEForwardSolution MNEForwardSolution::reduce_forward_solution ( qint32  p_iNumDipoles,
Eigen::MatrixXd &  p_D 
) const

reduces the forward solution and stores the result to p_fwdOut.

Parameters
[in]p_iNumDipolesDesired number of dipoles.
[out]p_DThe reduction operator.
Returns
reduced MNE forward solution.

Definition at line 702 of file mne_forwardsolution.cpp.

◆ restrict_gain_matrix()

void MNEForwardSolution::restrict_gain_matrix ( Eigen::MatrixXd &  G,
const FIFFLIB::FiffInfo info 
)
static

Restrict gain matrix entries for optimal depth weighting

Parameters
[in,out]GGain matrix to be restricted; result is stored in place.
[in]infoFiff information.

Definition at line 1772 of file mne_forwardsolution.cpp.

◆ to_fixed_ori()

void MNEForwardSolution::to_fixed_ori ( )

Helper to convert the forward solution to fixed ori from free

Definition at line 1817 of file mne_forwardsolution.cpp.

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const MNELIB::MNEForwardSolution p_MNEForwardSolution 
)
friend

overloading the stream out operator<<

Parameters
[in]outThe stream to which the MNE forward solution should be assigned to.
[in]p_MNEForwardSolutionMNE forward solution which should be assigned to the stream.
Returns
the stream with the attached fiff projector.

Definition at line 559 of file mne_forwardsolution.h.

◆ operator==

bool operator== ( const MNEForwardSolution a,
const MNEForwardSolution b 
)
friend

Overloaded == operator to compare an object to this instance.

Parameters
[in]objectThe object which should be compared to.
Returns
true if equal, false otherwise.

Definition at line 575 of file mne_forwardsolution.h.

Member Data Documentation

◆ coord_frame

FIFFLIB::fiff_int_t MNELIB::MNEForwardSolution::coord_frame

Coil coordinate system definition.

Definition at line 523 of file mne_forwardsolution.h.

◆ info

FIFFLIB::FiffInfoBase MNELIB::MNEForwardSolution::info

light weighted measurement info.

Definition at line 520 of file mne_forwardsolution.h.

◆ mri_head_t

FIFFLIB::FiffCoordTrans MNELIB::MNEForwardSolution::mri_head_t

MRI head coordinate transformation.

Definition at line 528 of file mne_forwardsolution.h.

◆ nchan

FIFFLIB::fiff_int_t MNELIB::MNEForwardSolution::nchan

Number of channels.

Definition at line 525 of file mne_forwardsolution.h.

◆ nsource

FIFFLIB::fiff_int_t MNELIB::MNEForwardSolution::nsource

Number of source dipoles.

Definition at line 524 of file mne_forwardsolution.h.

◆ sol

FIFFLIB::FiffNamedMatrix::SDPtr MNELIB::MNEForwardSolution::sol

Forward solution.

Definition at line 526 of file mne_forwardsolution.h.

◆ sol_grad

FIFFLIB::FiffNamedMatrix::SDPtr MNELIB::MNEForwardSolution::sol_grad

ToDo...

Definition at line 527 of file mne_forwardsolution.h.

◆ source_nn

Eigen::MatrixX3f MNELIB::MNEForwardSolution::source_nn

Source normals (number depends on fixed or free orientation).

Definition at line 531 of file mne_forwardsolution.h.

◆ source_ori

FIFFLIB::fiff_int_t MNELIB::MNEForwardSolution::source_ori

Source orientation: fixed or free.

Definition at line 521 of file mne_forwardsolution.h.

◆ source_rr

Eigen::MatrixX3f MNELIB::MNEForwardSolution::source_rr

Source locations.

Definition at line 530 of file mne_forwardsolution.h.

◆ src

MNESourceSpace MNELIB::MNEForwardSolution::src

Geometric description of the source spaces (hemispheres).

Definition at line 529 of file mne_forwardsolution.h.

◆ surf_ori

bool MNELIB::MNEForwardSolution::surf_ori

If surface oriented.

Definition at line 522 of file mne_forwardsolution.h.


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