MNE-CPP 0.1.9
A Framework for Electrophysiology
Loading...
Searching...
No Matches
mne_inverse_operator.h
Go to the documentation of this file.
1//=============================================================================================================
38#ifndef MNE_INVERSE_OPERATOR_H
39#define MNE_INVERSE_OPERATOR_H
40
41//=============================================================================================================
42// INCLUDES
43//=============================================================================================================
44
45#include "mne_global.h"
46#include "mne_sourcespace.h"
47#include "mne_forwardsolution.h"
48
49#include <fiff/fiff_types.h>
51#include <fiff/fiff_proj.h>
52#include <fiff/fiff_cov.h>
53#include <fiff/fiff_info.h>
54
55#include <utils/mnemath.h>
56
57//=============================================================================================================
58// EIGEN INCLUDES
59//=============================================================================================================
60
61#include <Eigen/Core>
62
63//=============================================================================================================
64// QT INCLUDES
65//=============================================================================================================
66
67#include <QList>
68
69//=============================================================================================================
70// FORWARD DECLARATIONS
71//=============================================================================================================
72
73namespace FSLIB
74{
75 class Label;
76}
77
78//=============================================================================================================
79// DEFINE NAMESPACE MNELIB
80//=============================================================================================================
81
82namespace MNELIB
83{
84
85//=========================================================================================================
90{
91 Eigen::VectorXi roiIdx;
92 Eigen::MatrixXd ctrs;
93 Eigen::VectorXd sumd;
94 Eigen::MatrixXd D;
96 qint32 iLabelIdxOut;
97};
98
99//=========================================================================================================
104{
105 Eigen::MatrixXd matRoiMT;
106 Eigen::MatrixXd matRoiMTOrig;
108 qint32 nClusters;
109 Eigen::VectorXi idcs;
110 qint32 iLabelIdxIn;
112 QString sDistMeasure;
114 RegionMTOut cluster() const
115 {
116 QString t_sDistMeasure;
117 if(sDistMeasure.isEmpty())
118 t_sDistMeasure = QString("cityblock");
119 else
120 t_sDistMeasure = sDistMeasure;
121
122 // Kmeans Reduction
123 RegionMTOut p_RegionMTOut;
124
125 UTILSLIB::KMeans t_kMeans(t_sDistMeasure, QString("sample"), 5);
126
127 t_kMeans.calculate(this->matRoiMT, this->nClusters, p_RegionMTOut.roiIdx, p_RegionMTOut.ctrs, p_RegionMTOut.sumd, p_RegionMTOut.D);
128
129 p_RegionMTOut.iLabelIdxOut = this->iLabelIdxIn;
130
131 return p_RegionMTOut;
132 }
133};
134
135//=============================================================================================================
142{
143public:
144 typedef QSharedPointer<MNEInverseOperator> SPtr;
145 typedef QSharedPointer<const MNEInverseOperator> ConstSPtr;
147 //=========================================================================================================
152
153 //=========================================================================================================
159 MNEInverseOperator(QIODevice& p_IODevice);
160
161 //=========================================================================================================
174 const MNEForwardSolution& forward,
175 const FIFFLIB::FiffCov& p_noise_cov,
176 float loose = 0.2f,
177 float depth = 0.8f,
178 bool fixed = false,
179 bool limit_depth_chs = true);
180
181 //=========================================================================================================
187 MNEInverseOperator(const MNEInverseOperator &p_MNEInverseOperator);
188
189 //=========================================================================================================
194
195 //=========================================================================================================
212 bool assemble_kernel(const FSLIB::Label &label,
213 QString method,
214 bool pick_normal,
215 Eigen::MatrixXd &K,
216 Eigen::SparseMatrix<double> &noise_norm,
217 QList<Eigen::VectorXi> &vertno);
218
219 //=========================================================================================================
227 bool check_ch_names(const FIFFLIB::FiffInfo &info) const;
228
229 //=========================================================================================================
240 Eigen::MatrixXd cluster_kernel(const FSLIB::AnnotationSet &p_AnnotationSet,
241 qint32 p_iClusterSize,
242 Eigen::MatrixXd& p_D,
243 QString p_sMethod = "cityblock") const;
244
245 //=========================================================================================================
251 inline Eigen::MatrixXd& getKernel();
252
253 //=========================================================================================================
259 inline Eigen::MatrixXd getKernel() const;
260
261 //=========================================================================================================
267 inline bool isFixedOrient() const;
268
269 //=========================================================================================================
283 static MNEInverseOperator make_inverse_operator(const FIFFLIB::FiffInfo &info,
284 MNEForwardSolution forward,
285 const FIFFLIB::FiffCov& p_noise_cov,
286 float loose = 0.2f,
287 float depth = 0.8f,
288 bool fixed = false,
289 bool limit_depth_chs = true);
290
291 //=========================================================================================================
306 MNEInverseOperator prepare_inverse_operator(qint32 nave,
307 float lambda2,
308 bool dSPM,
309 bool sLORETA = false) const;
310
311 //=========================================================================================================
324 static bool read_inverse_operator(QIODevice &p_IODevice, MNEInverseOperator& inv);
325
326 //=========================================================================================================
334 void write(QIODevice &p_IODevice);
335
336 //=========================================================================================================
342 void writeToStream(FIFFLIB::FiffStream* p_pStream);
343
344 //=========================================================================================================
353 friend std::ostream& operator<<(std::ostream& out, const MNELIB::MNEInverseOperator &p_MNEInverseOperator);
354
355public:
357 FIFFLIB::fiff_int_t methods;
358 FIFFLIB::fiff_int_t source_ori;
359 FIFFLIB::fiff_int_t nsource;
360 FIFFLIB::fiff_int_t nchan;
361 FIFFLIB::fiff_int_t coord_frame;
362 Eigen::MatrixXf source_nn;
363 Eigen::VectorXd sing;
374 FIFFLIB::fiff_int_t nave;
375 QList<FIFFLIB::FiffProj> projs;
376 Eigen::MatrixXd proj;
377 Eigen::MatrixXd whitener;
378 Eigen::VectorXd reginv;
379 Eigen::SparseMatrix<double> noisenorm;
381private:
382 Eigen::MatrixXd m_K;
383};
384
385//=============================================================================================================
386// INLINE DEFINITIONS
387//=============================================================================================================
388
389inline Eigen::MatrixXd& MNEInverseOperator::getKernel()
390{
391 return m_K;
392}
393
394//=============================================================================================================
395
396inline Eigen::MatrixXd MNEInverseOperator::getKernel() const
397{
398 return m_K;
399}
400
401//=============================================================================================================
402
404{
405 return this->source_ori == FIFFV_MNE_FIXED_ORI;
406}
407
408//=============================================================================================================
409
410inline std::ostream& operator<<(std::ostream& out, const MNELIB::MNEInverseOperator &p_MNEInverseOperator)
411{
412 out << "#### MNE Inverse Operator ####\n";
413
414 out << "\n methods: " << p_MNEInverseOperator.methods << std::endl;
415 out << "\n source_ori: " << p_MNEInverseOperator.source_ori << std::endl;
416 out << "\n nsource: " << p_MNEInverseOperator.nsource << std::endl;
417 out << "\n nchan: " << p_MNEInverseOperator.nchan << std::endl;
418 out << "\n coord_frame:\n\t" << p_MNEInverseOperator.coord_frame << std::endl;
419
420 out << "\n eigen_leads: " << *(p_MNEInverseOperator.eigen_leads) << std::endl;
421 out << "\n eigen_fields:\n\t" << *(p_MNEInverseOperator.eigen_fields) << std::endl;
422
423 return out;
424}
425} // NAMESPACE
426
427#ifndef metatype_mneinverseoperatorsptr
428#define metatype_mneinverseoperatorsptr
429Q_DECLARE_METATYPE(QSharedPointer<MNELIB::MNEInverseOperator>);
430#endif
431
432#ifndef metatype_mneinverseoperators
433#define metatype_mneinverseoperators
435#endif
436
437#endif // MNE_INVERSE_OPERATOR_H
FiffInfo class declaration.
FiffNamedMatrix class declaration.
FiffProj class declaration.
FiffCov class declaration.
Definitions for describing the objects in a FIFF file.
MNEMath class declaration.
mne library export/import macros.
#define MNESHARED_EXPORT
Definition mne_global.h:56
MNEForwardSolution class declaration, which provides the forward solution including the source space ...
MNESourceSpace class declaration.
Coordinate transformation description.
covariance data
Definition fiff_cov.h:78
QSharedDataPointer< FiffCov > SDPtr
Definition fiff_cov.h:82
FIFF measurement file information.
Definition fiff_info.h:85
light measurement info
QSharedDataPointer< FiffNamedMatrix > SDPtr
FIFF File I/O routines.
Annotation set.
Freesurfer/MNE label.
Definition label.h:81
Eigen::MatrixXd matRoiMTOrig
Eigen::MatrixXd matRoiMT
FIFFLIB::FiffCov::SDPtr fmri_prior
QSharedPointer< const MNEInverseOperator > ConstSPtr
QSharedPointer< MNEInverseOperator > SPtr
QList< FIFFLIB::FiffProj > projs
Eigen::SparseMatrix< double > noisenorm
FIFFLIB::FiffCov::SDPtr orient_prior
FIFFLIB::FiffNamedMatrix::SDPtr eigen_leads
FIFFLIB::FiffCoordTrans mri_head_t
FIFFLIB::FiffCov::SDPtr depth_prior
FIFFLIB::FiffCov::SDPtr noise_cov
FIFFLIB::FiffCov::SDPtr source_cov
FIFFLIB::FiffNamedMatrix::SDPtr eigen_fields
Source Space descritpion.
K-Means Clustering.
Definition kmeans.h:73
bool calculate(Eigen::MatrixXd X, qint32 kClusters, Eigen::VectorXi &idx, Eigen::MatrixXd &C, Eigen::VectorXd &sumD, Eigen::MatrixXd &D)
Definition kmeans.cpp:120
Q_DECLARE_METATYPE(QSharedPointer< MNELIB::MNEInverseOperator >)