v2.0.0
Loading...
Searching...
No Matches
mne_surface_or_volume.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
42#include "mne_surface.h"
43#include "mne_source_space.h"
44#include "mne_patch_info.h"
45//#include "fwd_bem_model.h"
46#include "mne_nearest.h"
47#include "filter_thread_arg.h"
48#include "mne_triangle.h"
50#include "mne_proj_data.h"
51#include "mne_vol_geom.h"
52#include "mne_mgh_tag_group.h"
53#include "mne_mgh_tag.h"
54
55#include <fiff/fiff_stream.h>
58#include <fiff/fiff_dig_point.h>
59
60#include <math/sphere.h>
61#include <utils/ioutils.h>
62
63#include <QFile>
64#include <QTextStream>
65#include <QCoreApplication>
66#include <QtConcurrent>
67#include <QDebug>
68
69#define _USE_MATH_DEFINES
70#include <math.h>
71
72#include <Eigen/Dense>
73#include <Eigen/Sparse>
74
75#include <algorithm>
76
77
78
79//=============================================================================================================
80// USED NAMESPACES
81//=============================================================================================================
82
83using namespace Eigen;
84using namespace FIFFLIB;
85using namespace MNELIB;
86
87//============================= dot.h =============================
88
89constexpr int FAIL = -1;
90constexpr int OK = 0;
91
92constexpr int X = 0;
93constexpr int Y = 1;
94constexpr int Z = 2;
95
96constexpr int NNEIGHBORS = 26;
97
98
99
100
101//============================= make_volume_source_space.c =============================
102
103static std::optional<FiffCoordTrans> make_voxel_ras_trans(const Eigen::Vector3f& r0,
104 const Eigen::Vector3f& x_ras,
105 const Eigen::Vector3f& y_ras,
106 const Eigen::Vector3f& z_ras,
107 const Eigen::Vector3f& voxel_size)
108{
109 Eigen::Matrix3f rot;
110 rot.row(0) = x_ras.transpose() * voxel_size[0];
111 rot.row(1) = y_ras.transpose() * voxel_size[1];
112 rot.row(2) = z_ras.transpose() * voxel_size[2];
113
115}
116
117//=============================================================================================================
118// DEFINE MEMBER METHODS
119//=============================================================================================================
120
124
125//=============================================================================================================
126
128{
129 if (curv.size() > 0)
130 return;
131 curv = Eigen::VectorXf::Ones(np);
132}
133
134//=============================================================================================================
135
137{
138 // rr, nn, itris, use_itris are Eigen matrices — auto-cleanup
139 // inuse, vertno are Eigen VectorXi — auto-cleanup
140 // tris, use_tris are std::vector<MNETriangle> — auto-cleanup
141 // neighbor_tri is std::vector<VectorXi> — auto-cleanup
142 // nneighbor_tri is Eigen VectorXi — auto-cleanup
143 // curv is Eigen VectorXf — auto-cleanup
144
145 // neighbor_vert is std::vector<VectorXi> — auto-cleanup
146 // nneighbor_vert is Eigen VectorXi — auto-cleanup
147 // vert_dist is std::vector<VectorXf> — auto-cleanup
148 // nearest is std::vector<MNENearest> — auto-cleanup
149 // patches is std::vector<optional<MNEPatchInfo>> — auto-cleanup
150 // dist is FiffSparseMatrix value, interpolator/vol_geom/mgh_tags are optional — auto-cleanup
151 // voxel_surf_RAS_t, MRI_voxel_surf_RAS_t, MRI_surf_RAS_RAS_t are optional — auto-cleanup
152 this->MRI_volume.clear();
153}
154
155//=============================================================================================================
156
158{
159 Eigen::VectorXi result(static_cast<int>(nearest.size()));
160 for (int i = 0; i < result.size(); ++i)
161 result[i] = nearest[i].nearest;
162 return result;
163}
164
165//=============================================================================================================
166
168{
169 Eigen::VectorXd result(static_cast<int>(nearest.size()));
170 for (int i = 0; i < result.size(); ++i)
171 result[i] = static_cast<double>(nearest[i].dist);
172 return result;
173}
174
175//=============================================================================================================
176
177void MNESurfaceOrVolume::setNearestData(const Eigen::VectorXi& nearestIdx, const Eigen::VectorXd& nearestDist)
178{
179 const int n = nearestIdx.size();
180 nearest.resize(n);
181 for (int i = 0; i < n; ++i) {
182 nearest[i].vert = i;
183 nearest[i].nearest = nearestIdx[i];
184 nearest[i].dist = static_cast<float>(nearestDist[i]);
185 nearest[i].patch = nullptr;
186 }
187}
188
189//=============================================================================================================
190
191double MNESurfaceOrVolume::solid_angle(const Eigen::Vector3f& from, const MNETriangle& tri) /* ...to this triangle */
192/*
193 * Compute the solid angle according to van Oosterom's
194 * formula
195 */
196{
197 Eigen::Vector3d v1 = (tri.r1 - from).cast<double>();
198 Eigen::Vector3d v2 = (tri.r2 - from).cast<double>();
199 Eigen::Vector3d v3 = (tri.r3 - from).cast<double>();
200
201 double triple = v1.cross(v2).dot(v3);
202
203 double l1 = v1.norm();
204 double l2 = v2.norm();
205 double l3 = v3.norm();
206 double s = (l1*l2*l3+v1.dot(v2)*l3+v1.dot(v3)*l2+v2.dot(v3)*l1);
207
208 return (2.0*atan2(triple,s));
209}
210
211//=============================================================================================================
212
214/*
215 * Add the triangle data structures
216 */
217{
218 int k;
219 MNETriangle* tri;
220
222 return;
223
224 tris.clear();
225 use_tris.clear();
226 /*
227 * Add information for the complete triangulation
228 */
229 if (itris.rows() > 0 && ntri > 0) {
230 tris.resize(ntri);
231 tot_area = 0.0;
232 for (k = 0, tri = tris.data(); k < ntri; k++, tri++) {
233 tri->vert = &itris(k,0);
234 tri->r1 = rr.row(tri->vert[0]).transpose();
235 tri->r2 = rr.row(tri->vert[1]).transpose();
236 tri->r3 = rr.row(tri->vert[2]).transpose();
237 tri->compute_data();
238 tot_area += tri->area;
239 }
240#ifdef TRIANGLE_SIZE_WARNING
241 for (k = 0, tri = tris.data(); k < ntri; k++, tri++)
242 if (tri->area < 1e-5*tot_area/ntri)
243 qWarning("Warning: Triangle area is only %g um^2 (%.5f %% of expected average)\n",
244 1e12*tri->area,100*ntri*tri->area/tot_area);
245#endif
246 }
247#ifdef DEBUG
248 qInfo("\ttotal area = %-.1f cm^2\n",1e4*tot_area);
249#endif
250 /*
251 * Add information for the selected subset if applicable
252 */
253 if (use_itris.rows() > 0 && nuse_tri > 0) {
254 use_tris.resize(nuse_tri);
255 for (k = 0, tri = use_tris.data(); k < nuse_tri; k++, tri++) {
256 tri->vert = &use_itris(k,0);
257 tri->r1 = rr.row(tri->vert[0]).transpose();
258 tri->r2 = rr.row(tri->vert[1]).transpose();
259 tri->r3 = rr.row(tri->vert[2]).transpose();
260 tri->compute_data();
261 }
262 }
263 return;
264}
265
266//=============================================================================================================
267
269/*
270 * Compute the center of mass of a set of points
271 */
272{
273 int q;
274 cm[0] = cm[1] = cm[2] = 0.0;
275 for (q = 0; q < np; q++) {
276 cm[0] += rr(q,0);
277 cm[1] += rr(q,1);
278 cm[2] += rr(q,2);
279 }
280 if (np > 0) {
281 cm[0] = cm[0]/np;
282 cm[1] = cm[1]/np;
283 cm[2] = cm[2]/np;
284 }
285 return;
286}
287
288//=============================================================================================================
289
291/*
292 * Compute the center of mass of a surface
293 */
294{
296 return;
297}
298
299//=============================================================================================================
300
302{
303 int k,p,ndist;
304 int nneigh;
305
306 if (neighbor_vert.empty() || nneighbor_vert.size() == 0)
307 return;
308
309 vert_dist.clear();
310 vert_dist.resize(np);
311 qInfo("\tDistances between neighboring vertices...");
312 for (k = 0, ndist = 0; k < np; k++) {
313 nneigh = nneighbor_vert[k];
314 vert_dist[k] = Eigen::VectorXf(nneigh);
315 const Eigen::VectorXi& neigh = neighbor_vert[k];
316 for (p = 0; p < nneigh; p++) {
317 if (neigh[p] >= 0) {
318 vert_dist[k][p] = (rr.row(neigh[p]) - rr.row(k)).norm();
319 }
320 else
321 vert_dist[k][p] = -1.0;
322 ndist++;
323 }
324 }
325 qInfo("[%d distances done]\n",ndist);
326 return;
327}
328
329//=============================================================================================================
330
332{
333 int k,c,p;
334 int *ii;
335 float w,size;
336 MNETriangle* tri;
337
339 return OK;
340 /*
341 * Reallocate the stuff and initialize
342 */
343 nn = MNESurfaceOrVolume::NormalsT::Zero(np,3);
344 /*
345 * One pass through the triangles will do it
346 */
348 for (p = 0, tri = tris.data(); p < ntri; p++, tri++) {
349 ii = tri->vert;
350 w = 1.0; /* This should be related to the triangle size */
351 /*
352 * Then the vertex normals
353 */
354 for (k = 0; k < 3; k++)
355 for (c = 0; c < 3; c++)
356 nn(ii[k],c) += w*tri->nn[c];
357 }
358 for (k = 0; k < np; k++) {
359 size = nn.row(k).norm();
360 if (size > 0.0)
361 nn.row(k) /= size;
362 }
364 return OK;
365}
366
367//=============================================================================================================
368
369int MNESurfaceOrVolume::add_geometry_info(bool do_normals, bool check_too_many_neighbors)
370/*
371 * Add vertex normals and neighbourhood information
372 */
373{
374 int k,c,p,q;
375 int vert;
376 int *ii;
377 int nneighbors;
378 float w,size;
379 int found;
380 int nfix_distinct,nfix_no_neighbors,nfix_defect;
381 MNETriangle* tri;
382
385 return OK;
386 }
388 return OK;
389 /*
390 * Reallocate the stuff and initialize
391 */
392 if (do_normals) {
393 nn = MNESurfaceOrVolume::NormalsT::Zero(np,3);
394 }
395 neighbor_tri.clear();
396 neighbor_tri.resize(np);
397 nneighbor_tri = Eigen::VectorXi::Zero(np);
398
399 /* nn is already zero-initialized above */
400 /*
401 * One pass through the triangles will do it
402 */
404 for (p = 0, tri = tris.data(); p < ntri; p++, tri++)
405 if (tri->area == 0)
406 qWarning("\tWarning : zero size triangle # %d\n",p);
407 qInfo("\tTriangle ");
408 if (do_normals)
409 qInfo("and vertex ");
410 qInfo("normals and neighboring triangles...");
411 for (p = 0, tri = tris.data(); p < ntri; p++, tri++) {
412 ii = tri->vert;
413 w = 1.0; /* This should be related to the triangle size */
414 for (k = 0; k < 3; k++) {
415 /*
416 * Then the vertex normals
417 */
418 if (do_normals)
419 for (c = 0; c < 3; c++)
420 nn(ii[k],c) += w*tri->nn[c];
421 /*
422 * Add to the list of neighbors
423 */
424 neighbor_tri[ii[k]].conservativeResize(nneighbor_tri[ii[k]]+1);
425 neighbor_tri[ii[k]][nneighbor_tri[ii[k]]] = p;
426 nneighbor_tri[ii[k]]++;
427 }
428 }
429 nfix_no_neighbors = 0;
430 nfix_defect = 0;
431 for (k = 0; k < np; k++) {
432 if (nneighbor_tri[k] <= 0) {
433#ifdef STRICT_ERROR
434 err_printf_set_error("Vertex %d does not have any neighboring triangles!",k);
435 return FAIL;
436#else
437#ifdef REPORT_WARNINGS
438 qWarning("Warning: Vertex %d does not have any neighboring triangles!\n",k);
439#endif
440#endif
441 nfix_no_neighbors++;
442 }
443 else if (nneighbor_tri[k] < 3) {
444#ifdef REPORT_WARNINGS
445 qWarning("\n\tTopological defect: Vertex %d has only %d neighboring triangle%s Vertex omitted.\n\t",
446 k,nneighbor_tri[k],nneighbor_tri[k] > 1 ? "s." : ".");
447#endif
448 nfix_defect++;
449 nneighbor_tri[k] = 0;
450 neighbor_tri[k].resize(0);
451 }
452 }
453 /*
454 * Scale the vertex normals to unit length
455 */
456 for (k = 0; k < np; k++)
457 if (nneighbor_tri[k] > 0) {
458 size = nn.row(k).norm();
459 if (size > 0.0)
460 nn.row(k) /= size;
461 }
462 qInfo("[done]\n");
463 /*
464 * Determine the neighboring vertices
465 */
466 qInfo("\tVertex neighbors...");
467 neighbor_vert.clear();
468 neighbor_vert.resize(np);
469 nneighbor_vert = VectorXi::Zero(np);
470 /*
471 * We know the number of neighbors beforehand
472 */
473 for (k = 0; k < np; k++) {
474 if (nneighbor_tri[k] > 0) {
475 neighbor_vert[k] = VectorXi(nneighbor_tri[k]);
477 }
478 else {
479 nneighbor_vert[k] = 0;
480 }
481 }
482 nfix_distinct = 0;
483 for (k = 0; k < np; k++) {
484 Eigen::VectorXi& neighbors = neighbor_vert[k];
485 nneighbors = 0;
486 for (p = 0; p < nneighbor_tri[k]; p++) {
487 /*
488 * Fit in the other vertices of the neighboring triangle
489 */
490 for (c = 0; c < 3; c++) {
491 vert = tris[neighbor_tri[k][p]].vert[c];
492 if (vert != k) {
493 for (q = 0, found = false; q < nneighbors; q++) {
494 if (neighbors[q] == vert) {
495 found = true;
496 break;
497 }
498 }
499 if (!found) {
500 if (nneighbors < nneighbor_vert[k])
501 neighbors[nneighbors++] = vert;
502 else {
503 if (check_too_many_neighbors) {
504 qCritical("Too many neighbors for vertex %d.",k);
505 return FAIL;
506 }
507 else
508 qWarning("\tWarning: Too many neighbors for vertex %d\n",k);
509 }
510 }
511 }
512 }
513 }
514 if (nneighbors != nneighbor_vert[k]) {
515#ifdef REPORT_WARNINGS
516 qWarning("\n\tIncorrect number of distinct neighbors for vertex %d (%d instead of %d) [fixed].",
517 k,nneighbors,nneighbor_vert[k]);
518#endif
519 nfix_distinct++;
520 nneighbor_vert[k] = nneighbors;
521 }
522 }
523 qInfo("[done]\n");
524 /*
525 * Distance calculation follows
526 */
529 /*
530 * Summarize the defects
531 */
532 if (nfix_defect > 0)
533 qWarning("\tWarning: %d topological defects were fixed.\n",nfix_defect);
534 if (nfix_distinct > 0)
535 qWarning("\tWarning: %d vertices had incorrect number of distinct neighbors (fixed).\n",nfix_distinct);
536 if (nfix_no_neighbors > 0)
537 qWarning("\tWarning: %d vertices did not have any neighboring triangles (fixed)\n",nfix_no_neighbors);
538#ifdef DEBUG
539 for (k = 0; k < np; k++) {
540 if (nneighbor_vert[k] <= 0)
541 qCritical("No neighbors for vertex %d\n",k);
542 if (nneighbor_tri[k] <= 0)
543 qCritical("No neighbor tris for vertex %d\n",k);
544 }
545#endif
546 return OK;
547}
548
549//=============================================================================================================
550
552{
553 return add_geometry_info(do_normals,true);
554}
555
556//=============================================================================================================
557
559
560{
561 return add_geometry_info(do_normals,false);
562}
563
564//=============================================================================================================
IOUtils class declaration.
Sphere class declaration.
constexpr int FAIL
constexpr int Y
constexpr int Z
constexpr int OK
constexpr int X
FiffDigitizerData class declaration.
FiffCoordTrans class declaration.
FiffStream class declaration.
#define FIFFV_MNE_COORD_MRI_VOXEL
#define FIFFV_COORD_MRI
FiffDigPoint class declaration.
#define MNE_SOURCE_SPACE_SURFACE
Definition mne_types.h:103
#define MNE_SOURCE_SPACE_VOLUME
Definition mne_types.h:104
MNESourceSpace class declaration.
constexpr int NNEIGHBORS
MNEProjData class declaration.
MNE Patch Information (MNEPatchInfo) class declaration.
MNENearest class declaration.
MNE FsSurface or Volume (MNESurfaceOrVolume) class declaration.
MNEMshDisplaySurface class declaration.
MNEMghTagGroup class declaration.
MNEMghTag class declaration.
MNETriangle class declaration.
Filter Thread Argument (FilterThreadArg) class declaration.
MNESurface class declaration.
MNEVolGeom class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
Coordinate transformation description.
std::vector< Eigen::VectorXi > neighbor_tri
void setNearestData(const Eigen::VectorXi &nearestIdx, const Eigen::VectorXd &nearestDist)
std::vector< Eigen::VectorXi > neighbor_vert
MNESurfaceOrVolume()
Constructs the MNE FsSurface or Volume.
Eigen::VectorXd nearestDistVec() const
int add_geometry_info(bool do_normals, bool check_too_many_neighbors)
FIFFLIB::FiffSparseMatrix dist
std::vector< MNENearest > nearest
static void compute_cm(const PointsT &rr, int np, float(&cm)[3])
static double solid_angle(const Eigen::Vector3f &from, const MNELIB::MNETriangle &tri)
virtual ~MNESurfaceOrVolume()
Destroys the MNE FsSurface or Volume description.
std::vector< MNETriangle > tris
std::vector< Eigen::VectorXf > vert_dist
Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > PointsT
int add_geometry_info2(bool do_normals)
Eigen::VectorXi nearestVertIdx() const
std::vector< MNETriangle > use_tris
Per-triangle geometric data for a cortical or BEM surface.
Eigen::Vector3f nn
Eigen::Vector3f r2
Eigen::Vector3f r1
Eigen::Vector3f r3