v2.0.0
Loading...
Searching...
No Matches
mne_hemisphere.h
Go to the documentation of this file.
1//=============================================================================================================
36
37#ifndef MNE_HEMISPHERE_H
38#define MNE_HEMISPHERE_H
39
40//=============================================================================================================
41// INCLUDES
42//=============================================================================================================
43
44#include "mne_global.h"
45#include "mne_source_space.h"
46#include "mne_cluster_info.h"
47
48#include <fiff/fiff_types.h>
49#include <fiff/fiff.h>
50
51//=============================================================================================================
52// EIGEN INCLUDES
53//=============================================================================================================
54
55#include <Eigen/Core>
56#include <Eigen/SparseCore>
57
58//=============================================================================================================
59// QT INCLUDES
60//=============================================================================================================
61
62#include <QList>
63#include <vector>
64
65//=============================================================================================================
66// DEFINE NAMESPACE MNELIB
67//=============================================================================================================
68
69namespace MNELIB
70{
71
72//=============================================================================================================
73// FORWARD DECLARATIONS
74//=============================================================================================================
75
76//=============================================================================================================
83{
84public:
85 using SPtr = std::shared_ptr<MNEHemisphere>;
86 using ConstSPtr = std::shared_ptr<const MNEHemisphere>;
87
88 //=========================================================================================================
93
94 //=========================================================================================================
100 MNEHemisphere(const MNEHemisphere& p_MNEHemisphere);
101
102 //=========================================================================================================
106 ~MNEHemisphere() override;
107
108 //=========================================================================================================
114 MNESourceSpace::SPtr clone() const override;
115
116 //=========================================================================================================
124 bool add_geometry_info();
125
126 //=========================================================================================================
134
135 //=========================================================================================================
142 bool compute_patch_info();
143
144 //=========================================================================================================
148 void clear();
149
150 //=========================================================================================================
158 Eigen::MatrixXf& getTriCoords(float p_fScaling = 1.0f);
159
160 //=========================================================================================================
166 inline bool isClustered() const;
167
168 //=========================================================================================================
183
184 //=========================================================================================================
194 void writeToStream(FIFFLIB::FiffStream* p_pStream);
195
196 //=========================================================================================================
205
206 //ToDo write(IODevice &)
207
215 friend bool operator== (const MNEHemisphere &a, const MNEHemisphere &b);
216
217public:
218 // --- Fields inherited from MNESurfaceOrVolume via MNESourceSpace ---
219 // type, id, np, ntri, coord_frame, rr, nn, nuse, inuse, vertno,
220 // nuse_tri, dist_limit, dist, nearest, neighbor_tri, neighbor_vert
221 // are all available through the base class.
222 // Use nearestVertIdx() / nearestDistVec() accessors for VectorXi/VectorXd views.
223
224 // --- Shadowing fields (different type from base — TODO: unify in Phase 2) ---
225 // tris / use_tris removed: now use inherited itris / use_itris (TrianglesT, row-major)
226 // nearest / nearest_dist removed: now use inherited vector<MNENearest> nearest
227 // + nearestVertIdx() / nearestDistVec() accessors
228 // dist removed: now use inherited FiffSparseMatrix dist + toEigenSparse() / fromEigenSparse()
229 QList<Eigen::VectorXi> pinfo;
230 Eigen::VectorXi patch_inds;
231
232 // --- MNEHemisphere-specific fields ---
233 Eigen::MatrixX3d tri_cent;
234 Eigen::MatrixX3d tri_nn;
235 Eigen::VectorXd tri_area;
236 Eigen::MatrixX3d use_tri_cent;
237 Eigen::MatrixX3d use_tri_nn;
238 Eigen::VectorXd use_tri_area;
239
241private:
242 // Newly added
243 Eigen::MatrixXf m_TriCoords;
244};
245
246//=============================================================================================================
247// INLINE DEFINITIONS
248//=============================================================================================================
249
250inline bool MNEHemisphere::isClustered() const
251{
252 return !cluster_info.isEmpty();
253}
254
255//=============================================================================================================
256
257inline bool operator== (const MNEHemisphere &a, const MNEHemisphere &b)
258{
259 if(a.pinfo.size() == b.pinfo.size()) {
260 for(int i = 0; i < a.pinfo.size(); ++i) {
261 if(!a.pinfo.at(i).isApprox(b.pinfo.at(i))) {
262 return false;
263 }
264 }
265 } else {
266 return false;
267 }
268
269 return (a.type == b.type &&
270 a.id == b.id &&
271 a.np == b.np &&
272 a.ntri == b.ntri &&
273 a.coord_frame == b.coord_frame &&
274 a.rr.isApprox(b.rr, 0.0001f) &&
275 a.nn.isApprox(b.nn, 0.0001f) &&
276 a.itris.isApprox(b.itris) &&
277 a.nuse == b.nuse &&
278 a.inuse.isApprox(b.inuse) &&
279 a.vertno.isApprox(b.vertno) &&
280 a.nuse_tri == b.nuse_tri &&
281 a.use_itris.isApprox(b.use_itris) &&
282 a.nearestVertIdx().isApprox(b.nearestVertIdx()) &&
283 a.nearestDistVec().isApprox(b.nearestDistVec(), 0.0001) &&
284 a.patch_inds.isApprox(b.patch_inds) &&
285 //a.dist_limit == b.dist_limit && //TODO: We still not sure if dist_limit can also be a matrix. This needs to be debugged
286 a.dist.toEigenSparse().isApprox(b.dist.toEigenSparse(), 0.0001) &&
287 a.tri_cent.isApprox(b.tri_cent, 0.0001) &&
288 a.tri_nn.isApprox(b.tri_nn, 0.0001) &&
289 a.tri_area.isApprox(b.tri_area, 0.0001) &&
290 a.use_tri_cent.isApprox(b.use_tri_cent, 0.0001) &&
291 a.use_tri_nn.isApprox(b.use_tri_nn, 0.0001) &&
292 a.use_tri_area.isApprox(b.use_tri_area, 0.0001) &&
293 [&a, &b]() {
294 if (a.neighbor_tri.size() != b.neighbor_tri.size()) return false;
295 for (size_t i = 0; i < a.neighbor_tri.size(); ++i) {
296 if (a.neighbor_tri[i].size() != b.neighbor_tri[i].size()) return false;
297 if (a.neighbor_tri[i].size() > 0 && !(a.neighbor_tri[i].array() == b.neighbor_tri[i].array()).all()) return false;
298 }
299 return true;
300 }() &&
301 [&a, &b]() {
302 if (a.neighbor_vert.size() != b.neighbor_vert.size()) return false;
303 for (size_t i = 0; i < a.neighbor_vert.size(); ++i) {
304 if (a.neighbor_vert[i].size() != b.neighbor_vert[i].size()) return false;
305 if (a.neighbor_vert[i].size() > 0 && !(a.neighbor_vert[i].array() == b.neighbor_vert[i].array()).all()) return false;
306 }
307 return true;
308 }() &&
309 a.cluster_info == b.cluster_info &&
310 a.m_TriCoords.isApprox(b.m_TriCoords, 0.0001f));
311}
312} // NAMESPACE
313
314#endif // MNE_HEMISPHERE_H
FIFF class declaration, which provides static wrapper functions to stay consistent with mne matlab to...
Old fiff_type declarations - replace them.
mne library export/import macros.
#define MNESHARED_EXPORT
Definition mne_global.h:52
MNEClusterInfo class declaration, which provides cluster information.
MNESourceSpace class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
bool operator==(const MNEClusterInfo &a, const MNEClusterInfo &b)
qint32 fiff_int_t
Definition fiff_types.h:89
Coordinate transformation description.
Eigen::SparseMatrix< double > toEigenSparse() const
FIFF File I/O routines.
cluster information
void writeToStream(FIFFLIB::FiffStream *p_pStream)
Eigen::MatrixX3d tri_cent
Eigen::VectorXd use_tri_area
Eigen::MatrixX3d use_tri_cent
Eigen::MatrixXf & getTriCoords(float p_fScaling=1.0f)
Eigen::MatrixX3d tri_nn
MNEHemisphere & operator=(const MNEHemisphere &other)
MNESourceSpace::SPtr clone() const override
bool transform_hemisphere_to(FIFFLIB::fiff_int_t dest, const FIFFLIB::FiffCoordTrans &p_Trans)
MNEClusterInfo cluster_info
Eigen::MatrixX3d use_tri_nn
Eigen::VectorXd tri_area
std::shared_ptr< const MNEHemisphere > ConstSPtr
QList< Eigen::VectorXi > pinfo
std::shared_ptr< MNEHemisphere > SPtr
Eigen::VectorXi patch_inds
std::shared_ptr< MNESourceSpace > SPtr
Eigen::VectorXd nearestDistVec() const
FIFFLIB::FiffSparseMatrix dist
Eigen::VectorXi nearestVertIdx() const