PolyFEM
Loading...
Searching...
No Matches
CMesh2D.cpp
Go to the documentation of this file.
3
9
10#include <geogram/basic/file_system.h>
11#include <geogram/mesh/mesh_io.h>
12#include <geogram/mesh/mesh_geometry.h>
13#include <geogram/mesh/mesh_repair.h>
14
15#include <cassert>
16#include <array>
17
18namespace polyfem
19{
20 using namespace io;
21 using namespace utils;
22
23 namespace mesh
24 {
25 void CMesh2D::refine(const int n_refinement, const double t)
26 {
27 // return;
28 if (n_refinement <= 0)
29 {
30 return;
31 }
32
33 orders_.resize(0, 0);
34
35 bool all_simplicial = true;
36 for (int e = 0; e < n_elements(); ++e)
37 {
38 all_simplicial &= is_simplex(e);
39 }
40
41 for (int i = 0; i < n_refinement; ++i)
42 {
43 GEO::Mesh mesh;
44 mesh.copy(mesh_);
45
46 c2e_.reset();
47 boundary_vertices_.reset();
48 boundary_edges_.reset();
49
50 mesh_.clear(false, false);
51
52 // TODO add tags to the refinement
53 if (all_simplicial)
54 {
56 }
57 else if (t <= 0)
58 {
60 }
61 else
62 {
64 }
65
67 c2e_ = std::make_unique<GEO::Attribute<GEO::index_t>>(mesh_.facet_corners.attributes(), "edge_id");
68 boundary_vertices_ = std::make_unique<GEO::Attribute<bool>>(mesh_.vertices.attributes(), "boundary_vertex");
69 boundary_edges_ = std::make_unique<GEO::Attribute<bool>>(mesh_.edges.attributes(), "boundary_edge");
70 }
71
74
75 in_ordered_vertices_ = Eigen::VectorXi::LinSpaced(mesh_.vertices.nb(), 0, mesh_.vertices.nb() - 1);
76 assert(in_ordered_vertices_[0] == 0);
77 assert(in_ordered_vertices_[1] == 1);
78 assert(in_ordered_vertices_[2] == 2);
79 assert(in_ordered_vertices_[in_ordered_vertices_.size() - 1] == mesh_.vertices.nb() - 1);
80
81 in_ordered_edges_.resize(mesh_.edges.nb(), 2);
82
83 for (int e = 0; e < (int)mesh_.edges.nb(); ++e)
84 {
85 for (int lv = 0; lv < 2; ++lv)
86 {
87 in_ordered_edges_(e, lv) = mesh_.edges.vertex(e, lv);
88 }
89 assert(in_ordered_edges_(e, 0) != in_ordered_edges_(e, 1));
90 }
91 assert(in_ordered_edges_.size() > 0);
92
93 in_ordered_faces_.resize(0, 0);
94 }
95
96 bool CMesh2D::load(const std::string &path)
97 {
98 // This method should be used for special loading, like hybrid in 3d
99
100 // edge_nodes_.clear();
101 // face_nodes_.clear();
102 // cell_nodes_.clear();
103 // order_ = 1;
104
105 // c2e_.reset();
106 // boundary_vertices_.reset();
107 // boundary_edges_.reset();
108
109 // mesh_.clear(false,false);
110
111 // if (!StringUtils::endswith(path, "msh"))
112 // {
113 // Eigen::MatrixXd vertices;
114 // Eigen::MatrixXi cells;
115 // std::vector<std::vector<int>> elements;
116 // std::vector<std::vector<double>> weights;
117
118 // if(!MshReader::load(path, vertices, cells, elements, weights))
119 // return false;
120
121 // build_from_matrices(vertices, cells);
122 // attach_higher_order_nodes(vertices, elements);
123 // cell_weights_ = weights;
124 // }
125 // else
126 // {
127 // if(!mesh_load(path, mesh_))
128 // return false;
129 // }
130
131 // orient_normals_2d(mesh_);
132 // Navigation::prepare_mesh(mesh_);
133 // c2e_ = std::make_unique<GEO::Attribute<GEO::index_t>>(mesh_.facet_corners.attributes(), "edge_id");
134 // boundary_vertices_ = std::make_unique<GEO::Attribute<bool>>(mesh_.vertices.attributes(), "boundary_vertex");
135 // boundary_edges_ = std::make_unique<GEO::Attribute<bool>>(mesh_.edges.attributes(), "boundary_edge");
136
137 // compute_elements_tag();
138 assert(false);
139 return false;
140 }
141
142 bool CMesh2D::load(const GEO::Mesh &mesh)
143 {
144 edge_nodes_.clear();
145 face_nodes_.clear();
146 cell_nodes_.clear();
147
148 c2e_.reset();
149 boundary_vertices_.reset();
150 boundary_edges_.reset();
151
152 mesh_.clear(false, false);
153 mesh_.copy(mesh);
154
157 c2e_ = std::make_unique<GEO::Attribute<GEO::index_t>>(mesh_.facet_corners.attributes(), "edge_id");
158 boundary_vertices_ = std::make_unique<GEO::Attribute<bool>>(mesh_.vertices.attributes(), "boundary_vertex");
159 boundary_edges_ = std::make_unique<GEO::Attribute<bool>>(mesh_.edges.attributes(), "boundary_edge");
160
162 return true;
163 }
164
165 bool CMesh2D::save(const std::string &path) const
166 {
167 if (!mesh_save(mesh_, path))
168 return false;
169
170 return true;
171 }
172
173 bool CMesh2D::build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
174 {
175 edge_nodes_.clear();
176 face_nodes_.clear();
177 cell_nodes_.clear();
178
179 c2e_.reset();
180 boundary_vertices_.reset();
181 boundary_edges_.reset();
182
183 mesh_.clear(false, false);
185
188
189 c2e_ = std::make_unique<GEO::Attribute<GEO::index_t>>(mesh_.facet_corners.attributes(), "edge_id");
190 boundary_vertices_ = std::make_unique<GEO::Attribute<bool>>(mesh_.vertices.attributes(), "boundary_vertex");
191 boundary_edges_ = std::make_unique<GEO::Attribute<bool>>(mesh_.edges.attributes(), "boundary_edge");
192
194 return true;
195 }
196
197 void CMesh2D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
198 {
199 edge_nodes_.clear();
200 face_nodes_.clear();
201 cell_nodes_.clear();
202
203 edge_nodes_.resize(n_edges());
204 face_nodes_.resize(n_faces());
205
206 orders_.resize(n_faces(), 1);
207
208 assert(nodes.size() == n_faces());
209
210 for (int f = 0; f < n_faces(); ++f)
211 {
212 auto index = get_index_from_face(f);
213
214 const auto &nodes_ids = nodes[f];
215
216 if (nodes_ids.size() == 3)
217 {
218 orders_(f) = 1;
219 continue;
220 }
221 // P2
222 else if (nodes_ids.size() == 6)
223 {
224 orders_(f) = 2;
225
226 for (int le = 0; le < 3; ++le)
227 {
228 auto &n = edge_nodes_[index.edge];
229
230 // nodes not aleardy created
231 if (n.nodes.size() <= 0)
232 {
233 n.v1 = index.vertex;
234 n.v2 = switch_vertex(index).vertex;
235
236 int node_index = 0;
237 if ((n.v1 == nodes_ids[0] && n.v2 == nodes_ids[1]) || (n.v2 == nodes_ids[0] && n.v1 == nodes_ids[1]))
238 node_index = 3;
239 else if ((n.v1 == nodes_ids[1] && n.v2 == nodes_ids[2]) || (n.v2 == nodes_ids[1] && n.v1 == nodes_ids[2]))
240 node_index = 4;
241 else
242 node_index = 5;
243
244 n.nodes.resize(1, 2);
245 n.nodes << V(nodes_ids[node_index], 0), V(nodes_ids[node_index], 1);
246 n.nodes_ids.push_back(nodes_ids[node_index]);
247 }
248 index = next_around_face(index);
249 }
250 }
251 // P3
252 else if (nodes_ids.size() == 10)
253 {
254 orders_(f) = 3;
255
256 for (int le = 0; le < 3; ++le)
257 {
258 auto &n = edge_nodes_[index.edge];
259
260 // nodes not aleardy created
261 if (n.nodes.size() <= 0)
262 {
263 n.v1 = index.vertex;
264 n.v2 = switch_vertex(index).vertex;
265
266 int node_index1 = 0;
267 int node_index2 = 0;
268 if (n.v1 == nodes_ids[0] && n.v2 == nodes_ids[1])
269 {
270 node_index1 = 3;
271 node_index2 = 4;
272 }
273 else if (n.v2 == nodes_ids[0] && n.v1 == nodes_ids[1])
274 {
275 node_index1 = 4;
276 node_index2 = 3;
277 }
278 else if (n.v1 == nodes_ids[1] && n.v2 == nodes_ids[2])
279 {
280 node_index1 = 5;
281 node_index2 = 6;
282 }
283 else if (n.v2 == nodes_ids[1] && n.v1 == nodes_ids[2])
284 {
285 node_index1 = 6;
286 node_index2 = 5;
287 }
288 else if (n.v1 == nodes_ids[2] && n.v2 == nodes_ids[0])
289 {
290 node_index1 = 7;
291 node_index2 = 8;
292 }
293 else
294 {
295 assert(n.v2 == nodes_ids[2] && n.v1 == nodes_ids[0]);
296 node_index1 = 8;
297 node_index2 = 7;
298 }
299
300 n.nodes.resize(2, 2);
301 n.nodes.row(0) << V(nodes_ids[node_index1], 0), V(nodes_ids[node_index1], 1);
302 n.nodes.row(1) << V(nodes_ids[node_index2], 0), V(nodes_ids[node_index2], 1);
303
304 n.nodes_ids.push_back(nodes_ids[node_index1]);
305 n.nodes_ids.push_back(nodes_ids[node_index2]);
306 }
307 index = next_around_face(index);
308 }
309
310 {
311 auto &n = face_nodes_[f];
312 n.v1 = mesh_.facets.vertex(f, 0);
313 n.v2 = mesh_.facets.vertex(f, 1);
314 n.v3 = mesh_.facets.vertex(f, 2);
315 n.nodes.resize(1, 2);
316 n.nodes << V(nodes_ids[9], 0), V(nodes_ids[9], 1);
317 n.nodes_ids.push_back(nodes_ids[9]);
318 }
319 }
320 // P4
321 else if (nodes_ids.size() == 15)
322 {
323 orders_(f) = 4;
324 assert(false);
325 // unsupported P4 for geometry, need meshes for testing
326 }
327 // unsupported
328 else
329 {
330 assert(false);
331 }
332 }
333
334 if (orders_.maxCoeff() == 1)
335 orders_.resize(0, 0);
336 }
337
338 std::pair<RowVectorNd, int> CMesh2D::edge_node(const Navigation::Index &index, const int n_new_nodes, const int i) const
339 {
340 if (orders_.size() <= 0 || orders_(index.face) == 1 || edge_nodes_.empty() || edge_nodes_[index.edge].nodes.rows() != n_new_nodes)
341 {
342 const auto v1 = point(index.vertex);
343 const auto v2 = point(switch_vertex(index).vertex);
344
345 const double t = i / (n_new_nodes + 1.0);
346
347 return std::make_pair((1 - t) * v1 + t * v2, -1);
348 }
349
350 const auto &n = edge_nodes_[index.edge];
351 if (n.v1 == index.vertex)
352 return std::make_pair(n.nodes.row(i - 1), n.nodes_ids[i - 1]);
353 else
354 {
355 assert(n.v2 == index.vertex);
356 return std::make_pair(n.nodes.row(n.nodes.rows() - i), n.nodes_ids[n.nodes_ids.size() - i]);
357 }
358 }
359
360 std::pair<RowVectorNd, int> CMesh2D::face_node(const Navigation::Index &index, const int n_new_nodes, const int i, const int j) const
361 {
362 if (is_simplex(index.face))
363 {
364 if (orders_.size() <= 0 || orders_(index.face) == 1 || orders_(index.face) == 2 || face_nodes_.empty() || face_nodes_[index.face].nodes.rows() != n_new_nodes)
365 {
366 const auto v1 = point(index.vertex);
367 const auto v2 = point(switch_vertex(index).vertex);
368 const auto v3 = point(switch_vertex(switch_edge(index)).vertex);
369
370 const double b2 = i / (n_new_nodes + 2.0);
371 const double b3 = j / (n_new_nodes + 2.0);
372 const double b1 = 1 - b3 - b2;
373 assert(b3 < 1);
374 assert(b3 > 0);
375
376 return std::make_pair(b1 * v1 + b2 * v2 + b3 * v3, -1);
377 }
378
379 assert(orders_(index.face) == 3);
380 // unsupported P4 for geometry
381 const auto &n = face_nodes_[index.face];
382 return std::make_pair(n.nodes.row(0), n.nodes_ids[0]);
383 }
384 else if (is_cube(index.face))
385 {
386 // supports only blilinear quads
387 assert(orders_.size() <= 0 || orders_(index.face) == 1);
388
389 const auto v1 = point(index.vertex);
390 const auto v2 = point(switch_vertex(index).vertex);
391 const auto v3 = point(switch_vertex(switch_edge(switch_vertex(index))).vertex);
392 const auto v4 = point(switch_vertex(switch_edge(index)).vertex);
393
394 const double b1 = i / (n_new_nodes + 1.0);
395 const double b2 = j / (n_new_nodes + 1.0);
396
397 return std::make_pair(v1 * (1 - b1) * (1 - b2) + v2 * b1 * (1 - b2) + v3 * b1 * b2 + v4 * (1 - b1) * b2, -1);
398 }
399
400 assert(false);
401 return std::make_pair(RowVectorNd(2, 1), -1);
402 }
403
405 {
406 GEO::vec3 min_corner, max_corner;
407 GEO::get_bbox(mesh_, &min_corner[0], &max_corner[0]);
408 min.resize(2);
409 max.resize(2);
410
411 min(0) = min_corner.x;
412 min(1) = min_corner.y;
413
414 max(0) = max_corner.x;
415 max(1) = max_corner.y;
416 }
417
419 {
420
421 GEO::vec3 min_corner, max_corner;
422 GEO::get_bbox(mesh_, &min_corner[0], &max_corner[0]);
423 GEO::vec3 extent = max_corner - min_corner;
424 double scaling = std::max(extent[0], std::max(extent[1], extent[2]));
425 // const GEO::vec3 origin = 0.5 * (min_corner + max_corner);
426 const GEO::vec3 origin = min_corner;
427 for (GEO::index_t v = 0; v < mesh_.vertices.nb(); ++v)
428 {
429 mesh_.vertices.point(v) = (mesh_.vertices.point(v) - origin) / scaling;
430 }
431 Eigen::RowVector2d shift;
432 shift << origin[0], origin[1];
433 for (auto &n : edge_nodes_)
434 {
435 if (n.nodes.size() > 0)
436 n.nodes = (n.nodes.rowwise() - shift) / scaling;
437 }
438 for (auto &n : face_nodes_)
439 {
440 if (n.nodes.size() > 0)
441 n.nodes = (n.nodes.rowwise() - shift) / scaling;
442 }
443
444 logger().debug("-- bbox before normalization:");
445 logger().debug(" min : {} {}", min_corner[0], min_corner[1]);
446 logger().debug(" max : {} {}", max_corner[0], max_corner[1]);
447 logger().debug(" extent: {} {}", max_corner[0] - min_corner[0], max_corner[1] - min_corner[1]);
448 GEO::get_bbox(mesh_, &min_corner[0], &max_corner[0]);
449 logger().debug("-- bbox after normalization:");
450 logger().debug(" min : {} {}", min_corner[0], min_corner[1]);
451 logger().debug(" max : {} {}", max_corner[0], max_corner[1]);
452 logger().debug(" extent: {} {}", max_corner[0] - min_corner[0], max_corner[1] - min_corner[1]);
453
454 Eigen::MatrixXd p0, p1, p;
455 get_edges(p0, p1);
456 p = p0 - p1;
457 logger().debug("-- edge length after normalization:");
458 logger().debug(" min: {}", p.rowwise().norm().minCoeff());
459 logger().debug(" max: {}", p.rowwise().norm().maxCoeff());
460 logger().debug(" avg: {}", p.rowwise().norm().mean());
461 }
462
463 double CMesh2D::edge_length(const int gid) const
464 {
465 const int v0 = mesh_.edges.vertex(gid, 0);
466 const int v1 = mesh_.edges.vertex(gid, 1);
467
468 return (point(v0) - point(v1)).norm();
469 }
470
471 void CMesh2D::set_point(const int global_index, const RowVectorNd &p)
472 {
473 mesh_.vertices.point(global_index).x = p(0);
474 mesh_.vertices.point(global_index).y = p(1);
475 }
476
477 RowVectorNd CMesh2D::point(const int global_index) const
478 {
479 const double *ptr = mesh_.vertices.point_ptr(global_index);
480 RowVectorNd p(2);
481 p(0) = ptr[0];
482 p(1) = ptr[1];
483 return p;
484 }
485
486 bool CMesh2D::is_boundary_element(const int element_global_id) const
487 {
488 auto index = get_index_from_face(element_global_id);
489
490 for (int i = 0; i < n_face_vertices(element_global_id); ++i)
491 {
492 if (is_boundary_edge(index.edge))
493 return true;
494
495 index = next_around_face(index);
496 }
497
498 return false;
499 }
500
501 // void CMesh2D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
502 // {
503 // ranges.clear();
504
505 // std::vector<Eigen::MatrixXi> local_tris(mesh_.facets.nb());
506 // std::vector<Eigen::MatrixXd> local_pts(mesh_.facets.nb());
507
508 // int total_tris = 0;
509 // int total_pts = 0;
510
511 // ranges.push_back(0);
512
513 // for (GEO::index_t f = 0; f < mesh_.facets.nb(); ++f)
514 // {
515 // const int n_vertices = mesh_.facets.nb_vertices(f);
516
517 // Eigen::MatrixXd face_pts(n_vertices, 2);
518 // // Eigen::MatrixXi edges(n_vertices,2);
519 // local_tris[f].resize(n_vertices - 2, 3);
520
521 // for (int i = 0; i < n_vertices; ++i)
522 // {
523 // const int vertex = mesh_.facets.vertex(f, i);
524 // const double *pt = mesh_.vertices.point_ptr(vertex);
525 // face_pts(i, 0) = pt[0];
526 // face_pts(i, 1) = pt[1];
527
528 // // edges(i, 0) = i;
529 // // edges(i, 1) = (i+1) % n_vertices;
530 // }
531
532 // for (int i = 1; i < n_vertices - 1; ++i)
533 // {
534 // local_tris[f].row(i - 1) << 0, i, i + 1;
535 // }
536
537 // local_pts[f] = face_pts;
538
539 // total_tris += local_tris[f].rows();
540 // total_pts += local_pts[f].rows();
541
542 // ranges.push_back(total_tris);
543
544 // assert(local_pts[f].rows() == face_pts.rows());
545 // }
546
547 // tris.resize(total_tris, 3);
548 // pts.resize(total_pts, 2);
549
550 // int tri_index = 0;
551 // int pts_index = 0;
552 // for (std::size_t i = 0; i < local_tris.size(); ++i)
553 // {
554 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
555 // tri_index += local_tris[i].rows();
556
557 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
558 // pts_index += local_pts[i].rows();
559 // }
560 // }
561
567
572
574 {
575 const int v0 = mesh_.edges.vertex(index, 0);
576 const int v1 = mesh_.edges.vertex(index, 1);
577
578 return 0.5 * (point(v0) + point(v1));
579 }
580
581 void CMesh2D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
582 {
583 body_ids_.resize(n_elements());
584 std::fill(body_ids_.begin(), body_ids_.end(), -1);
585
586 for (int e = 0; e < n_elements(); ++e)
587 {
588 const auto bary = face_barycenter(e);
589 body_ids_[e] = marker(e, element_vertices(e), bary);
590 }
591 }
592
593 void CMesh2D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
594 {
595 boundary_ids_.resize(n_edges());
596
597 for (int e = 0; e < n_edges(); ++e)
598 {
599 bool is_boundary = is_boundary_edge(e);
600 const auto p = edge_barycenter(e);
601 std::vector<int> vs = {edge_vertex(e, 0), edge_vertex(e, 1)};
602 std::sort(vs.begin(), vs.end());
603 boundary_ids_[e] = marker(e, vs, p, is_boundary);
604 }
605 }
606
607 void CMesh2D::append(const Mesh &mesh)
608 {
609 assert(typeid(mesh) == typeid(CMesh2D));
610 Mesh::append(mesh);
611
612 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
613
614 const int n_v = n_vertices();
615 const int n_f = n_faces();
616
617 mesh_.vertices.create_vertices(mesh2d.n_vertices());
618 for (int i = n_v; i < (int)mesh_.vertices.nb(); ++i)
619 {
620 GEO::vec3 &p = mesh_.vertices.point(i);
621 set_point(i, mesh2d.point(i - n_v));
622 }
623
624 std::vector<GEO::index_t> indices;
625 for (int i = 0; i < mesh2d.n_faces(); ++i)
626 {
627 indices.clear();
628 for (int j = 0; j < mesh2d.mesh_.facets.nb_vertices(i); ++j)
629 indices.push_back(mesh2d.mesh_.facets.vertex(i, j) + n_v);
630
631 mesh_.facets.create_polygon(indices.size(), &indices[0]);
632 }
633
634 assert(n_vertices() == n_v + mesh2d.n_vertices());
635 assert(n_faces() == n_f + mesh2d.n_faces());
636
637 c2e_.reset();
638 boundary_vertices_.reset();
639 boundary_edges_.reset();
641 c2e_ = std::make_unique<GEO::Attribute<GEO::index_t>>(mesh_.facet_corners.attributes(), "edge_id");
642 boundary_vertices_ = std::make_unique<GEO::Attribute<bool>>(mesh_.vertices.attributes(), "boundary_vertex");
643 boundary_edges_ = std::make_unique<GEO::Attribute<bool>>(mesh_.edges.attributes(), "boundary_edge");
644 }
645
646 std::unique_ptr<Mesh> CMesh2D::copy() const
647 {
648 std::unique_ptr<CMesh2D> copy_mesh = std::make_unique<CMesh2D>();
649 copy_mesh->load(this->mesh_);
650
651 // Manually copy parent's data
652 copy_mesh->elements_tag_ = this->elements_tag_;
653 copy_mesh->node_ids_ = this->node_ids_;
654 copy_mesh->boundary_ids_ = this->boundary_ids_;
655 copy_mesh->body_ids_ = this->body_ids_;
656 copy_mesh->orders_ = this->orders_;
657 copy_mesh->is_rational_ = this->is_rational_;
658 copy_mesh->edge_nodes_ = this->edge_nodes_;
659 copy_mesh->face_nodes_ = this->face_nodes_;
660 copy_mesh->cell_nodes_ = this->cell_nodes_;
661 copy_mesh->cell_weights_ = this->cell_weights_;
662 copy_mesh->in_ordered_vertices_ = this->in_ordered_vertices_;
663 copy_mesh->in_ordered_edges_ = this->in_ordered_edges_;
664 copy_mesh->in_ordered_faces_ = this->in_ordered_faces_;
665
666 return copy_mesh;
667 }
668 } // namespace mesh
669} // namespace polyfem
int V
virtual RowVectorNd edge_barycenter(const int index) const override
edge barycenter
Definition CMesh2D.cpp:573
Navigation::Index switch_edge(Navigation::Index idx) const override
Definition CMesh2D.hpp:87
bool load(const std::string &path) override
loads a mesh from the path
Definition CMesh2D.cpp:96
void normalize() override
normalize the mesh
Definition CMesh2D.cpp:418
void compute_elements_tag() override
compute element types, see ElementType
Definition CMesh2D.cpp:562
int n_faces() const override
number of faces
Definition CMesh2D.hpp:33
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition CMesh2D.hpp:37
int n_vertices() const override
number of vertices
Definition CMesh2D.hpp:35
int n_edges() const override
number of edges
Definition CMesh2D.hpp:34
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition CMesh2D.cpp:486
double edge_length(const int gid) const override
edge length
Definition CMesh2D.cpp:463
std::unique_ptr< GEO::Attribute< bool > > boundary_edges_
Definition CMesh2D.hpp:104
std::pair< RowVectorNd, int > edge_node(const Navigation::Index &index, const int n_new_nodes, const int i) const override
Definition CMesh2D.cpp:338
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
Definition CMesh2D.cpp:646
Navigation::Index switch_vertex(Navigation::Index idx) const override
Definition CMesh2D.hpp:86
void compute_boundary_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &, bool)> &marker) override
computes boundary selections based on a function
Definition CMesh2D.cpp:593
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition CMesh2D.cpp:471
bool is_boundary_edge(const int edge_global_id) const override
is edge boundary
Definition CMesh2D.hpp:48
std::unique_ptr< GEO::Attribute< GEO::index_t > > c2e_
Definition CMesh2D.hpp:102
std::unique_ptr< GEO::Attribute< bool > > boundary_vertices_
Definition CMesh2D.hpp:103
std::pair< RowVectorNd, int > face_node(const Navigation::Index &index, const int n_new_nodes, const int i, const int j) const override
Definition CMesh2D.cpp:360
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition CMesh2D.cpp:607
int edge_vertex(const int e_id, const int lv_id) const override
id of the edge vertex
Definition CMesh2D.hpp:40
Navigation::Index get_index_from_face(int f, int lv=0) const override
Definition CMesh2D.hpp:83
virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition CMesh2D.cpp:404
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition CMesh2D.cpp:197
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
Definition CMesh2D.cpp:173
void refine(const int n_refinement, const double t) override
refine the mesh
Definition CMesh2D.cpp:25
bool save(const std::string &path) const override
Definition CMesh2D.cpp:165
void compute_body_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &)> &marker) override
computes boundary selections based on a function
Definition CMesh2D.cpp:581
virtual RowVectorNd point(const int global_index) const override
point coordinates
Definition CMesh2D.cpp:477
virtual void update_elements_tag() override
Update elements types.
Definition CMesh2D.cpp:568
RowVectorNd face_barycenter(const int index) const override
face barycenter
Definition Mesh2D.cpp:69
Navigation::Index next_around_face(Navigation::Index idx) const
Definition Mesh2D.hpp:61
void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const override
Get all the edges.
Definition Mesh2D.cpp:21
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:161
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:664
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:656
bool is_rational_
stores if the mesh is rational
Definition Mesh.hpp:666
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:352
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:682
bool is_simplex(const int el_id) const
checks if element is simples compatible
Definition Mesh.cpp:422
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:660
std::vector< int > node_ids_
list of node labels
Definition Mesh.hpp:658
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:673
std::vector< std::vector< double > > cell_weights_
weights associates to cells for rational polynomail meshes
Definition Mesh.hpp:675
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
Definition Mesh.hpp:229
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:662
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:671
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:669
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:680
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:494
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:678
void prepare_mesh(GEO::Mesh &M)
SplitFunction catmul_clark_split_func()
SplitFunction polar_split_func(double t)
Helper function.
void refine_triangle_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out)
Refine a triangle mesh.
void orient_normals_2d(GEO::Mesh &M)
Orient facets of a 2D mesh so that each connected component has positive volume.
void to_geogram_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, GEO::Mesh &M)
Converts a triangle mesh to a Geogram mesh.
void generate_edges(GEO::Mesh &M)
assing edges to M
void refine_polygonal_mesh(const GEO::Mesh &M_in, GEO::Mesh &M_out, Polygons::SplitFunction split_func)
Refine a polygonal 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
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13