v2.0.0
Loading...
Searching...
No Matches
mne_hemisphere.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "mne_hemisphere.h"
42#include "mne_nearest.h"
43
44#include <utils/mnemath.h>
45
46#include <algorithm>
47
48//=============================================================================================================
49// USED NAMESPACES
50//=============================================================================================================
51
52using namespace MNELIB;
53using namespace UTILSLIB;
54using namespace Eigen;
55using namespace FIFFLIB;
56
57//=============================================================================================================
58// DEFINE MEMBER METHODS
59//=============================================================================================================
60
63, patch_inds(VectorXi::Zero(0))
64, tri_cent(MatrixX3d::Zero(0,3))
65, tri_nn(MatrixX3d::Zero(0,3))
66, tri_area(VectorXd::Zero(0))
67, use_tri_cent(MatrixX3d::Zero(0,3))
68, use_tri_nn(MatrixX3d::Zero(0,3))
69, use_tri_area(VectorXd::Zero(0))
70//, m_TriCoords()
71//, m_pGeometryData(NULL)
72{
73 // Override some base class defaults to match MNEHemisphere semantics
74 this->type = 1;
75 this->id = -1;
76 this->np = -1;
77 this->ntri = -1;
78 this->coord_frame = -1;
79 this->nuse = -1;
80 this->nuse_tri = -1;
81 this->dist_limit = -1;
82}
83
84//=============================================================================================================
85
88, pinfo(p_MNEHemisphere.pinfo)
89, patch_inds(p_MNEHemisphere.patch_inds)
90, tri_cent(p_MNEHemisphere.tri_cent)
91, tri_nn(p_MNEHemisphere.tri_nn)
92, tri_area(p_MNEHemisphere.tri_area)
93, use_tri_cent(p_MNEHemisphere.use_tri_cent)
94, use_tri_nn(p_MNEHemisphere.use_tri_nn)
95, use_tri_area(p_MNEHemisphere.use_tri_area)
96, cluster_info(p_MNEHemisphere.cluster_info)
97, m_TriCoords(p_MNEHemisphere.m_TriCoords)
98{
99 // Copy base class (MNESurfaceOrVolume) fields that MNEHemisphere uses
100 this->type = p_MNEHemisphere.type;
101 this->id = p_MNEHemisphere.id;
102 this->np = p_MNEHemisphere.np;
103 this->ntri = p_MNEHemisphere.ntri;
104 this->coord_frame = p_MNEHemisphere.coord_frame;
105 this->rr = p_MNEHemisphere.rr;
106 this->nn = p_MNEHemisphere.nn;
107 this->nuse = p_MNEHemisphere.nuse;
108 this->inuse = p_MNEHemisphere.inuse;
109 this->vertno = p_MNEHemisphere.vertno;
110 this->itris = p_MNEHemisphere.itris;
111 this->use_itris = p_MNEHemisphere.use_itris;
112 this->nuse_tri = p_MNEHemisphere.nuse_tri;
113 this->dist_limit = p_MNEHemisphere.dist_limit;
114 this->dist = p_MNEHemisphere.dist;
115 this->nearest = p_MNEHemisphere.nearest;
116 this->neighbor_tri = p_MNEHemisphere.neighbor_tri;
117 this->neighbor_vert = p_MNEHemisphere.neighbor_vert;
118}
119
120//=============================================================================================================
121
123{
124 if (this != &other) {
125 // Copy base class (MNESurfaceOrVolume) fields
126 this->type = other.type;
127 this->id = other.id;
128 this->np = other.np;
129 this->ntri = other.ntri;
130 this->coord_frame = other.coord_frame;
131 this->rr = other.rr;
132 this->nn = other.nn;
133 this->nuse = other.nuse;
134 this->inuse = other.inuse;
135 this->vertno = other.vertno;
136 this->nuse_tri = other.nuse_tri;
137 this->dist_limit = other.dist_limit;
138 this->neighbor_tri = other.neighbor_tri;
139 this->neighbor_vert = other.neighbor_vert;
140
141 // Copy triangle index fields (inherited)
142 this->itris = other.itris;
143 this->use_itris = other.use_itris;
144
145 // Copy inherited fields with value semantics
146 this->dist = other.dist;
147 this->nearest = other.nearest;
148
149 // Copy MNEHemisphere fields
150 this->pinfo = other.pinfo;
151 this->patch_inds = other.patch_inds;
152 this->tri_cent = other.tri_cent;
153 this->tri_nn = other.tri_nn;
154 this->tri_area = other.tri_area;
155 this->use_tri_cent = other.use_tri_cent;
156 this->use_tri_nn = other.use_tri_nn;
157 this->use_tri_area = other.use_tri_area;
158 this->cluster_info = other.cluster_info;
159 this->m_TriCoords = other.m_TriCoords;
160 }
161 return *this;
162}
163
164//=============================================================================================================
165
169
170//=============================================================================================================
171
173{
174 return std::make_shared<MNEHemisphere>(*this);
175}
176
177//=============================================================================================================
178
180{
181 //
182 // Main triangulation
183 //
184 printf("\tCompleting triangulation info...");
185 tri_cent = MatrixX3d::Zero(ntri,3);
186 tri_nn = MatrixX3d::Zero(ntri,3);
187 tri_area = VectorXd::Zero(ntri);
188
189 Matrix3d r;
190 Vector3d a, b;
191 int k = 0;
192 float size = 0;
193 for (qint32 i = 0; i < ntri; ++i)
194 {
195 for ( qint32 j = 0; j < 3; ++j)
196 {
197 k = itris(i, j);
198
199 r(j,0) = rr(k, 0);
200 r(j,1) = rr(k, 1);
201 r(j,2) = rr(k, 2);
202
203 tri_cent(i, 0) += rr(k, 0);
204 tri_cent(i, 1) += rr(k, 1);
205 tri_cent(i, 2) += rr(k, 2);
206 }
207 tri_cent.row(i) /= 3.0f;
208
209 //cross product {cross((r2-r1),(r3-r1))}
210 a = r.row(1) - r.row(0 );
211 b = r.row(2) - r.row(0);
212 tri_nn(i,0) = a(1)*b(2)-a(2)*b(1);
213 tri_nn(i,1) = a(2)*b(0)-a(0)*b(2);
214 tri_nn(i,2) = a(0)*b(1)-a(1)*b(0);
215
216 //area
217 size = tri_nn.row(i)*tri_nn.row(i).transpose();
218 size = std::pow(size, 0.5f );
219
220 tri_area(i) = size/2.0f;
221 tri_nn.row(i) /= size;
222
223 }
224 printf("[done]\n");
225
226 //
227 // Selected triangles
228 //
229 printf("\tCompleting selection triangulation info...");
230 if (nuse_tri > 0)
231 {
232 use_tri_cent = MatrixX3d::Zero(nuse_tri,3);
233 use_tri_nn = MatrixX3d::Zero(nuse_tri,3);
234 use_tri_area = VectorXd::Zero(nuse_tri);
235
236 for (qint32 i = 0; i < nuse_tri; ++i)
237 {
238 for ( qint32 j = 0; j < 3; ++j)
239 {
240 k = use_itris(i, j);
241
242 r(j,0) = rr(k, 0);
243 r(j,1) = rr(k, 1);
244 r(j,2) = rr(k, 2);
245
246 use_tri_cent(i, 0) += rr(k, 0);
247 use_tri_cent(i, 1) += rr(k, 1);
248 use_tri_cent(i, 2) += rr(k, 2);
249 }
250 use_tri_cent.row(i) /= 3.0f;
251
252 //cross product {cross((r2-r1),(r3-r1))}
253 a = r.row(1) - r.row(0 );
254 b = r.row(2) - r.row(0);
255 use_tri_nn(i,0) = a(1)*b(2)-a(2)*b(1);
256 use_tri_nn(i,1) = a(2)*b(0)-a(0)*b(2);
257 use_tri_nn(i,2) = a(0)*b(1)-a(1)*b(0);
258
259 //area
260 size = use_tri_nn.row(i)*use_tri_nn.row(i).transpose();
261 size = std::pow(size, 0.5f );
262
263 use_tri_area(i) = size/2.0f;
264 }
265
266 }
267 printf("[done]\n");
268
269 printf("\tCompleting triangle and vertex neighboring info...");
271 printf("[done]\n");
272
273 return true;
274}
275
276//=============================================================================================================
277
279{
280 if (nearest.empty())
281 {
282 pinfo.clear();
283 patch_inds = VectorXi();
284 return false;
285 }
286
287 printf("\tComputing patch statistics...");
288
289 std::vector< std::pair<int,int> > t_vIndn;
290
291 for(size_t i = 0; i < nearest.size(); ++i)
292 {
293 std::pair<int,int> t_pair(static_cast<int>(i), nearest[i].nearest);
294 t_vIndn.push_back(t_pair);
295 }
296 std::sort(t_vIndn.begin(),t_vIndn.end(), MNEMath::compareIdxValuePairSmallerThan<int> );
297
298 VectorXi nearest_sorted(t_vIndn.size());
299
300 qint32 current = 0;
301 std::vector<qint32> t_vfirsti;
302 t_vfirsti.push_back(current);
303 std::vector<qint32> t_vlasti;
304
305 for(quint32 i = 0; i < t_vIndn.size(); ++i)
306 {
307 nearest_sorted[i] = t_vIndn[i].second;
308 if (t_vIndn[current].second != t_vIndn[i].second)
309 {
310 current = i;
311 t_vlasti.push_back(i-1);
312 t_vfirsti.push_back(current);
313 }
314 }
315 t_vlasti.push_back(static_cast<int>(t_vIndn.size()-1));
316
317 for(quint32 k = 0; k < t_vfirsti.size(); ++k)
318 {
319 std::vector<int> t_vIndex;
320
321 for(int l = t_vfirsti[k]; l <= t_vlasti[k]; ++l)
322 t_vIndex.push_back(t_vIndn[l].first);
323
324 std::sort(t_vIndex.begin(),t_vIndex.end());
325
326 int* t_pV = &t_vIndex[0];
327 Eigen::Map<Eigen::VectorXi> t_vPInfo(t_pV, t_vIndex.size());
328
329 pinfo.append(t_vPInfo);
330 }
331
332 // compute patch indices of the in-use source space vertices
333 std::vector<qint32> patch_verts;
334 patch_verts.reserve(t_vlasti.size());
335 for(quint32 i = 0; i < t_vlasti.size(); ++i)
336 patch_verts.push_back(nearest_sorted[t_vlasti[i]]);
337
338 patch_inds.resize(vertno.size());
339 std::vector<qint32>::iterator it;
340 for(qint32 i = 0; i < vertno.size(); ++i)
341 {
342 it = std::find(patch_verts.begin(), patch_verts.end(), vertno[i]);
343 patch_inds[i] = it-patch_verts.begin();
344 }
345
346 return true;
347}
348
349//=============================================================================================================
350
352{
353 int k,c,p,q;
354 bool found;
355
356 //Create neighboring triangle vector using temporary std::vector for efficient appending
357 {
358 std::vector<std::vector<int>> temp_ntri(this->itris.rows());
359 for (p = 0; p < this->itris.rows(); p++) {
360 for (k = 0; k < 3; k++) {
361 temp_ntri[this->itris(p,k)].push_back(p);
362 }
363 }
364 neighbor_tri.resize(this->itris.rows());
365 for (k = 0; k < static_cast<int>(temp_ntri.size()); k++) {
366 neighbor_tri[k] = Eigen::Map<Eigen::VectorXi>(temp_ntri[k].data(), temp_ntri[k].size());
367 }
368 }
369
370 //Create the neighboring vertices vector using temporary std::vector
371 {
372 std::vector<std::vector<int>> temp_nvert(this->np);
373 for (k = 0; k < this->np; k++) {
374 for (p = 0; p < neighbor_tri[k].size(); p++) {
375 //Fit in the other vertices of the neighboring triangle
376 for (c = 0; c < 3; c++) {
377 int vert = this->itris(neighbor_tri[k][p], c);
378
379 if (vert != k) {
380 found = false;
381
382 for (q = 0; q < static_cast<int>(temp_nvert[k].size()); q++) {
383 if (temp_nvert[k][q] == vert) {
384 found = true;
385 break;
386 }
387 }
388
389 if(!found) {
390 temp_nvert[k].push_back(vert);
391 }
392 }
393 }
394 }
395 }
396 neighbor_vert.resize(this->np);
397 for (k = 0; k < this->np; k++) {
398 neighbor_vert[k] = Eigen::Map<Eigen::VectorXi>(temp_nvert[k].data(), temp_nvert[k].size());
399 }
400 }
401
402 return true;
403}
404
405//=============================================================================================================
406
408{
409 // Reset base class fields
410 type = 1;
411 id = -1;
412 np = -1;
413 ntri = -1;
414 coord_frame = -1;
415 rr = PointsT::Zero(0,3);
416 nn = NormalsT::Zero(0,3);
417 nuse = -1;
418 inuse = VectorXi::Zero(0);
419 vertno = VectorXi::Zero(0);
420 nuse_tri = -1;
421 dist_limit = -1;
422 neighbor_tri.clear();
423 neighbor_vert.clear();
424
425 // Reset inherited value-semantic fields
426 itris = TrianglesT::Zero(0,3);
427 use_itris = TrianglesT::Zero(0,3);
429 nearest.clear();
430
431 // Reset MNEHemisphere fields
432 pinfo.clear();
433 patch_inds = VectorXi::Zero(0);
434 tri_cent = MatrixX3d::Zero(0,3);
435 tri_nn = MatrixX3d::Zero(0,3);
436 tri_area = VectorXd::Zero(0);
437 use_tri_cent = MatrixX3d::Zero(0,3);
438 use_tri_nn = MatrixX3d::Zero(0,3);
439 use_tri_area = VectorXd::Zero(0);
440
441 cluster_info.clear();
442
443 m_TriCoords = MatrixXf();
444}
445
446//=============================================================================================================
447
448MatrixXf& MNEHemisphere::getTriCoords(float p_fScaling)
449{
450 if(m_TriCoords.size() == 0)
451 {
452 m_TriCoords = MatrixXf(3,3*itris.rows());
453 for(qint32 i = 0; i < itris.rows(); ++i)
454 {
455 m_TriCoords.col(i*3) = rr.row( itris(i,0) ).transpose().cast<float>();
456 m_TriCoords.col(i*3+1) = rr.row( itris(i,1) ).transpose().cast<float>();
457 m_TriCoords.col(i*3+2) = rr.row( itris(i,2) ).transpose().cast<float>();
458 }
459 }
460
461 m_TriCoords *= p_fScaling;
462
463 return m_TriCoords;
464}
465
466//=============================================================================================================
467
469{
470 FiffCoordTrans trans(p_Trans);
471
472 if (this->coord_frame == dest)
473 {
474// res = src;
475 return true;
476 }
477
478 if (trans.to == this->coord_frame && trans.from == dest)
479 trans.invert_transform();
480 else if(trans.from != this->coord_frame || trans.to != dest)
481 {
482 printf("Cannot transform the source space using this coordinate transformation");//Consider throw
483 return false;
484 }
485
486 MatrixXf t = trans.trans.block(0,0,3,4);
487// res = src;
488 this->coord_frame = dest;
489 MatrixXf t_rr = MatrixXf::Ones(this->np, 4);
490 t_rr.block(0, 0, this->np, 3) = this->rr;
491 MatrixXf t_nn = MatrixXf::Zero(this->np, 4);
492 t_nn.block(0, 0, this->np, 3) = this->nn;
493
494 this->rr = (t*t_rr.transpose()).transpose();
495 this->nn = (t*t_nn.transpose()).transpose();
496
497 return true;
498}
499
500//=============================================================================================================
501//ToDo
503{
504 if(this->type == 1 || this->type == 2)
505 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_TYPE, &this->type);
506 else
507 printf("Unknown source space type (%d)\n", this->type);
508 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_ID, &this->id);
509
510// data = this.get('subject_his_id', None)
511// if data:
512// write_string(fid, FIFF.FIFF_SUBJ_HIS_ID, data)
513 p_pStream->write_int(FIFF_MNE_COORD_FRAME, &this->coord_frame);
514
515 if(this->type == 2) //2 = Vol
516 {
517 qDebug() << "ToDo: Write Volume not implemented yet!!!!!!!!";
518// p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_VOXEL_DIMS, this->shape)
519// p_pStream->write_coord_trans(this->src_mri_t);
520
522// write_coord_trans(fid, this['vox_mri_t'])
523
524// write_coord_trans(fid, this['mri_ras_t'])
525
526// write_float_sparse_rcs(fid, FIFF.FIFF_MNE_SOURCE_SPACE_INTERPOLATOR,
527// this['interpolator'])
528
529// if 'mri_file' in this and this['mri_file'] is not None:
530// write_string(fid, FIFF.FIFF_MNE_SOURCE_SPACE_MRI_FILE,
531// this['mri_file'])
532
533// write_int(fid, FIFF.FIFF_MRI_WIDTH, this['mri_width'])
534// write_int(fid, FIFF.FIFF_MRI_HEIGHT, this['mri_height'])
535// write_int(fid, FIFF.FIFF_MRI_DEPTH, this['mri_depth'])
536
538 }
539
540 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_NPOINTS, &this->np);
543
544 // Which vertices are active
545 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_SELECTION, this->inuse.data(), this->inuse.size());
546 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_NUSE, &this->nuse);
547
548 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_NTRI, &this->ntri);
549 if (this->ntri > 0)
550 p_pStream->write_int_matrix(FIFF_MNE_SOURCE_SPACE_TRIANGLES, (this->itris.array() + 1).matrix());
551
552 if (this->type != 2 && this->use_itris.rows() > 0)
553 {
554 // Use triangulation
556 p_pStream->write_int_matrix(FIFF_MNE_SOURCE_SPACE_USE_TRIANGLES, (this->use_itris.array() + 1).matrix());
557 }
558
559 // Patch-related information
560 if (!this->nearest.empty())
561 {
562 Eigen::VectorXi nearestIdx = this->nearestVertIdx();
563 Eigen::VectorXf nearestDistF = this->nearestDistVec().cast<float>();
564 p_pStream->write_int(FIFF_MNE_SOURCE_SPACE_NEAREST, nearestIdx.data(), nearestIdx.size());
566 }
567
568 // Distances
569 if (!this->dist.is_empty())
570 {
571 // Convert FiffSparseMatrix to Eigen and save only upper triangular portion
572 Eigen::SparseMatrix<double> eigenDist = this->dist.toEigenSparse();
573 typedef Eigen::Triplet<float> T;
574 std::vector<T> tripletList;
575 tripletList.reserve(eigenDist.nonZeros());
576 for (int k=0; k < eigenDist.outerSize(); ++k)
577 for (Eigen::SparseMatrix<double>::InnerIterator it(eigenDist,k); it; ++it)
578 if(it.col() >= it.row())//only upper triangle -> todo iteration can be optimized
579 tripletList.push_back(T(it.row(), it.col(), (float)it.value()));
580 Eigen::SparseMatrix<float> dists(eigenDist.rows(), eigenDist.cols());
581 dists.setFromTriplets(tripletList.begin(), tripletList.end());
582
584 //ToDo check if write_float_matrix or write float is okay
585 p_pStream->write_float(FIFF_MNE_SOURCE_SPACE_DIST_LIMIT, &this->dist_limit); //p_pStream->write_float_matrix(FIFF_MNE_SOURCE_SPACE_DIST_LIMIT, this->dist_limit);
586 }
587}
588
590
591//QGeometryData* MNEHemisphere::getGeometryData(float p_fScaling)
592//{
593// if(m_pGeometryData == NULL)
594// {
595// m_pGeometryData = new QGeometryData();
596
597// MatrixXd* triCoords = getTriCoords(p_fScaling);
598
599// m_pGeometryData->appendVertexArray(QArray<QVector3D>::fromRawData( reinterpret_cast<const QVector3D*>(triCoords->data()), triCoords->cols() ));
600// }
601
602// return m_pGeometryData;
603//}
#define FIFF_MNE_SOURCE_SPACE_ID
#define FIFF_MNE_SOURCE_SPACE_DIST_LIMIT
#define FIFF_MNE_COORD_FRAME
#define FIFF_MNE_SOURCE_SPACE_TYPE
#define FIFF_MNE_SOURCE_SPACE_SELECTION
#define FIFF_MNE_SOURCE_SPACE_USE_TRIANGLES
#define FIFF_MNE_SOURCE_SPACE_NUSE_TRI
#define FIFF_MNE_SOURCE_SPACE_NORMALS
#define FIFF_MNE_SOURCE_SPACE_NEAREST_DIST
#define FIFF_MNE_SOURCE_SPACE_DIST
#define FIFF_MNE_SOURCE_SPACE_POINTS
#define FIFF_MNE_SOURCE_SPACE_NTRI
#define FIFF_MNE_SOURCE_SPACE_NPOINTS
#define FIFF_MNE_SOURCE_SPACE_TRIANGLES
#define FIFF_MNE_SOURCE_SPACE_NUSE
#define FIFF_MNE_SOURCE_SPACE_NEAREST
#define FIFFB_MNE_PARENT_MRI_FILE
MNEMath class declaration.
MNEHemisphere class declaration.
MNENearest class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
qint32 fiff_int_t
Definition fiff_types.h:89
Shared utilities (I/O helpers, spectral analysis, layout management, warp algorithms).
Definition buildinfo.h:45
Coordinate transformation description.
Eigen::Matrix< float, 4, 4, Eigen::DontAlign > trans
FIFF sparse matrix storage.
FIFF File I/O routines.
fiff_long_t start_block(fiff_int_t kind)
fiff_long_t write_float_matrix(fiff_int_t kind, const Eigen::MatrixXf &mat)
fiff_long_t write_int_matrix(fiff_int_t kind, const Eigen::MatrixXi &mat)
fiff_long_t write_int(fiff_int_t kind, const fiff_int_t *data, fiff_int_t nel=1, fiff_int_t next=FIFFV_NEXT_SEQ)
fiff_long_t write_float(fiff_int_t kind, const float *data, fiff_int_t nel=1)
fiff_long_t write_float_sparse_rcs(fiff_int_t kind, const Eigen::SparseMatrix< float > &mat)
fiff_long_t end_block(fiff_int_t kind, fiff_int_t next=FIFFV_NEXT_SEQ)
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
QList< Eigen::VectorXi > pinfo
Eigen::VectorXi patch_inds
std::shared_ptr< MNESourceSpace > SPtr
std::vector< Eigen::VectorXi > neighbor_tri
std::vector< Eigen::VectorXi > neighbor_vert
Eigen::VectorXd nearestDistVec() const
FIFFLIB::FiffSparseMatrix dist
std::vector< MNENearest > nearest
Eigen::VectorXi nearestVertIdx() const
static bool compareIdxValuePairSmallerThan(const std::pair< int, T > &lhs, const std::pair< int, T > &rhs)
Definition mnemath.h:580