v2.0.0
Loading...
Searching...
No Matches
mne_surface.cpp
Go to the documentation of this file.
1//=============================================================================================================
36
37//=============================================================================================================
38// INCLUDES
39//=============================================================================================================
40
41#include "mne_surface.h"
42#include "mne_source_space.h"
43#include "mne_triangle.h"
44#include "mne_proj_data.h"
45
46#include <fiff/fiff_stream.h>
48
49#include <QFile>
50#include <QDebug>
51
52#define _USE_MATH_DEFINES
53#include <math.h>
54
55#include <Eigen/Dense>
56
57//=============================================================================================================
58// USED NAMESPACES
59//=============================================================================================================
60
61using namespace Eigen;
62using namespace FIFFLIB;
63using namespace MNELIB;
64
65//=============================================================================================================
66
67constexpr int FAIL = -1;
68constexpr int OK = 0;
69
70constexpr int X = 0;
71constexpr int Y = 1;
72constexpr int Z = 2;
73
74//=============================================================================================================
75// DEFINE MEMBER METHODS
76//=============================================================================================================
77
81
82//=============================================================================================================
83
87
88//=============================================================================================================
89// Const geometry methods
90//=============================================================================================================
91
92double MNESurface::sum_solids(const Eigen::Vector3f& from) const
93{
94 int k;
95 double tot_angle, angle;
96 for (k = 0, tot_angle = 0.0; k < ntri; k++) {
97 angle = solid_angle(from, tris[k]);
98 tot_angle += angle;
99 }
100 return tot_angle;
101}
102
103//=============================================================================================================
104
105void MNESurface::triangle_coords(const Eigen::Vector3f& r, int tri, float &x, float &y, float &z) const
106{
107 double a, b, c, v1, v2, det;
108 const MNETriangle* this_tri;
109
110 this_tri = &tris[tri];
111
112 Eigen::Vector3d rr = (r - this_tri->r1).cast<double>();
113 z = rr.dot(this_tri->nn.cast<double>());
114
115 a = this_tri->r12.cast<double>().squaredNorm();
116 b = this_tri->r13.cast<double>().squaredNorm();
117 c = this_tri->r12.cast<double>().dot(this_tri->r13.cast<double>());
118
119 v1 = rr.dot(this_tri->r12.cast<double>());
120 v2 = rr.dot(this_tri->r13.cast<double>());
121
122 det = a * b - c * c;
123
124 x = (b * v1 - c * v2) / det;
125 y = (a * v2 - c * v1) / det;
126}
127
128//=============================================================================================================
129
130int MNESurface::nearest_triangle_point(const Eigen::Vector3f& r, const MNEProjData *user, int tri, float &x, float &y, float &z) const
131{
132 double p, q, p0, q0, t0;
133 double a, b, c, v1, v2, det;
134 double best, dist, dist0;
135 const MNEProjData* pd = user;
136 const MNETriangle* this_tri;
137
138 this_tri = &tris[tri];
139 Eigen::Vector3d rr = (r - this_tri->r1).cast<double>();
140 dist = rr.dot(this_tri->nn.cast<double>());
141
142 if (pd) {
143 if (!pd->act[tri])
144 return false;
145 a = pd->a[tri];
146 b = pd->b[tri];
147 c = pd->c[tri];
148 }
149 else {
150 a = this_tri->r12.cast<double>().squaredNorm();
151 b = this_tri->r13.cast<double>().squaredNorm();
152 c = this_tri->r12.cast<double>().dot(this_tri->r13.cast<double>());
153 }
154
155 v1 = rr.dot(this_tri->r12.cast<double>());
156 v2 = rr.dot(this_tri->r13.cast<double>());
157
158 det = a * b - c * c;
159
160 p = (b * v1 - c * v2) / det;
161 q = (a * v2 - c * v1) / det;
162
163 if (p >= 0.0 && p <= 1.0 &&
164 q >= 0.0 && q <= 1.0 &&
165 q <= 1.0 - p) {
166 x = p;
167 y = q;
168 z = dist;
169 return true;
170 }
171 /*
172 * Side 1 -> 2
173 */
174 p0 = p + 0.5 * (q * c) / a;
175 if (p0 < 0.0) p0 = 0.0;
176 else if (p0 > 1.0) p0 = 1.0;
177 q0 = 0.0;
178 dist0 = sqrt((p - p0) * (p - p0) * a +
179 (q - q0) * (q - q0) * b +
180 (p - p0) * (q - q0) * c +
181 dist * dist);
182 best = dist0;
183 x = p0;
184 y = q0;
185 z = dist0;
186 /*
187 * Side 2 -> 3
188 */
189 t0 = 0.5 * ((2.0 * a - c) * (1.0 - p) + (2.0 * b - c) * q) / (a + b - c);
190 if (t0 < 0.0) t0 = 0.0;
191 else if (t0 > 1.0) t0 = 1.0;
192 p0 = 1.0 - t0;
193 q0 = t0;
194 dist0 = sqrt((p - p0) * (p - p0) * a +
195 (q - q0) * (q - q0) * b +
196 (p - p0) * (q - q0) * c +
197 dist * dist);
198 if (dist0 < best) {
199 best = dist0;
200 x = p0;
201 y = q0;
202 z = dist0;
203 }
204 /*
205 * Side 1 -> 3
206 */
207 p0 = 0.0;
208 q0 = q + 0.5 * (p * c) / b;
209 if (q0 < 0.0) q0 = 0.0;
210 else if (q0 > 1.0) q0 = 1.0;
211 dist0 = sqrt((p - p0) * (p - p0) * a +
212 (q - q0) * (q - q0) * b +
213 (p - p0) * (q - q0) * c +
214 dist * dist);
215 if (dist0 < best) {
216 best = dist0;
217 x = p0;
218 y = q0;
219 z = dist0;
220 }
221 return true;
222}
223
224//=============================================================================================================
225
226int MNESurface::nearest_triangle_point(const Eigen::Vector3f& r, int tri, float &x, float &y, float &z) const
227{
228 return nearest_triangle_point(r, nullptr, tri, x, y, z);
229}
230
231//=============================================================================================================
232
233Eigen::Vector3f MNESurface::project_to_triangle(int tri, float p, float q) const
234{
235 const MNETriangle* this_tri = &tris[tri];
236
237 return Eigen::Vector3f(
238 this_tri->r1 + p * this_tri->r12 + q * this_tri->r13
239 );
240}
241
242//=============================================================================================================
243
244Eigen::Vector3f MNESurface::project_to_triangle(int best, const Eigen::Vector3f& r) const
245{
246 float p, q, dist;
247 nearest_triangle_point(r, best, p, q, dist);
248 return project_to_triangle(best, p, q);
249}
250
251//=============================================================================================================
252
253int MNESurface::project_to_surface(const MNEProjData *proj_data, const Eigen::Vector3f& r, float &distp) const
254{
255 float dist;
256 float p, q;
257 float p0, q0, dist0;
258 int best;
259 int k;
260
261 p0 = q0 = 0.0;
262 dist0 = 0.0;
263 for (best = -1, k = 0; k < ntri; k++) {
264 if (nearest_triangle_point(r, proj_data, k, p, q, dist)) {
265 if (best < 0 || std::fabs(dist) < std::fabs(dist0)) {
266 dist0 = dist;
267 best = k;
268 p0 = p;
269 q0 = q;
270 }
271 }
272 }
273 distp = dist0;
274 return best;
275}
276
277//=============================================================================================================
278
280 Eigen::VectorXi& nearest_tri,
281 Eigen::VectorXf& dist, int nstep) const
282{
283 auto p = std::make_unique<MNEProjData>(this);
284 int k, was;
285
286 qInfo("%s for %d points %d steps...", nearest_tri[0] < 0 ? "Closest" : "Approx closest", np_points, nstep);
287
288 for (k = 0; k < np_points; k++) {
289 was = nearest_tri[k];
290 Eigen::Vector3f pt = Eigen::Map<const Eigen::Vector3f>(r.row(k).data());
291 decide_search_restriction(*p, nearest_tri[k], nstep, pt);
292 nearest_tri[k] = project_to_surface(p.get(), pt, dist[k]);
293 if (nearest_tri[k] < 0) {
294 decide_search_restriction(*p, -1, nstep, pt);
295 nearest_tri[k] = project_to_surface(p.get(), pt, dist[k]);
296 }
297 }
298 (void)was;
299
300 qInfo("[done]\n");
301}
302
303//=============================================================================================================
304
306 int approx_best,
307 int nstep,
308 const Eigen::Vector3f& r) const
309{
310 int k;
311 Eigen::Vector3f diff_vec;
312 float dist_val, mindist;
313 int minvert;
314
315 for (k = 0; k < ntri; k++)
316 p.act[k] = 0;
317
318 if (approx_best < 0) {
319 mindist = 1000.0;
320 minvert = 0;
321 for (k = 0; k < np; k++) {
322 diff_vec = rr.row(k).transpose() - r;
323 dist_val = diff_vec.norm();
324 if (dist_val < mindist && nneighbor_tri[k] > 0) {
325 mindist = dist_val;
326 minvert = k;
327 }
328 }
329 }
330 else {
331 const MNETriangle* this_tri = &tris[approx_best];
332 diff_vec = this_tri->r1 - r;
333 mindist = diff_vec.norm();
334 minvert = this_tri->vert[0];
335
336 diff_vec = this_tri->r2 - r;
337 dist_val = diff_vec.norm();
338 if (dist_val < mindist) {
339 mindist = dist_val;
340 minvert = this_tri->vert[1];
341 }
342 diff_vec = this_tri->r3 - r;
343 dist_val = diff_vec.norm();
344 if (dist_val < mindist) {
345 mindist = dist_val;
346 minvert = this_tri->vert[2];
347 }
348 }
349
350 activate_neighbors(minvert, p.act, nstep);
351
352 for (k = 0, p.nactive = 0; k < ntri; k++)
353 if (p.act[k])
354 p.nactive++;
355}
356
357//=============================================================================================================
358
359void MNESurface::activate_neighbors(int start, Eigen::VectorXi &act, int nstep) const
360{
361 int k;
362
363 if (nstep == 0)
364 return;
365
366 for (k = 0; k < nneighbor_tri[start]; k++)
367 act[neighbor_tri[start][k]] = 1;
368 for (k = 0; k < nneighbor_vert[start]; k++)
369 activate_neighbors(neighbor_vert[start][k], act, nstep - 1);
370}
371
372//=============================================================================================================
373// Static factory methods
374//=============================================================================================================
375
376std::unique_ptr<MNESurface> MNESurface::read_bem_surface(const QString& name, int which, bool add_geometry)
377{
378 float sigma;
379 return read_bem_surface(name, which, add_geometry, sigma, true);
380}
381
382//=============================================================================================================
383
384std::unique_ptr<MNESurface> MNESurface::read_bem_surface(const QString& name, int which, bool add_geometry, float& sigma)
385{
386 return read_bem_surface(name, which, add_geometry, sigma, true);
387}
388
389//=============================================================================================================
390
391std::unique_ptr<MNESurface> MNESurface::read_bem_surface2(const QString& name, int which, bool add_geometry)
392{
393 float sigma;
394 return read_bem_surface(name, which, add_geometry, sigma, false);
395}
396
397//=============================================================================================================
398
399std::unique_ptr<MNESurface> MNESurface::read_bem_surface(const QString& name, int which, bool add_geometry, float& sigma, bool check_too_many_neighbors)
400{
401 QFile file(name);
402 FiffStream::SPtr stream(new FiffStream(&file));
403
404 QList<FiffDirNode::SPtr> surfs;
405 QList<FiffDirNode::SPtr> bems;
407 FiffTag::UPtr t_pTag;
408
409 int id = -1;
410 int nnode, ntri_count;
411 MNESurface* s = nullptr;
412 std::unique_ptr<MNESurface> s_ptr;
413 int k;
414 MatrixXf tmp_node_normals;
416 float sigmaLocal = -1.0;
417 MatrixXf tmp_nodes;
418 MatrixXi tmp_triangles;
419
420 if (!stream->open())
421 return nullptr;
422
423 bems = stream->dirtree()->dir_tree_find(FIFFB_BEM);
424 if (bems.size() > 0) {
425 node = bems[0];
426 if (node->find_tag(stream, FIFF_BEM_COORD_FRAME, t_pTag)) {
427 coord_frame = *t_pTag->toInt();
428 }
429 }
430 surfs = stream->dirtree()->dir_tree_find(FIFFB_BEM_SURF);
431 if (surfs.size() == 0) {
432 qCritical("No BEM surfaces found in %s", name.toUtf8().constData());
433 stream->close(); return nullptr;
434 }
435 if (which >= 0) {
436 for (k = 0; k < surfs.size(); ++k) {
437 node = surfs[k];
438 if (node->find_tag(stream, FIFF_BEM_SURF_ID, t_pTag)) {
439 id = *t_pTag->toInt();
440 if (id == which)
441 break;
442 }
443 }
444 if (id != which) {
445 qCritical("Desired surface not found in %s", name.toUtf8().constData());
446 stream->close(); return nullptr;
447 }
448 }
449 else
450 node = surfs[0];
451
452 if (!node->find_tag(stream, FIFF_BEM_SURF_NNODE, t_pTag)) {
453 stream->close(); return nullptr;
454 }
455 nnode = *t_pTag->toInt();
456
457 if (!node->find_tag(stream, FIFF_BEM_SURF_NTRI, t_pTag)) {
458 stream->close(); return nullptr;
459 }
460 ntri_count = *t_pTag->toInt();
461
462 if (!node->find_tag(stream, FIFF_BEM_SURF_NODES, t_pTag)) {
463 stream->close(); return nullptr;
464 }
465 tmp_nodes = t_pTag->toFloatMatrix().transpose();
466
467 if (node->find_tag(stream, FIFF_BEM_SURF_NORMALS, t_pTag)) {
468 tmp_node_normals = t_pTag->toFloatMatrix().transpose();
469 }
470
471 if (!node->find_tag(stream, FIFF_BEM_SURF_TRIANGLES, t_pTag)) {
472 stream->close(); return nullptr;
473 }
474 tmp_triangles = t_pTag->toIntMatrix().transpose();
475
476 if (node->find_tag(stream, FIFF_MNE_COORD_FRAME, t_pTag)) {
477 coord_frame = *t_pTag->toInt();
478 }
479 else if (node->find_tag(stream, FIFF_BEM_COORD_FRAME, t_pTag)) {
480 coord_frame = *t_pTag->toInt();
481 }
482 if (node->find_tag(stream, FIFF_BEM_SIGMA, t_pTag)) {
483 sigmaLocal = *t_pTag->toFloat();
484 }
485
486 stream->close();
487
488 s_ptr = std::make_unique<MNESurface>();
489 s = s_ptr.get();
490 tmp_triangles.array() -= 1;
491 s->itris = tmp_triangles;
492 s->id = which;
493 s->sigma = sigmaLocal;
495 s->rr = tmp_nodes;
496 if (tmp_node_normals.rows() > 0)
497 s->nn = tmp_node_normals;
498 s->ntri = ntri_count;
499 s->np = nnode;
501 s->nuse_tri = 0;
502 s->tot_area = 0.0;
503 s->dist_limit = -1.0;
504 s->vol_dims[0] = s->vol_dims[1] = s->vol_dims[2] = 0;
505 s->MRI_vol_dims[0] = s->MRI_vol_dims[1] = s->MRI_vol_dims[2] = 0;
506 s->cm[0] = s->cm[1] = s->cm[2] = 0.0;
507
508 if (add_geometry) {
509 if (check_too_many_neighbors) {
510 if (s->add_geometry_info(s->nn.rows() == 0) != OK) {
511 return nullptr;
512 }
513 }
514 else {
515 if (s->add_geometry_info2(s->nn.rows() == 0) != OK) {
516 return nullptr;
517 }
518 }
519 }
520 else if (s->nn.rows() == 0) {
521 if (s->add_vertex_normals() != OK) {
522 return nullptr;
523 }
524 }
525 else
527
528 s->nuse = s->np;
529 s->inuse = Eigen::VectorXi::Ones(s->np);
530 s->vertno = Eigen::VectorXi::LinSpaced(s->np, 0, s->np - 1);
531 sigma = sigmaLocal;
532
533 return s_ptr;
534}
constexpr int FAIL
constexpr int Y
constexpr int Z
constexpr int OK
constexpr int X
FiffCoordTrans class declaration.
#define FIFFB_BEM_SURF
Definition fiff_file.h:404
#define FIFF_BEM_SURF_NNODE
Definition fiff_file.h:733
#define FIFF_BEM_SURF_NODES
Definition fiff_file.h:735
#define FIFF_BEM_SURF_TRIANGLES
Definition fiff_file.h:736
#define FIFF_BEM_COORD_FRAME
Definition fiff_file.h:743
#define FIFFB_BEM
Definition fiff_file.h:403
#define FIFF_BEM_SURF_ID
Definition fiff_file.h:731
#define FIFF_BEM_SURF_NTRI
Definition fiff_file.h:734
#define FIFF_BEM_SIGMA
Definition fiff_file.h:744
#define FIFF_BEM_SURF_NORMALS
Definition fiff_file.h:737
FiffStream class declaration.
#define FIFF_MNE_COORD_FRAME
#define FIFFV_MNE_SPACE_SURFACE
#define FIFFV_COORD_MRI
MNESourceSpace class declaration.
MNEProjData class declaration.
MNETriangle class declaration.
MNESurface class declaration.
Core MNE data structures (source spaces, source estimates, hemispheres).
FIFF file I/O and data structures (raw, epochs, evoked, covariance, forward).
QSharedPointer< FiffDirNode > SPtr
FIFF File I/O routines.
QSharedPointer< FiffStream > SPtr
std::unique_ptr< FiffTag > UPtr
Definition fiff_tag.h:158
Auxiliary projection data computed from MNEProjOp for efficient repeated application.
Eigen::VectorXf b
Eigen::VectorXf c
Eigen::VectorXi act
Eigen::VectorXf a
double sum_solids(const Eigen::Vector3f &from) const
void find_closest_on_surface_approx(const PointsT &r, int np, Eigen::VectorXi &nearest_tri, Eigen::VectorXf &dist, int nstep) const
Eigen::Vector3f project_to_triangle(int tri, float p, float q) const
static std::unique_ptr< MNESurface > read_bem_surface2(const QString &name, int which, bool add_geometry)
void triangle_coords(const Eigen::Vector3f &r, int tri, float &x, float &y, float &z) const
static std::unique_ptr< MNESurface > read_bem_surface(const QString &name, int which, bool add_geometry)
int project_to_surface(const MNEProjData *proj_data, const Eigen::Vector3f &r, float &distp) const
void decide_search_restriction(MNEProjData &p, int approx_best, int nstep, const Eigen::Vector3f &r) const
int nearest_triangle_point(const Eigen::Vector3f &r, const MNEProjData *user, int tri, float &x, float &y, float &z) const
void activate_neighbors(int start, Eigen::VectorXi &act, int nstep) const
std::vector< Eigen::VectorXi > neighbor_tri
std::vector< Eigen::VectorXi > neighbor_vert
int add_geometry_info(bool do_normals, bool check_too_many_neighbors)
FIFFLIB::FiffSparseMatrix dist
static double solid_angle(const Eigen::Vector3f &from, const MNELIB::MNETriangle &tri)
std::vector< MNETriangle > tris
Eigen::Matrix< float, Eigen::Dynamic, 3, Eigen::RowMajor > PointsT
int add_geometry_info2(bool do_normals)
Per-triangle geometric data for a cortical or BEM surface.
Eigen::Vector3f r13
Eigen::Vector3f nn
Eigen::Vector3f r2
Eigen::Vector3f r1
Eigen::Vector3f r3
Eigen::Vector3f r12