PolyFEM
Loading...
Searching...
No Matches
MeshUtils.cpp
Go to the documentation of this file.
1
2#include "MeshUtils.hpp"
3
10
11#include <unordered_set>
12
13#include <igl/PI.h>
14#include <igl/read_triangle_mesh.h>
15
16#include <geogram/mesh/mesh_io.h>
17#include <geogram/mesh/mesh_geometry.h>
18#include <geogram/basic/geometry.h>
19#include <geogram/mesh/mesh_preprocessing.h>
20#include <geogram/mesh/mesh_topology.h>
21#include <geogram/mesh/mesh_geometry.h>
22#include <geogram/mesh/mesh_repair.h>
23#include <geogram/mesh/mesh_AABB.h>
24#include <geogram/voronoi/CVT.h>
25#include <geogram/basic/logger.h>
27
28using namespace polyfem::io;
29using namespace polyfem::utils;
30
31bool polyfem::mesh::is_planar(const GEO::Mesh &M, const double tol)
32{
33 if (M.vertices.dimension() == 2)
34 return true;
35
36 assert(M.vertices.dimension() == 3);
37 GEO::vec3 min_corner, max_corner;
38 GEO::get_bbox(M, &min_corner[0], &max_corner[0]);
39 const double diff = (max_corner[2] - min_corner[2]);
40
41 return fabs(diff) < tol;
42}
43
44GEO::vec3 polyfem::mesh::mesh_vertex(const GEO::Mesh &M, GEO::index_t v)
45{
46 using GEO::index_t;
47 GEO::vec3 p(0, 0, 0);
48 for (index_t d = 0; d < std::min(3u, (index_t)M.vertices.dimension()); ++d)
49 {
50 if (M.vertices.double_precision())
51 {
52 p[d] = M.vertices.point_ptr(v)[d];
53 }
54 else
55 {
56 p[d] = M.vertices.single_precision_point_ptr(v)[d];
57 }
58 }
59 return p;
60}
61
62// -----------------------------------------------------------------------------
63
64GEO::vec3 polyfem::mesh::facet_barycenter(const GEO::Mesh &M, GEO::index_t f)
65{
66 using GEO::index_t;
67 GEO::vec3 p(0, 0, 0);
68 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
69 {
70 p += polyfem::mesh::mesh_vertex(M, M.facets.vertex(f, lv));
71 }
72 return p / M.facets.nb_vertices(f);
73}
74
75// -----------------------------------------------------------------------------
76
77GEO::index_t polyfem::mesh::mesh_create_vertex(GEO::Mesh &M, const GEO::vec3 &p)
78{
79 using GEO::index_t;
80 auto v = M.vertices.create_vertex();
81 for (index_t d = 0; d < std::min(3u, (index_t)M.vertices.dimension()); ++d)
82 {
83 if (M.vertices.double_precision())
84 {
85 M.vertices.point_ptr(v)[d] = p[d];
86 }
87 else
88 {
89 M.vertices.single_precision_point_ptr(v)[d] = (float)p[d];
90 }
91 }
92 return v;
93}
94
96
97void polyfem::mesh::compute_element_tags(const GEO::Mesh &M, std::vector<ElementType> &element_tags)
98{
99 using GEO::index_t;
100
101 std::vector<ElementType> old_tags = element_tags;
102
103 element_tags.resize(M.facets.nb());
104
105 // Step 0: Compute boundary vertices as true boundary + vertices incident to a polygon
106 GEO::Attribute<bool> is_boundary_vertex(M.vertices.attributes(), "boundary_vertex");
107 std::vector<bool> is_interface_vertex(M.vertices.nb(), false);
108 {
109 for (index_t f = 0; f < M.facets.nb(); ++f)
110 {
111 if (M.facets.nb_vertices(f) != 4 || (!old_tags.empty() && old_tags[f] == ElementType::INTERIOR_POLYTOPE))
112 {
113 // Vertices incident to polygonal facets (triangles or > 4 vertices) are marked as interface
114 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
115 {
116 is_interface_vertex[M.facets.vertex(f, lv)] = true;
117 }
118 }
119 }
120 }
121
122 // Step 1: Determine which vertices are regular or not
123 //
124 // Interior vertices are regular if they are incident to exactly 4 quads
125 // Boundary vertices are regular if they are incident to at most 2 quads, and no other facets
126 std::vector<int> degree(M.vertices.nb(), 0);
127 std::vector<bool> is_regular_vertex(M.vertices.nb());
128 for (index_t f = 0; f < M.facets.nb(); ++f)
129 {
130 if (M.facets.nb_vertices(f) == 4)
131 {
132 // Only count incident quads for the degree
133 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
134 {
135 index_t v = M.facets.vertex(f, lv);
136 degree[v]++;
137 }
138 }
139 }
140 for (index_t v = 0; v < M.vertices.nb(); ++v)
141 {
142 // assert(degree[v] > 0); // We assume there are no isolated vertices here
143 if (is_boundary_vertex[v] || is_interface_vertex[v])
144 {
145 is_regular_vertex[v] = (degree[v] <= 2);
146 }
147 else
148 {
149 is_regular_vertex[v] = (degree[v] == 4);
150 }
151 }
152
153 // Step 2: Iterate over the facets and determine the type
154 for (index_t f = 0; f < M.facets.nb(); ++f)
155 {
156 assert(M.facets.nb_vertices(f) > 2);
157 if (!old_tags.empty() && old_tags[f] == ElementType::INTERIOR_POLYTOPE)
158 continue;
159
160 if (M.facets.nb_vertices(f) == 4)
161 {
162 // Quad facet
163
164 // a) Determine if it is on the mesh boundary
165 bool is_boundary_facet = false;
166 bool is_interface_facet = false;
167 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
168 {
169 if (is_boundary_vertex[M.facets.vertex(f, lv)])
170 {
171 is_boundary_facet = true;
172 }
173 if (is_interface_vertex[M.facets.vertex(f, lv)])
174 {
175 is_interface_facet = true;
176 }
177 }
178
179 // b) Determine if it is regular or not
180 if (is_boundary_facet || is_interface_facet)
181 {
182 // A boundary quad is regular iff all its vertices are incident to at most 2 other quads
183 // We assume that non-boundary vertices of a boundary quads are always regular
184 bool is_singular = false;
185 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
186 {
187 index_t v = M.facets.vertex(f, lv);
188 if (is_boundary_vertex[v] || is_interface_vertex[v])
189 {
190 if (!is_regular_vertex[v])
191 {
192 is_singular = true;
193 break;
194 }
195 }
196 else
197 {
198 if (!is_regular_vertex[v])
199 {
200 element_tags[f] = ElementType::UNDEFINED;
201 break;
202 }
203 }
204 }
205
206 if (is_interface_facet)
207 {
208 element_tags[f] = ElementType::INTERFACE_CUBE;
209 }
210 else if (is_singular)
211 {
212 element_tags[f] = ElementType::SIMPLE_SINGULAR_BOUNDARY_CUBE;
213 }
214 else
215 {
216 element_tags[f] = ElementType::REGULAR_BOUNDARY_CUBE;
217 }
218 }
219 else
220 {
221 // An interior quad is regular if all its vertices are singular
222 int nb_singulars = 0;
223 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
224 {
225 if (!is_regular_vertex[M.facets.vertex(f, lv)])
226 {
227 ++nb_singulars;
228 }
229 }
230
231 if (nb_singulars == 0)
232 {
233 element_tags[f] = ElementType::REGULAR_INTERIOR_CUBE;
234 }
235 else if (nb_singulars == 1)
236 {
237 element_tags[f] = ElementType::SIMPLE_SINGULAR_INTERIOR_CUBE;
238 }
239 else
240 {
241 element_tags[f] = ElementType::MULTI_SINGULAR_INTERIOR_CUBE;
242 }
243 }
244 }
245 else
246 {
247 // Polygonal facet
248
249 // Note: In this function, we consider triangles as polygonal facets
250 ElementType tag = ElementType::INTERIOR_POLYTOPE;
251 GEO::Attribute<bool> boundary_vertices(M.vertices.attributes(), "boundary_vertex");
252 for (index_t lv = 0; lv < M.facets.nb_vertices(f); ++lv)
253 {
254 if (boundary_vertices[M.facets.vertex(f, lv)])
255 {
256 tag = ElementType::BOUNDARY_POLYTOPE;
257 // std::cout << "foo" << std::endl;
258 break;
259 }
260 }
261
262 element_tags[f] = tag;
263 }
264 }
265
266 // TODO what happens at the neighs?
267 // Override for simplices
268 for (index_t f = 0; f < M.facets.nb(); ++f)
269 {
270 if (M.facets.nb_vertices(f) == 3)
271 {
272 element_tags[f] = ElementType::SIMPLEX;
273 }
274 }
275}
276
278
279namespace
280{
281
282 // Signed area of a polygonal facet
283 double signed_area(const GEO::Mesh &M, GEO::index_t f)
284 {
285 using namespace GEO;
286 double result = 0;
287 index_t v0 = M.facet_corners.vertex(M.facets.corners_begin(f));
288 const vec3 &p0 = Geom::mesh_vertex(M, v0);
289 for (index_t c =
290 M.facets.corners_begin(f) + 1;
291 c + 1 < M.facets.corners_end(f); c++)
292 {
293 index_t v1 = M.facet_corners.vertex(c);
294 const vec3 &p1 = polyfem::mesh::mesh_vertex(M, v1);
295 index_t v2 = M.facet_corners.vertex(c + 1);
296 const vec3 &p2 = polyfem::mesh::mesh_vertex(M, v2);
297 result += Geom::triangle_signed_area(vec2(&p0[0]), vec2(&p1[0]), vec2(&p2[0]));
298 }
299 return result;
300 }
301
302} // anonymous namespace
303
305{
306 using namespace GEO;
307 vector<index_t> component;
308 index_t nb_components = get_connected_components(M, component);
309 vector<double> comp_signed_volume(nb_components, 0.0);
310 for (index_t f = 0; f < M.facets.nb(); ++f)
311 {
312 comp_signed_volume[component[f]] += signed_area(M, f);
313 }
314 for (index_t f = 0; f < M.facets.nb(); ++f)
315 {
316 if (comp_signed_volume[component[f]] < 0.0)
317 {
318 M.facets.flip(f);
319 }
320 }
321}
322
324
325void polyfem::mesh::reorder_mesh(Eigen::MatrixXd &V, Eigen::MatrixXi &F, const Eigen::VectorXi &C, Eigen::VectorXi &R)
326{
327 assert(V.rows() == C.size());
328 int num_colors = C.maxCoeff() + 1;
329 Eigen::VectorXi count(num_colors);
330 count.setZero();
331 for (int i = 0; i < C.size(); ++i)
332 {
333 ++count[C(i)];
334 }
335 R.resize(num_colors + 1);
336 R(0) = 0;
337 for (int c = 0; c < num_colors; ++c)
338 {
339 R(c + 1) = R(c) + count(c);
340 }
341 count.setZero();
342 Eigen::VectorXi remap(C.size());
343 for (int i = 0; i < C.size(); ++i)
344 {
345 remap[i] = R(C(i)) + count[C(i)];
346 ++count[C(i)];
347 }
348 // Remap vertices
349 Eigen::MatrixXd NV(V.rows(), V.cols());
350 for (int v = 0; v < V.rows(); ++v)
351 {
352 NV.row(remap(v)) = V.row(v);
353 }
354 V = NV;
355 // Remap face indices
356 for (int f = 0; f < F.rows(); ++f)
357 {
358 for (int lv = 0; lv < F.cols(); ++lv)
359 {
360 F(f, lv) = remap(F(f, lv));
361 }
362 }
363}
364
366
367namespace
368{
369
370 void compute_unsigned_distance_field(const GEO::Mesh &M,
371 const GEO::MeshFacetsAABB &aabb_tree, const Eigen::MatrixXd &P, Eigen::VectorXd &D)
372 {
373 assert(P.cols() == 3);
374 D.resize(P.rows());
375#pragma omp parallel for
376 for (int i = 0; i < P.rows(); ++i)
377 {
378 GEO::vec3 pos(P(i, 0), P(i, 1), P(i, 2));
379 double sq_dist = aabb_tree.squared_distance(pos);
380 D(i) = sq_dist;
381 }
382 }
383
384 // calculate twice signed area of triangle (0,0)-(x1,y1)-(x2,y2)
385 // return an SOS-determined sign (-1, +1, or 0 only if it's a truly degenerate triangle)
386 int orientation(
387 double x1, double y1, double x2, double y2, double &twice_signed_area)
388 {
389 twice_signed_area = y1 * x2 - x1 * y2;
390 if (twice_signed_area > 0)
391 return 1;
392 else if (twice_signed_area < 0)
393 return -1;
394 else if (y2 > y1)
395 return 1;
396 else if (y2 < y1)
397 return -1;
398 else if (x1 > x2)
399 return 1;
400 else if (x1 < x2)
401 return -1;
402 else
403 return 0; // only true when x1==x2 and y1==y2
404 }
405
406 // robust test of (x0,y0) in the triangle (x1,y1)-(x2,y2)-(x3,y3)
407 // if true is returned, the barycentric coordinates are set in a,b,c.
408 //
409 // Note: This function comes from SDFGen by Christopher Batty.
410 // https://github.com/christopherbatty/SDFGen/blob/master/makelevelset3.cpp
411 bool point_in_triangle_2d(
412 double x0, double y0, double x1, double y1,
413 double x2, double y2, double x3, double y3,
414 double &a, double &b, double &c)
415 {
416 x1 -= x0;
417 x2 -= x0;
418 x3 -= x0;
419 y1 -= y0;
420 y2 -= y0;
421 y3 -= y0;
422 int signa = orientation(x2, y2, x3, y3, a);
423 if (signa == 0)
424 return false;
425 int signb = orientation(x3, y3, x1, y1, b);
426 if (signb != signa)
427 return false;
428 int signc = orientation(x1, y1, x2, y2, c);
429 if (signc != signa)
430 return false;
431 double sum = a + b + c;
432 geo_assert(sum != 0); // if the SOS signs match and are nonzero, there's no way all of a, b, and c are zero.
433 a /= sum;
434 b /= sum;
435 c /= sum;
436 return true;
437 }
438
439 // -----------------------------------------------------------------------------
440
441 // \brief Computes the (approximate) orientation predicate in 2d.
442 // \details Computes the sign of the (approximate) signed volume of
443 // the triangle p0, p1, p2
444 // \param[in] p0 first vertex of the triangle
445 // \param[in] p1 second vertex of the triangle
446 // \param[in] p2 third vertex of the triangle
447 // \retval POSITIVE if the triangle is oriented positively
448 // \retval ZERO if the triangle is flat
449 // \retval NEGATIVE if the triangle is oriented negatively
450 // \todo check whether orientation is inverted as compared to
451 // Shewchuk's version.
452 // Taken from geogram/src/lib/geogram/delaunay/delaunay_2d.cpp
453 inline GEO::Sign orient_2d_inexact(GEO::vec2 p0, GEO::vec2 p1, GEO::vec2 p2)
454 {
455 double a11 = p1[0] - p0[0];
456 double a12 = p1[1] - p0[1];
457
458 double a21 = p2[0] - p0[0];
459 double a22 = p2[1] - p0[1];
460
461 double Delta = GEO::det2x2(
462 a11, a12,
463 a21, a22);
464
465 return GEO::geo_sgn(Delta);
466 }
467
468 // -----------------------------------------------------------------------------
469
480 template <int X = 0, int Y = 1, int Z = 2>
481 int intersect_ray_z(const GEO::Mesh &M, GEO::index_t f, const GEO::vec3 &q, double &z)
482 {
483 using namespace GEO;
484
485 index_t c = M.facets.corners_begin(f);
486 const vec3 &p1 = Geom::mesh_vertex(M, M.facet_corners.vertex(c++));
487 const vec3 &p2 = Geom::mesh_vertex(M, M.facet_corners.vertex(c++));
488 const vec3 &p3 = Geom::mesh_vertex(M, M.facet_corners.vertex(c));
489
490 double u, v, w;
491 if (point_in_triangle_2d(
492 q[X], q[Y], p1[X], p1[Y], p2[X], p2[Y], p3[X], p3[Y], u, v, w))
493 {
494 z = u * p1[Z] + v * p2[Z] + w * p3[Z];
495 auto sign = orient_2d_inexact(vec2(p1[X], p1[Y]), vec2(p2[X], p2[Y]), vec2(p3[X], p3[Y]));
496 switch (sign)
497 {
498 case GEO::POSITIVE:
499 return 1;
500 case GEO::NEGATIVE:
501 return -1;
502 case GEO::ZERO:
503 default:
504 return 0;
505 }
506 }
507
508 return 0;
509 }
510
511 // -----------------------------------------------------------------------------
512
513 void compute_sign(const GEO::Mesh &M, const GEO::MeshFacetsAABB &aabb_tree,
514 const Eigen::MatrixXd &P, Eigen::VectorXd &D)
515 {
516 assert(P.cols() == 3);
517 assert(D.size() == P.rows());
518
519 GEO::vec3 min_corner, max_corner;
520 GEO::get_bbox(M, &min_corner[0], &max_corner[0]);
521
522#pragma omp parallel for
523 for (int k = 0; k < P.rows(); ++k)
524 {
525 GEO::vec3 center(P(k, 0), P(k, 1), P(k, 2));
526
527 GEO::Box box;
528 box.xyz_min[0] = box.xyz_max[0] = center[0];
529 box.xyz_min[1] = box.xyz_max[1] = center[1];
530 box.xyz_min[2] = min_corner[2];
531 box.xyz_max[2] = max_corner[2];
532
533 std::vector<std::pair<double, int>> inter;
534 auto action = [&M, &inter, &center](GEO::index_t f) {
535 double z;
536 if (int s = intersect_ray_z(M, f, center, z))
537 {
538 inter.emplace_back(z, s);
539 }
540 };
541 aabb_tree.compute_bbox_facet_bbox_intersections(box, action);
542 std::sort(inter.begin(), inter.end());
543
544 std::vector<double> reduced;
545 for (int i = 0, s = 0; i < (int)inter.size(); ++i)
546 {
547 const int ds = inter[i].second;
548 s += ds;
549 if ((s == -1 && ds < 0) || (s == 0 && ds > 0))
550 {
551 reduced.push_back(inter[i].first);
552 }
553 }
554
555 int num_before = 0;
556 for (double z : reduced)
557 {
558 if (z < center[2])
559 {
560 ++num_before;
561 }
562 }
563 if (num_before % 2 == 1)
564 {
565 // Point is inside
566 D(k) *= -1.0;
567 }
568 }
569 }
570
571} // anonymous namespace
572
574
575void polyfem::mesh::to_geogram_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, GEO::Mesh &M)
576{
577 M.clear();
578 // Setup vertices
579 M.vertices.create_vertices((int)V.rows());
580 for (int i = 0; i < (int)M.vertices.nb(); ++i)
581 {
582 GEO::vec3 &p = M.vertices.point(i);
583 p[0] = V(i, 0);
584 p[1] = V(i, 1);
585 p[2] = V.cols() >= 3 ? V(i, 2) : 0;
586 }
587 // Setup faces
588 if (F.cols() == 3)
589 {
590 M.facets.create_triangles((int)F.rows());
591 }
592 else if (F.cols() == 4)
593 {
594 M.facets.create_quads((int)F.rows());
595 }
596 else
597 {
598 throw std::runtime_error("Mesh format not supported");
599 }
600 for (int c = 0; c < (int)M.facets.nb(); ++c)
601 {
602 for (int lv = 0; lv < F.cols(); ++lv)
603 {
604 M.facets.set_vertex(c, lv, F(c, lv));
605 }
606 }
607}
608
609// void polyfem::mesh::to_geogram_mesh_3d(const Eigen::MatrixXd &V, const Eigen::MatrixXi &C, GEO::Mesh &M) {
610// M.clear();
611// // Setup vertices
612// M.vertices.create_vertices((int) V.rows());
613// assert(V.cols() == 3);
614// for (int i = 0; i < (int) M.vertices.nb(); ++i) {
615// GEO::vec3 &p = M.vertices.point(i);
616// p[0] = V(i, 0);
617// p[1] = V(i, 1);
618// p[2] = V(i, 2);
619// }
620
621// if(C.cols() == 4)
622// M.cells.create_tets((int) C.rows());
623// else if(C.cols() == 8)
624// M.cells.create_hexes((int) C.rows());
625// else
626// assert(false);
627
628// for (int c = 0; c < (int) M.cells.nb(); ++c) {
629// for (int lv = 0; lv < C.cols(); ++lv) {
630// M.cells.set_vertex(c, lv, C(c, lv));
631// }
632// }
633// M.cells.connect();
634// // GEO::mesh_reorient(M);
635// }
636
637// -----------------------------------------------------------------------------
638
639void polyfem::mesh::from_geogram_mesh(const GEO::Mesh &M, Eigen::MatrixXd &V, Eigen::MatrixXi &F, Eigen::MatrixXi &T)
640{
641 V.resize(M.vertices.nb(), 3);
642 for (int i = 0; i < (int)M.vertices.nb(); ++i)
643 {
644 GEO::vec3 p = M.vertices.point(i);
645 V.row(i) << p[0], p[1], p[2];
646 }
647 assert(M.facets.are_simplices());
648 F.resize(M.facets.nb(), 3);
649 for (int c = 0; c < (int)M.facets.nb(); ++c)
650 {
651 for (int lv = 0; lv < 3; ++lv)
652 {
653 F(c, lv) = M.facets.vertex(c, lv);
654 }
655 }
656 assert(M.cells.are_simplices());
657 T.resize(M.cells.nb(), 4);
658 for (int c = 0; c < (int)M.cells.nb(); ++c)
659 {
660 for (int lv = 0; lv < 4; ++lv)
661 {
662 T(c, lv) = M.cells.vertex(c, lv);
663 }
664 }
665}
666
667// -----------------------------------------------------------------------------
668
669void polyfem::mesh::signed_squared_distances(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F,
670 const Eigen::MatrixXd &P, Eigen::VectorXd &D)
671{
672 GEO::Mesh M;
673 to_geogram_mesh(V, F, M);
674 GEO::MeshFacetsAABB aabb_tree(M);
675 compute_unsigned_distance_field(M, aabb_tree, P, D);
676 compute_sign(M, aabb_tree, P, D);
677}
678
679// -----------------------------------------------------------------------------
680
681double polyfem::mesh::signed_volume(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
682{
683 assert(F.cols() == 3);
684 assert(V.cols() == 3);
685 std::array<Eigen::RowVector3d, 4> t;
686 t[3] = Eigen::RowVector3d::Zero(V.cols());
687 double volume_total = 0;
688 for (int f = 0; f < F.rows(); ++f)
689 {
690 for (int lv = 0; lv < F.cols(); ++lv)
691 {
692 t[lv] = V.row(F(f, lv));
693 }
694 double vol = GEO::Geom::tetra_signed_volume(t[0].data(), t[1].data(), t[2].data(), t[3].data());
695 volume_total += vol;
696 }
697 return -volume_total;
698}
699
700// -----------------------------------------------------------------------------
701
702void polyfem::mesh::orient_closed_surface(const Eigen::MatrixXd &V, Eigen::MatrixXi &F, bool positive)
703{
704 if ((positive ? 1 : -1) * signed_volume(V, F) < 0)
705 {
706 for (int f = 0; f < F.rows(); ++f)
707 {
708 F.row(f) = F.row(f).reverse().eval();
709 }
710 }
711}
712
713// -----------------------------------------------------------------------------
714
715void polyfem::mesh::extract_polyhedra(const Mesh3D &mesh, std::vector<std::unique_ptr<GEO::Mesh>> &polys, bool triangulated)
716{
717 std::vector<int> vertex_g2l(mesh.n_vertices() + mesh.n_faces(), -1);
718 std::vector<int> vertex_l2g;
719 for (int c = 0; c < mesh.n_cells(); ++c)
720 {
721 if (!mesh.is_polytope(c))
722 {
723 continue;
724 }
725 auto poly = std::make_unique<GEO::Mesh>();
726 int nv = mesh.n_cell_vertices(c);
727 int nf = mesh.n_cell_faces(c);
728 poly->vertices.create_vertices((triangulated ? nv + nf : nv) + 1);
729 vertex_l2g.clear();
730 vertex_l2g.reserve(nv);
731 for (int lf = 0; lf < nf; ++lf)
732 {
733 GEO::vector<GEO::index_t> facet_vertices;
734 auto index = mesh.get_index_from_element(c, lf, 0);
735 for (int lv = 0; lv < mesh.n_face_vertices(index.face); ++lv)
736 {
737 Eigen::RowVector3d p = mesh.point(index.vertex);
738 if (vertex_g2l[index.vertex] < 0)
739 {
740 vertex_g2l[index.vertex] = vertex_l2g.size();
741 vertex_l2g.push_back(index.vertex);
742 }
743 int v1 = vertex_g2l[index.vertex];
744 facet_vertices.push_back(v1);
745 poly->vertices.point(v1) = GEO::vec3(p.data());
746 index = mesh.next_around_face(index);
747 }
748 if (triangulated)
749 {
750 GEO::vec3 p(0, 0, 0);
751 for (GEO::index_t lv = 0; lv < facet_vertices.size(); ++lv)
752 {
753 p += poly->vertices.point(facet_vertices[lv]);
754 }
755 p /= facet_vertices.size();
756 int v0 = vertex_l2g.size();
757 vertex_l2g.push_back(0);
758 poly->vertices.point(v0) = p;
759 for (GEO::index_t lv = 0; lv < facet_vertices.size(); ++lv)
760 {
761 int v1 = facet_vertices[lv];
762 int v2 = facet_vertices[(lv + 1) % facet_vertices.size()];
763 poly->facets.create_triangle(v0, v1, v2);
764 }
765 }
766 else
767 {
768 poly->facets.create_polygon(facet_vertices);
769 }
770 }
771 {
772 Eigen::RowVector3d p = mesh.kernel(c);
773 poly->vertices.point(nv) = GEO::vec3(p.data());
774 }
775 assert(vertex_l2g.size() == size_t(triangulated ? nv + nf : nv));
776
777 for (int v : vertex_l2g)
778 {
779 vertex_g2l[v] = -1;
780 }
781
782 poly->facets.compute_borders();
783 poly->facets.connect();
784
785 polys.emplace_back(std::move(poly));
786 }
787}
788
789// -----------------------------------------------------------------------------
790
791// In Geogram, local vertices of a hex are numbered as follows:
792//
793// v5────v7
794// ╱┆ ╱│
795// v1─┼──v3 │
796// │v4┄┄┄┼v6
797// │╱ │╱
798// v0────v2
799//
800// However, `get_ordered_vertices_from_hex()` retrieves the local vertices in
801// this order:
802//
803// v7────v6
804// ╱┆ ╱│
805// v4─┼──v5 │
806// │v3┄┄┄┼v2
807// │╱ │╱
808// v0────v1
809//
810
811void polyfem::mesh::to_geogram_mesh(const Mesh3D &mesh, GEO::Mesh &M)
812{
813 M.clear();
814 // Convert vertices
815 M.vertices.create_vertices((int)mesh.n_vertices());
816 for (int i = 0; i < (int)M.vertices.nb(); ++i)
817 {
818 auto pt = mesh.point(i);
819 GEO::vec3 &p = M.vertices.point(i);
820 p[0] = pt[0];
821 p[1] = pt[1];
822 p[2] = pt[2];
823 }
824 // Convert faces
825 for (int f = 0, lf = 0; f < mesh.n_faces(); ++f)
826 {
827 if (mesh.is_boundary_face(f))
828 {
829 int nv = mesh.n_face_vertices(f);
830 M.facets.create_polygon(nv);
831 for (int lv = 0; lv < nv; ++lv)
832 {
833 M.facets.set_vertex(lf, lv, mesh.face_vertex(f, lv));
834 }
835 ++lf;
836 }
837 }
838 // Convert cells
839 typedef std::array<int, 8> Vector8i;
840 Vector8i g2p = {{0, 4, 1, 5, 3, 7, 2, 6}};
841 for (int c = 0; c < mesh.n_cells(); ++c)
842 {
843 if (mesh.is_cube(c))
844 {
845 Vector8i lvp = mesh.get_ordered_vertices_from_hex(c);
846 Vector8i lvg;
847 for (size_t k = 0; k < 8; ++k)
848 {
849 lvg[k] = lvp[g2p[k]];
850 }
851 std::reverse(lvg.begin(), lvg.end());
852 M.cells.create_hex(
853 lvg[0], lvg[1], lvg[2], lvg[3],
854 lvg[4], lvg[5], lvg[6], lvg[7]);
855 }
856 else
857 {
858 // TODO: Support conversion of tets as well!
859 }
860 }
861 M.facets.connect();
862 M.cells.connect();
863 // M.cells.compute_borders();
864 GEO::mesh_reorient(M);
865}
866
868
869void polyfem::mesh::tertrahedralize_star_shaped_surface(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F,
870 const Eigen::RowVector3d &kernel, Eigen::MatrixXd &OV, Eigen::MatrixXi &OF, Eigen::MatrixXi &OT)
871{
872 assert(V.cols() == 3);
873 OV.resize(V.rows() + 1, V.cols());
874 OV.topRows(V.rows()) = V;
875 OV.bottomRows(1) = kernel;
876 OF = F;
877 OT.resize(OF.rows(), 4);
878 OT.col(0).setConstant(V.rows());
879 OT.rightCols(3) = F;
880}
881
883
884void polyfem::mesh::sample_surface(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, int num_samples,
885 Eigen::MatrixXd &P, Eigen::MatrixXd *N, int num_lloyd, int num_newton)
886{
887 assert(num_samples > 3);
888 GEO::Mesh M;
889 to_geogram_mesh(V, F, M);
890 GEO::CentroidalVoronoiTesselation CVT(&M);
891 // GEO::mesh_save(M, "foo.obj");
892 bool was_quiet = GEO::Logger::instance()->is_quiet();
893 GEO::Logger::instance()->set_quiet(true);
894 CVT.compute_initial_sampling(num_samples);
895 GEO::Logger::instance()->set_quiet(was_quiet);
896
897 if (num_lloyd > 0)
898 {
899 CVT.Lloyd_iterations(num_lloyd);
900 }
901
902 if (num_newton > 0)
903 {
904 CVT.Newton_iterations(num_newton);
905 }
906
907 P.resize(3, num_samples);
908 std::copy_n(CVT.embedding(0), 3 * num_samples, P.data());
909 P.transposeInPlace();
910
911 if (N)
912 {
913 GEO::MeshFacetsAABB aabb(M);
914 N->resizeLike(P);
915 for (int i = 0; i < num_samples; ++i)
916 {
917 GEO::vec3 p(P(i, 0), P(i, 1), P(i, 2));
918 GEO::vec3 nearest_point;
919 double sq_dist;
920 auto f = aabb.nearest_facet(p, nearest_point, sq_dist);
921 GEO::vec3 n = normalize(GEO::Geom::mesh_facet_normal(M, f));
922 N->row(i) << n[0], n[1], n[2];
923 }
924 }
925}
926
928
929namespace
930{
931
932 bool approx_aligned(const double *a_, const double *b_, const double *p_, const double *q_, double tol = 1e-6)
933 {
934 using namespace GEO;
935 vec3 a(a_), b(b_), p(p_), q(q_);
936 double da = std::sqrt(Geom::point_segment_squared_distance(a, p, q));
937 double db = std::sqrt(Geom::point_segment_squared_distance(b, p, q));
938 double cos_theta = Geom::cos_angle(b - a, p - q);
939 return (da < tol && db < tol && std::abs(std::abs(cos_theta) - 1.0) < tol);
940 }
941
942} // anonymous namespace
943
944// -----------------------------------------------------------------------------
945
946void polyfem::mesh::extract_parent_edges(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IE,
947 const Eigen::MatrixXd &BV, const Eigen::MatrixXi &BE, Eigen::MatrixXi &OE)
948{
949 assert(IV.cols() == 2 || IV.cols() == 3);
950 assert(BV.cols() == 2 || BV.cols() == 3);
951 typedef std::pair<int, int> Edge;
952 std::vector<Edge> selected;
953 for (int e1 = 0; e1 < IE.rows(); ++e1)
954 {
955 Eigen::RowVector3d a;
956 a.setZero();
957 a.head(IV.cols()) = IV.row(IE(e1, 0));
958 Eigen::RowVector3d b;
959 b.setZero();
960 b.head(IV.cols()) = IV.row(IE(e1, 1));
961 for (int e2 = 0; e2 < BE.rows(); ++e2)
962 {
963 Eigen::RowVector3d p;
964 p.setZero();
965 p.head(BV.cols()) = BV.row(BE(e2, 0));
966 Eigen::RowVector3d q;
967 q.setZero();
968 q.head(BV.cols()) = BV.row(BE(e2, 1));
969 if (approx_aligned(a.data(), b.data(), p.data(), q.data()))
970 {
971 selected.emplace_back(IE(e1, 0), IE(e1, 1));
972 break;
973 }
974 }
975 }
976
977 OE.resize(selected.size(), 2);
978 for (int e = 0; e < OE.rows(); ++e)
979 {
980 OE.row(e) << selected[e].first, selected[e].second;
981 }
982}
983
985
987 const Eigen::MatrixXd &vertices,
988 const Eigen::MatrixXi &codim_edges,
989 const Eigen::MatrixXi &faces,
990 Eigen::VectorXi &codim_vertices)
991{
992 std::vector<bool> is_vertex_codim(vertices.rows(), true);
993 for (int i = 0; i < codim_edges.rows(); i++)
994 {
995 for (int j = 0; j < codim_edges.cols(); j++)
996 {
997 is_vertex_codim[codim_edges(i, j)] = false;
998 }
999 }
1000 for (int i = 0; i < faces.rows(); i++)
1001 {
1002 for (int j = 0; j < faces.cols(); j++)
1003 {
1004 is_vertex_codim[faces(i, j)] = false;
1005 }
1006 }
1007 const auto n_codim_vertices = std::count(is_vertex_codim.begin(), is_vertex_codim.end(), true);
1008 codim_vertices.resize(n_codim_vertices);
1009 for (int i = 0, ci = 0; i < vertices.rows(); i++)
1010 {
1011 if (is_vertex_codim[i])
1012 {
1013 codim_vertices[ci++] = i;
1014 }
1015 }
1016}
1017
1019 const Eigen::MatrixXi &tets,
1020 Eigen::MatrixXi &faces)
1021{
1022 std::unordered_set<Eigen::Vector3i, HashMatrix> tri_to_tet(4 * tets.rows());
1023 for (int i = 0; i < tets.rows(); i++)
1024 {
1025 tri_to_tet.emplace(tets(i, 0), tets(i, 2), tets(i, 1));
1026 tri_to_tet.emplace(tets(i, 0), tets(i, 3), tets(i, 2));
1027 tri_to_tet.emplace(tets(i, 0), tets(i, 1), tets(i, 3));
1028 tri_to_tet.emplace(tets(i, 1), tets(i, 2), tets(i, 3));
1029 }
1030
1031 std::vector<Eigen::RowVector3i> faces_vector;
1032 for (const auto &tri : tri_to_tet)
1033 {
1034 // find dual triangle with reversed indices:
1035 bool is_surface_triangle =
1036 tri_to_tet.find(Eigen::Vector3i(tri[2], tri[1], tri[0])) == tri_to_tet.end()
1037 && tri_to_tet.find(Eigen::Vector3i(tri[1], tri[0], tri[2])) == tri_to_tet.end()
1038 && tri_to_tet.find(Eigen::Vector3i(tri[0], tri[2], tri[1])) == tri_to_tet.end();
1039 if (is_surface_triangle)
1040 {
1041 faces_vector.emplace_back(tri[0], tri[1], tri[2]);
1042 }
1043 }
1044
1045 faces.resize(faces_vector.size(), 3);
1046 for (int i = 0; i < faces.rows(); i++)
1047 {
1048 faces.row(i) = faces_vector[i];
1049 }
1050}
1051
1053 const Eigen::MatrixXd &vertices,
1054 const Eigen::MatrixXi &tets,
1055 Eigen::MatrixXd &surface_vertices,
1056 Eigen::MatrixXi &tris)
1057{
1058 Eigen::MatrixXi full_tris;
1059 find_triangle_surface_from_tets(tets, full_tris);
1060
1061 std::unordered_map<int, int> full_to_surface;
1062 std::vector<size_t> surface_to_full;
1063 for (int i = 0; i < full_tris.rows(); i++)
1064 {
1065 for (int j = 0; j < full_tris.cols(); j++)
1066 {
1067 if (full_to_surface.find(full_tris(i, j)) == full_to_surface.end())
1068 {
1069 full_to_surface[full_tris(i, j)] = surface_to_full.size();
1070 surface_to_full.push_back(full_tris(i, j));
1071 }
1072 }
1073 }
1074
1075 surface_vertices.resize(surface_to_full.size(), 3);
1076 for (int i = 0; i < surface_to_full.size(); i++)
1077 {
1078 surface_vertices.row(i) = vertices.row(surface_to_full[i]);
1079 }
1080
1081 tris.resize(full_tris.rows(), full_tris.cols());
1082 for (int i = 0; i < tris.rows(); i++)
1083 {
1084 for (int j = 0; j < tris.cols(); j++)
1085 {
1086 tris(i, j) = full_to_surface[full_tris(i, j)];
1087 }
1088 }
1089}
1090
1092 const std::string &mesh_path,
1093 Eigen::MatrixXd &vertices,
1094 Eigen::VectorXi &codim_vertices,
1095 Eigen::MatrixXi &codim_edges,
1096 Eigen::MatrixXi &faces)
1097{
1098 vertices.resize(0, 0);
1099 codim_vertices.resize(0);
1100 codim_edges.resize(0, 0);
1101 faces.resize(0, 0);
1102
1103 std::string lowername = mesh_path;
1104 std::transform(
1105 lowername.begin(), lowername.end(), lowername.begin(), ::tolower);
1106
1107 if (StringUtils::endswith(lowername, ".msh"))
1108 {
1109 Eigen::MatrixXi cells;
1110 std::vector<std::vector<int>> elements;
1111 std::vector<std::vector<double>> weights;
1112 std::vector<int> body_ids;
1113 if (!MshReader::load(mesh_path, vertices, cells, elements, weights, body_ids))
1114 {
1115 logger().error("Unable to load mesh: {}", mesh_path);
1116 return false;
1117 }
1118
1119 if (cells.cols() == 1)
1120 codim_vertices = cells;
1121 else if (cells.cols() == 2)
1122 codim_edges = cells;
1123 else if (cells.cols() == 3)
1124 faces = cells;
1125 else if (cells.cols() == 4)
1126 {
1127 if (vertices.cols() == 2)
1128 {
1129 logger().error("read_surface_mesh not implemented for 2D quad meshes");
1130 return false;
1131 }
1132 else
1133 {
1134 // TODO: how to distinguish between 3D tet mesh and 3D surface quad mesh?
1135 assert(vertices.cols() == 3);
1136 Eigen::MatrixXd surface_vertices;
1137 extract_triangle_surface_from_tets(vertices, cells, surface_vertices, faces);
1138 vertices = surface_vertices;
1139 }
1140 }
1141 else
1142 {
1143 logger().error("read_surface_mesh not implemented for hexahedral and polygonal/polyhedral meshes");
1144 return false;
1145 }
1146 }
1147 else if (StringUtils::endswith(lowername, ".obj")) // Use specialized OBJ reader function with polyline support
1148 {
1149 if (!OBJReader::read(mesh_path, vertices, codim_edges, faces))
1150 {
1151 logger().error("Unable to load mesh: {}", mesh_path);
1152 return false;
1153 }
1154 }
1155 else if (!igl::read_triangle_mesh(mesh_path, vertices, faces))
1156 {
1157 GEO::Mesh mesh;
1158 if (!GEO::mesh_load(mesh_path, mesh))
1159 {
1160 logger().error("Unable to load mesh: {}", mesh_path);
1161 return false;
1162 }
1163
1164 int dim = is_planar(mesh) ? 2 : 3;
1165 vertices.resize(mesh.vertices.nb(), dim);
1166 for (int vi = 0; vi < mesh.vertices.nb(); vi++)
1167 {
1168 const auto &v = mesh.vertices.point(vi);
1169 for (int vj = 0; vj < dim; vj++)
1170 {
1171 vertices(vi, vj) = v[vj];
1172 }
1173 }
1174
1175 // TODO: Check that this works even for a volumetric mesh
1176 assert(mesh.facets.nb());
1177 int face_cols = mesh.facets.nb_vertices(0);
1178 faces.resize(mesh.facets.nb(), face_cols);
1179 for (int fi = 0; fi < mesh.facets.nb(); fi++)
1180 {
1181 assert(face_cols == mesh.facets.nb_vertices(fi));
1182 for (int fj = 0; fj < mesh.facets.nb_vertices(fi); fj++)
1183 {
1184 faces(fi, fj) = mesh.facets.vertex(fi, fj);
1185 }
1186 }
1187 }
1188
1189 find_codim_vertices(vertices, codim_edges, faces, codim_vertices);
1190
1191 return true;
1192}
1193
1194int polyfem::mesh::count_faces(const int dim, const Eigen::MatrixXi &cells)
1195{
1196 std::unordered_set<std::vector<int>, HashVector> boundaries;
1197
1198 auto insert = [&](std::vector<int> v) {
1199 std::sort(v.begin(), v.end());
1200 boundaries.insert(v);
1201 };
1202
1203 for (int i = 0; i < cells.rows(); i++)
1204 {
1205 const auto &cell = cells.row(i);
1206 if (cells.cols() == 3) // triangle
1207 {
1208 insert({{cell(0), cell(1)}});
1209 insert({{cell(1), cell(2)}});
1210 insert({{cell(2), cell(0)}});
1211 }
1212 else if (cells.cols() == 4 && dim == 2) // quadralateral
1213 {
1214 insert({{cell(0), cell(1)}});
1215 insert({{cell(1), cell(2)}});
1216 insert({{cell(2), cell(3)}});
1217 insert({{cell(3), cell(0)}});
1218 }
1219 else if (cells.cols() == 4 && dim == 3) // tetrahedron
1220 {
1221 insert({{cell(0), cell(2), cell(1)}});
1222 insert({{cell(0), cell(3), cell(2)}});
1223 insert({{cell(0), cell(1), cell(3)}});
1224 insert({{cell(1), cell(2), cell(3)}});
1225 }
1226 else if (cells.cols() == 8) // hexahedron
1227 {
1228 insert({{cell(0), cell(1), cell(2), cell(3)}});
1229 insert({{cell(1), cell(5), cell(6), cell(3)}});
1230 insert({{cell(5), cell(4), cell(7), cell(6)}});
1231 insert({{cell(0), cell(4), cell(7), cell(3)}});
1232 insert({{cell(0), cell(4), cell(5), cell(1)}});
1233 insert({{cell(2), cell(6), cell(7), cell(3)}});
1234 }
1235 else
1236 {
1237 throw std::runtime_error("count_boundary_elements not implemented for polygons");
1238 }
1239 }
1240
1241 return boundaries.size();
1242}
1243
1245{
1246 using namespace GEO;
1247 typedef std::pair<index_t, index_t> Edge;
1248
1249 if (M.edges.nb() > 0)
1250 return;
1251
1252 M.facets.connect();
1253 M.cells.connect();
1254 if (M.cells.nb() != 0 && M.facets.nb() == 0)
1255 {
1256 M.cells.compute_borders();
1257 }
1258
1259 // Compute a list of all the edges, and store edge index as a corner attribute
1260 std::vector<std::pair<Edge, index_t>> e2c; // edge to corner id
1261 for (index_t f = 0; f < M.facets.nb(); ++f)
1262 {
1263 for (index_t c = M.facets.corners_begin(f); c < M.facets.corners_end(f); ++c)
1264 {
1265 index_t v = M.facet_corners.vertex(c);
1266 index_t c2 = M.facets.next_corner_around_facet(f, c);
1267 index_t v2 = M.facet_corners.vertex(c2);
1268 e2c.emplace_back(std::make_pair(std::min(v, v2), std::max(v, v2)), c);
1269 }
1270 }
1271 std::sort(e2c.begin(), e2c.end());
1272
1273 // Assign unique id to edges
1274 M.edges.clear();
1275 Edge prev_e(-1, -1);
1276 for (const auto &kv : e2c)
1277 {
1278 Edge e = kv.first;
1279 if (e != prev_e)
1280 {
1281 M.edges.create_edge(e.first, e.second);
1282 prev_e = e;
1283 }
1284 }
1285}
int V
QuadratureVector da
Definition Assembler.cpp:23
void find_codim_vertices(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &codim_edges, const Eigen::MatrixXi &faces, Eigen::VectorXi &codim_vertices)
void find_triangle_surface_from_tets(const Eigen::MatrixXi &tets, Eigen::MatrixXi &faces)
std::vector< Eigen::VectorXi > faces
int z
static bool load(const std::string &path, Eigen::MatrixXd &vertices, Eigen::MatrixXi &cells, std::vector< std::vector< int > > &elements, std::vector< std::vector< double > > &weights, std::vector< int > &body_ids)
Definition MshReader.cpp:29
static bool read(const std::string obj_file_name, std::vector< std::vector< double > > &V, std::vector< std::vector< double > > &TC, std::vector< std::vector< double > > &N, std::vector< std::vector< int > > &F, std::vector< std::vector< int > > &FTC, std::vector< std::vector< int > > &FN, std::vector< std::vector< int > > &L)
Read a mesh from an ascii obj file.
Definition OBJReader.cpp:32
virtual Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const =0
virtual RowVectorNd kernel(const int cell_id) const =0
std::array< int, 8 > get_ordered_vertices_from_hex(const int element_index) const
Definition Mesh3D.cpp:322
virtual int n_cell_faces(const int c_id) const =0
virtual Navigation3D::Index next_around_face(Navigation3D::Index idx) const =0
virtual int n_vertices() const =0
number of vertices
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:365
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:352
virtual RowVectorNd point(const int global_index) const =0
point coordinates
virtual bool is_boundary_face(const int face_global_id) const =0
is face boundary
virtual int n_cells() const =0
number of cells
virtual int n_faces() const =0
number of faces
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
virtual int n_cell_vertices(const int c_id) const =0
number of vertices of a cell
virtual int face_vertex(const int f_id, const int lv_id) const =0
id of the face vertex
M
Definition eigs.py:94
Eigen::ArrayXd P(const int m, const int p, const Eigen::ArrayXd &z)
Definition p_n_bases.cpp:42
bool is_planar(const GEO::Mesh &M, const double tol=1e-5)
Determine if the given mesh is planar (2D or tiny z-range).
Definition MeshUtils.cpp:31
void sample_surface(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, int num_samples, Eigen::MatrixXd &P, Eigen::MatrixXd *N=nullptr, int num_lloyd=10, int num_newton=10)
Samples points on a surface.
void orient_closed_surface(const Eigen::MatrixXd &V, Eigen::MatrixXi &F, bool positive=true)
Orient a triangulated surface to have positive volume.
double signed_volume(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
Compute the signed volume of a surface mesh.
void reorder_mesh(Eigen::MatrixXd &V, Eigen::MatrixXi &F, const Eigen::VectorXi &C, Eigen::VectorXi &R)
Reorder vertices of a mesh using color tags, so that vertices are ordered by increasing colors.
void orient_normals_2d(GEO::Mesh &M)
Orient facets of a 2D mesh so that each connected component has positive volume.
void extract_parent_edges(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IE, const Eigen::MatrixXd &BV, const Eigen::MatrixXi &BE, Eigen::MatrixXi &OE)
Extract a set of edges that are overlap with a set given set of parent edges, using vertices position...
ElementType
Type of Element, check [Poly-Spline Finite Element Method] for a complete description.
Definition Mesh.hpp:23
void extract_triangle_surface_from_tets(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &tets, Eigen::MatrixXd &surface_vertices, Eigen::MatrixXi &tris)
Extract triangular surface from a tetmesh.
GEO::vec3 mesh_vertex(const GEO::Mesh &M, GEO::index_t v)
Retrieve a 3D vector with the position of a given vertex.
Definition MeshUtils.cpp:44
GEO::index_t mesh_create_vertex(GEO::Mesh &M, const GEO::vec3 &p)
Definition MeshUtils.cpp:77
void to_geogram_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, GEO::Mesh &M)
Converts a triangle mesh to a Geogram mesh.
void signed_squared_distances(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXd &P, Eigen::VectorXd &D)
Computes the signed squared distance from a list of points to a triangle mesh.
GEO::vec3 facet_barycenter(const GEO::Mesh &M, GEO::index_t f)
Definition MeshUtils.cpp:64
void from_geogram_mesh(const GEO::Mesh &M, Eigen::MatrixXd &V, Eigen::MatrixXi &F, Eigen::MatrixXi &T)
Extract simplices from a Geogram mesh.
void tertrahedralize_star_shaped_surface(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::RowVector3d &kernel, Eigen::MatrixXd &OV, Eigen::MatrixXi &OF, Eigen::MatrixXi &OT)
Tetrahedralize a star-shaped mesh, with a given point in its kernel.
void generate_edges(GEO::Mesh &M)
assing edges to M
bool read_surface_mesh(const std::string &mesh_path, Eigen::MatrixXd &vertices, Eigen::VectorXi &codim_vertices, Eigen::MatrixXi &codim_edges, Eigen::MatrixXi &faces)
read a surface mesh
void extract_polyhedra(const Mesh3D &mesh, std::vector< std::unique_ptr< GEO::Mesh > > &polys, bool triangulated=false)
Extract polyhedra from a 3D volumetric mesh.
int count_faces(const int dim, const Eigen::MatrixXi &cells)
Count the number of boundary elements (triangles for tetmesh and edges for triangle mesh)
void compute_element_tags(const GEO::Mesh &M, std::vector< ElementType > &element_tags)
Compute the type of each facet in a surface mesh.
Definition MeshUtils.cpp:97
bool endswith(const std::string &str, const std::string &suffix)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42