PolyFEM
Loading...
Searching...
No Matches
NCMesh2D.cpp
Go to the documentation of this file.
2
4
5#include <igl/writeOBJ.h>
6
8
9namespace polyfem
10{
11 namespace mesh
12 {
13 bool NCMesh2D::is_boundary_element(const int element_global_id) const
14 {
15 assert(index_prepared);
16 for (int le = 0; le < n_face_vertices(element_global_id); le++)
17 if (is_boundary_edge(face_edge(element_global_id, le)))
18 return true;
19
20 return false;
21 }
22
23 void NCMesh2D::refine(const int n_refinement, const double t)
24 {
25 if (n_refinement <= 0)
26 return;
27 std::vector<bool> refine_mask(elements.size(), false);
28 for (int i = 0; i < elements.size(); i++)
29 if (elements[i].is_valid())
30 refine_mask[i] = true;
31
32 for (int i = 0; i < refine_mask.size(); i++)
33 if (refine_mask[i])
35
36 refine(n_refinement - 1, t);
37 }
38
39 bool NCMesh2D::build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
40 {
41 GEO::Mesh mesh_;
42 mesh_.clear(false, false);
43 to_geogram_mesh(V, F, mesh_);
44 orient_normals_2d(mesh_);
45
46 n_elements = 0;
47 elements.clear();
48 vertices.clear();
49 edges.clear();
50 midpointMap.clear();
51 edgeMap.clear();
52 refineHistory.clear();
53 index_prepared = false;
54 adj_prepared = false;
55
56 vertices.reserve(V.rows());
57 for (int i = 0; i < V.rows(); i++)
58 {
59 vertices.emplace_back(V.row(i));
60 }
61 for (int i = 0; i < F.rows(); i++)
62 {
63 add_element(F.row(i), -1);
64 }
65
67
68 return true;
69 }
70
71 void NCMesh2D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
72 {
73 for (int f = 0; f < n_faces(); ++f)
74 if (nodes[f].size() != 3)
75 throw std::runtime_error("NCMesh doesn't support high order mesh!");
76 }
77 std::pair<RowVectorNd, int> NCMesh2D::edge_node(const Navigation::Index &index, const int n_new_nodes, const int i) const
78 {
79 const auto v1 = point(index.vertex);
80 const auto v2 = point(switch_vertex(index).vertex);
81
82 const double t = i / (n_new_nodes + 1.0);
83
84 return std::make_pair((1 - t) * v1 + t * v2, -1);
85 }
86 std::pair<RowVectorNd, int> NCMesh2D::face_node(const Navigation::Index &index, const int n_new_nodes, const int i, const int j) const
87 {
88 const auto v1 = point(index.vertex);
89 const auto v2 = point(switch_vertex(index).vertex);
90 const auto v3 = point(switch_vertex(switch_edge(index)).vertex);
91
92 const double b2 = i / (n_new_nodes + 2.0);
93 const double b3 = j / (n_new_nodes + 2.0);
94 const double b1 = 1 - b3 - b2;
95 assert(b3 < 1);
96 assert(b3 > 0);
97
98 return std::make_pair(b1 * v1 + b2 * v2 + b3 * v3, -1);
99 }
100
101 int NCMesh2D::find_vertex(Eigen::Vector2i v) const
102 {
103 std::sort(v.data(), v.data() + v.size());
104 auto search = midpointMap.find(v);
105 if (search != midpointMap.end())
106 return search->second;
107 else
108 return -1;
109 }
110
111 int NCMesh2D::get_vertex(Eigen::Vector2i v)
112 {
113 std::sort(v.data(), v.data() + v.size());
114 int id = find_vertex(v);
115 if (id < 0)
116 {
117 Eigen::VectorXd v_mid = (vertices[v[0]].pos + vertices[v[1]].pos) / 2.;
118 id = vertices.size();
119 vertices.emplace_back(v_mid);
120 midpointMap.emplace(v, id);
121 }
122 return id;
123 }
124
125 int NCMesh2D::find_edge(Eigen::Vector2i v) const
126 {
127 std::sort(v.data(), v.data() + v.size());
128 auto search = edgeMap.find(v);
129 if (search != edgeMap.end())
130 return search->second;
131 else
132 return -1;
133 }
134
135 int NCMesh2D::get_edge(Eigen::Vector2i v)
136 {
137 std::sort(v.data(), v.data() + v.size());
138 int id = find_edge(v);
139 if (id < 0)
140 {
141 edges.emplace_back(v);
142 id = edges.size() - 1;
143 edgeMap.emplace(v, id);
144 }
145 return id;
146 }
147
148 int NCMesh2D::add_element(Eigen::Vector3i v, int parent)
149 {
150 const int id = elements.size();
151 const int level = (parent < 0) ? 0 : elements[parent].level + 1;
152 elements.emplace_back(2, v, level, parent);
153
154 if (parent >= 0)
155 elements[id].body_id = elements[parent].body_id;
156
157 for (int i = 0; i < v.size(); i++)
158 vertices[v(i)].n_elem++;
159
160 // add edges if not exist
161 int edge01 = get_edge(Eigen::Vector2i(v[0], v[1]));
162 int edge12 = get_edge(Eigen::Vector2i(v[2], v[1]));
163 int edge20 = get_edge(Eigen::Vector2i(v[0], v[2]));
164
165 elements[id].edges << edge01, edge12, edge20;
166
167 edges[edge01].add_element(id);
168 edges[edge12].add_element(id);
169 edges[edge20].add_element(id);
170
171 n_elements++;
172 index_prepared = false;
173 adj_prepared = false;
174
175 return id;
176 }
177
178 void NCMesh2D::refine_element(int id_full)
179 {
180 auto &elem = elements[id_full];
181 if (elem.is_not_valid())
182 throw std::runtime_error("Cannot refine an invalid element!");
183
184 const auto v = elem.vertices;
185 elem.is_refined = true;
186 n_elements--;
187
188 // remove the old element from edge reference
189 for (int e = 0; e < 3; e++)
190 edges[elem.edges(e)].remove_element(id_full);
191
192 for (int i = 0; i < v.size(); i++)
193 vertices[v(i)].n_elem--;
194
195 if (elem.children(0) >= 0)
196 {
197 for (int c = 0; c < elem.children.size(); c++)
198 {
199 auto &child = elements[elem.children(c)];
200 child.is_ghost = false;
201 n_elements++;
202 for (int le = 0; le < child.edges.size(); le++)
203 edges[child.edges(le)].add_element(child.children(c));
204 for (int i = 0; i < child.vertices.size(); i++)
205 vertices[child.vertices(i)].n_elem++;
206 }
207 }
208 else
209 {
210 // create mid-points if not exist
211 const int v01 = get_vertex(Eigen::Vector2i(v[0], v[1]));
212 const int v12 = get_vertex(Eigen::Vector2i(v[2], v[1]));
213 const int v20 = get_vertex(Eigen::Vector2i(v[0], v[2]));
214
215 // inherite line singularity flag from parent edge
216 for (int i = 0; i < v.size(); i++)
217 for (int j = 0; j < i; j++)
218 {
219 int mid_id = find_vertex(v[i], v[j]);
220 int edge_id = find_edge(v[i], v[j]);
221 int edge1 = get_edge(v[i], mid_id);
222 int edge2 = get_edge(v[j], mid_id);
223 edges[edge1].boundary_id = edges[edge_id].boundary_id;
224 edges[edge2].boundary_id = edges[edge_id].boundary_id;
225 }
226
227 // create and insert child elements
228 elements[id_full].children(0) = elements.size();
229 add_element(Eigen::Vector3i(v[0], v01, v20), id_full);
230 elements[id_full].children(1) = elements.size();
231 add_element(Eigen::Vector3i(v[1], v12, v01), id_full);
232 elements[id_full].children(2) = elements.size();
233 add_element(Eigen::Vector3i(v[2], v20, v12), id_full);
234 elements[id_full].children(3) = elements.size();
235 add_element(Eigen::Vector3i(v12, v20, v01), id_full);
236 }
237
238 refineHistory.push_back(id_full);
239
240 index_prepared = false;
241 adj_prepared = false;
242 }
243
244 void NCMesh2D::refine_elements(const std::vector<int> &ids)
245 {
246 std::vector<int> full_ids(ids.size());
247 for (int i = 0; i < ids.size(); i++)
248 full_ids[i] = valid_to_all_elem(ids[i]);
249
250 for (int i : full_ids)
252 }
253
255 {
256 const int parent_id = elements[id_full].parent;
257 auto &parent = elements[parent_id];
258
259 for (int i = 0; i < parent.children.size(); i++)
260 if (elements[parent.children(i)].is_not_valid())
261 throw std::runtime_error("Coarsen operation invalid!");
262
263 // remove elements
264 for (int i = 0; i < parent.children.size(); i++)
265 {
266 auto &elem = elements[parent.children(i)];
267 elem.is_ghost = true;
268 n_elements--;
269 for (int le = 0; le < elem.edges.size(); le++)
270 edges[elem.edges(le)].remove_element(parent.children(i));
271 for (int v = 0; v < elem.vertices.size(); v++)
272 vertices[elem.vertices(v)].n_elem--;
273 }
274
275 // add element
276 parent.is_refined = false;
277 n_elements++;
278 for (int le = 0; le < parent.edges.size(); le++)
279 edges[parent.edges(le)].add_element(parent_id);
280 for (int v = 0; v < parent.vertices.size(); v++)
281 vertices[parent.vertices(v)].n_elem++;
282
283 refineHistory.push_back(parent_id);
284
285 index_prepared = false;
286 adj_prepared = false;
287 }
288
289 int find(const Eigen::VectorXi &vec, int x)
290 {
291 for (int i = 0; i < vec.size(); i++)
292 {
293 if (x == vec[i])
294 return i;
295 }
296 return -1;
297 }
298
300 {
301 for (auto &edge : edges)
302 {
303 edge.leader = -1;
304 edge.followers.clear();
305 edge.weights.setConstant(-1);
306 }
307
308 Eigen::Vector2i v;
309 std::vector<follower_edge> followers;
310 for (int e_id = 0; e_id < elements.size(); e_id++)
311 {
312 const auto &element = elements[e_id];
313 if (element.is_not_valid())
314 continue;
315 for (int edge_local = 0; edge_local < 3; edge_local++)
316 {
317 v << element.vertices[edge_local], element.vertices[(edge_local + 1) % 3]; // order is important here!
318 int edge_global = element.edges[edge_local];
319 assert(edge_global >= 0);
320 traverse_edge(v, 0, 1, 0, followers);
321 for (auto &s : followers)
322 {
323 edges[s.id].leader = edge_global;
324 edges[edge_global].followers.push_back(s.id);
325 edges[s.id].weights << s.p1, s.p2;
326 }
327 followers.clear();
328 }
329 }
330 }
331
333 {
334 for (auto &edge : edges)
335 {
336 if (edge.n_elem() == 1)
337 edge.isboundary = true;
338 else
339 edge.isboundary = false;
340 }
341
342 for (auto &edge : edges)
343 {
344 if (edge.leader >= 0)
345 {
346 edge.isboundary = false;
347 edges[edge.leader].isboundary = false;
348 }
349 }
350
351 for (auto &vert : vertices)
352 vert.isboundary = false;
353
354 for (auto &edge : edges)
355 {
356 if (edge.isboundary && edge.n_elem())
357 {
358 for (int j = 0; j < 2; j++)
359 vertices[edge.vertices(j)].isboundary = true;
360 }
361 }
362 }
363
364 double line_weight(Eigen::Matrix<double, 2, 2> &e, Eigen::VectorXd &v)
365 {
366 assert(v.size() == 2);
367 double w1 = (v(0) - e(0, 0)) / (e(1, 0) - e(0, 0));
368 double w2 = (v(1) - e(0, 1)) / (e(1, 1) - e(0, 1));
369 if (0 <= w1 && w1 <= 1)
370 return w1;
371 else
372 return w2;
373 }
374
376 {
377 for (auto &vert : vertices)
378 {
379 vert.edge = -1;
380 vert.weight = -1;
381 }
382
383 Eigen::VectorXi vertexEdgeAdjacency;
384 vertexEdgeAdjacency.setConstant(vertices.size(), 1, -1);
385
386 for (int small_edge = 0; small_edge < edges.size(); small_edge++)
387 {
388 if (edges[small_edge].n_elem() == 0)
389 continue;
390
391 int large_edge = edges[small_edge].leader;
392 if (large_edge < 0)
393 continue;
394
395 int large_elem = edges[large_edge].get_element();
396 for (int j = 0; j < 2; j++)
397 {
398 int v_id = edges[small_edge].vertices(j);
399 if (find(elements[large_elem].vertices, v_id) < 0) // or maybe 0 < weights(large_edge, v_id) < 1
400 vertexEdgeAdjacency[v_id] = large_edge;
401 }
402 }
403
404 for (auto &element : elements)
405 {
406 if (element.is_not_valid())
407 continue;
408 for (int v_local = 0; v_local < 3; v_local++)
409 {
410 int v_global = element.vertices[v_local];
411 if (vertexEdgeAdjacency[v_global] < 0)
412 continue;
413
414 auto &large_edge = edges[vertexEdgeAdjacency[v_global]];
415 auto &large_element = elements[large_edge.get_element()];
416 vertices[v_global].edge = vertexEdgeAdjacency[v_global];
417
418 int e_local = find(large_element.edges, vertices[v_global].edge);
419 Eigen::Matrix<double, 2, 2> edge;
420 edge.row(0) = vertices[large_element.vertices[e_local]].pos;
421 edge.row(1) = vertices[large_element.vertices[(e_local + 1) % 3]].pos;
422 vertices[v_global].weight = line_weight(edge, vertices[v_global].pos);
423 }
424 }
425 }
426
427 double NCMesh2D::element_weight_to_edge_weight(const int l, const Eigen::Vector2d &pos)
428 {
429 double w = -1;
430 switch (l)
431 {
432 case 0:
433 w = pos(0);
434 assert(fabs(pos(1)) < 1e-12);
435 break;
436 case 1:
437 w = pos(1);
438 assert(fabs(pos(0) + pos(1) - 1) < 1e-12);
439 break;
440 case 2:
441 w = 1 - pos(1);
442 assert(fabs(pos(0)) < 1e-12);
443 break;
444 default:
445 assert(false);
446 }
447 return w;
448 }
449
451 {
452 all_to_valid_elemMap.assign(elements.size(), -1);
454
455 for (int i = 0, e = 0; i < elements.size(); i++)
456 {
457 if (elements[i].is_not_valid())
458 continue;
461 e++;
462 }
463
464 const int n_verts = n_vertices();
465
466 all_to_valid_vertexMap.assign(vertices.size(), -1);
467 valid_to_all_vertexMap.resize(n_verts);
468
469 for (int i = 0, j = 0; i < vertices.size(); i++)
470 {
471 if (vertices[i].n_elem == 0)
472 continue;
475 j++;
476 }
477
478 all_to_valid_edgeMap.assign(edges.size(), -1);
480
481 for (int i = 0, j = 0; i < edges.size(); i++)
482 {
483 if (edges[i].n_elem() == 0)
484 continue;
487 j++;
488 }
489 index_prepared = true;
490 }
491
492 void NCMesh2D::append(const Mesh &mesh)
493 {
494 assert(typeid(mesh) == typeid(NCMesh2D));
495 Mesh::append(mesh);
496
497 const NCMesh2D &mesh2d = dynamic_cast<const NCMesh2D &>(mesh);
498
499 const int n_v = n_vertices();
500 const int n_f = n_faces();
501
502 vertices.reserve(n_v + mesh2d.n_vertices());
503 for (int i = 0; i < mesh2d.n_vertices(); i++)
504 {
505 vertices.emplace_back(mesh2d.vertices[i].pos);
506 }
507 for (int i = 0; i < mesh2d.n_faces(); i++)
508 {
509 Eigen::Vector3i face = mesh2d.elements[i].vertices;
510 face = face.array() + n_v;
511 add_element(face, -1);
512 }
513
514 prepare_mesh();
515 }
516
517 std::unique_ptr<Mesh> NCMesh2D::copy() const
518 {
519 return std::make_unique<NCMesh2D>(*this);
520 }
521
522 void NCMesh2D::traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector<follower_edge> &list) const
523 {
524 int v_mid = find_vertex(v);
525 std::vector<follower_edge> list1, list2;
526 if (v_mid >= 0)
527 {
528 double p_mid = (p1 + p2) / 2;
529 traverse_edge(Eigen::Vector2i(v[0], v_mid), p1, p_mid, depth + 1, list1);
530 list.insert(
531 list.end(),
532 std::make_move_iterator(list1.begin()),
533 std::make_move_iterator(list1.end()));
534 traverse_edge(Eigen::Vector2i(v_mid, v[1]), p_mid, p2, depth + 1, list2);
535 list.insert(
536 list.end(),
537 std::make_move_iterator(list2.begin()),
538 std::make_move_iterator(list2.end()));
539 }
540 if (depth > 0)
541 {
542 int follower_id = find_edge(v);
543 if (follower_id >= 0 && edges[follower_id].n_elem() > 0)
544 list.emplace_back(follower_id, p1, p2);
545 }
546 }
547
548 bool NCMesh2D::load(const std::string &path)
549 {
550 assert(false);
551 return false;
552 }
553
554 bool NCMesh2D::load(const GEO::Mesh &mesh)
555 {
556 GEO::Mesh mesh_;
557 mesh_.clear(false, false);
558 mesh_.copy(mesh);
559 orient_normals_2d(mesh_);
560
561 Eigen::MatrixXd V(mesh_.vertices.nb(), 2);
562 Eigen::MatrixXi F(mesh_.facets.nb(), 3);
563
564 for (int v = 0; v < V.rows(); v++)
565 {
566 const double *ptr = mesh_.vertices.point_ptr(v);
567 V.row(v) << ptr[0], ptr[1];
568 }
569
570 for (int f = 0; f < F.rows(); f++)
571 for (int i = 0; i < F.cols(); i++)
572 F(f, i) = mesh_.facets.vertex(f, i);
573
574 n_elements = 0;
575 vertices.reserve(V.rows());
576 for (int i = 0; i < V.rows(); i++)
577 {
578 vertices.emplace_back(V.row(i));
579 }
580 for (int i = 0; i < F.rows(); i++)
581 {
582 add_element(F.row(i), -1);
583 }
584
585 prepare_mesh();
586
587 return true;
588 }
589
591 {
592 min = vertices[0].pos;
593 max = vertices[0].pos;
594
595 for (const auto &v : vertices)
596 {
597 for (int d = 0; d < 2; d++)
598 {
599 if (v.pos[d] > max[d])
600 max[d] = v.pos[d];
601 if (v.pos[d] < min[d])
602 min[d] = v.pos[d];
603 }
604 }
605 }
606
608 {
609 const auto &elem = elements[valid_to_all_elem(f)];
610
612 idx2.face = f;
613 idx2.vertex = all_to_valid_vertex(elem.vertices(lv));
614 idx2.edge = all_to_valid_edge(elem.edges(lv));
615 idx2.face_corner = -1;
616
617 return idx2;
618 }
619
621 {
622 const auto &elem = elements[valid_to_all_elem(idx.face)];
623 const auto &edge = edges[valid_to_all_edge(idx.edge)];
624
626 idx2.face = idx.face;
627 idx2.edge = idx.edge;
628
629 int v1 = valid_to_all_vertex(idx.vertex);
630 int v2 = -1;
631 for (int i = 0; i < edge.vertices.size(); i++)
632 if (edge.vertices(i) != v1)
633 {
634 v2 = edge.vertices(i);
635 break;
636 }
637
638 idx2.vertex = all_to_valid_vertex(v2);
639 idx2.face_corner = -1;
640
641 return idx2;
642 }
643
645 {
646 const auto &elem = elements[valid_to_all_elem(idx.face)];
647
649 idx2.face = idx.face;
650 idx2.vertex = idx.vertex;
651 idx2.face_corner = -1;
652
653 for (int i = 0; i < elem.edges.size(); i++)
654 {
655 const auto &edge = edges[elem.edges(i)];
656 const int valid_edge_id = all_to_valid_edge(elem.edges(i));
657 if (valid_edge_id != idx.edge && find(edge.vertices, idx.vertex) >= 0)
658 {
659 idx2.edge = valid_edge_id;
660 break;
661 }
662 }
663
664 return idx2;
665 }
666
668 {
670 idx2.edge = idx.edge;
671 idx2.vertex = idx.vertex;
672 idx2.face_corner = -1;
673
674 const auto &edge = edges[valid_to_all_edge(idx.edge)];
675 if (edge.n_elem() == 2)
676 idx2.face = (edge.find_opposite_element(valid_to_all_elem(idx.face)));
677 else
678 idx2.face = -1;
679
680 return idx2;
681 }
682
684 {
685 polyfem::RowVectorNd min, max;
686 bounding_box(min, max);
687
688 auto extent = max - min;
689 double scale = std::max(extent(0), extent(1));
690
691 for (auto &v : vertices)
692 v.pos = (v.pos - min.transpose()) / scale;
693 }
694
695 double NCMesh2D::edge_length(const int gid) const
696 {
697 const int v1 = edge_vertex(gid, 0);
698 const int v2 = edge_vertex(gid, 1);
699
700 return (point(v1) - point(v2)).norm();
701 }
702
711
712 void NCMesh2D::set_point(const int global_index, const RowVectorNd &p)
713 {
714 vertices[valid_to_all_vertex(global_index)].pos = p;
715 }
716
718 {
719 const int v1 = edge_vertex(index, 0);
720 const int v2 = edge_vertex(index, 1);
721
722 return 0.5 * (point(v1) + point(v2));
723 }
724
725 // void NCMesh2D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
726 // {
727 // ranges.clear();
728
729 // std::vector<Eigen::MatrixXi> local_tris(n_faces());
730 // std::vector<Eigen::MatrixXd> local_pts(n_faces());
731
732 // int total_tris = 0;
733 // int total_pts = 0;
734
735 // ranges.push_back(0);
736
737 // for (int f = 0; f < n_faces(); ++f)
738 // {
739 // const int n_vertices = n_face_vertices(f);
740
741 // Eigen::MatrixXd face_pts(n_vertices, 2);
742 // local_tris[f].resize(n_vertices - 2, 3);
743
744 // for (int i = 0; i < n_vertices; ++i)
745 // {
746 // const int vertex = face_vertex(f, i);
747 // auto pt = point(vertex);
748 // face_pts(i, 0) = pt[0];
749 // face_pts(i, 1) = pt[1];
750 // }
751
752 // for (int i = 1; i < n_vertices - 1; ++i)
753 // {
754 // local_tris[f].row(i - 1) << 0, i, i + 1;
755 // }
756
757 // local_pts[f] = face_pts;
758
759 // total_tris += local_tris[f].rows();
760 // total_pts += local_pts[f].rows();
761
762 // ranges.push_back(total_tris);
763
764 // assert(local_pts[f].rows() == face_pts.rows());
765 // }
766
767 // tris.resize(total_tris, 3);
768 // pts.resize(total_pts, 2);
769
770 // int tri_index = 0;
771 // int pts_index = 0;
772 // for (std::size_t i = 0; i < local_tris.size(); ++i)
773 // {
774 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
775 // tri_index += local_tris[i].rows();
776
777 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
778 // pts_index += local_pts[i].rows();
779 // }
780 // }
781
782 void NCMesh2D::set_body_ids(const std::vector<int> &body_ids)
783 {
784 assert(body_ids.size() == n_faces());
785 for (int i = 0; i < body_ids.size(); i++)
786 {
787 elements[valid_to_all_elem(i)].body_id = body_ids[i];
788 }
789 }
790
791 void NCMesh2D::set_boundary_ids(const std::vector<int> &boundary_ids)
792 {
793 assert(boundary_ids.size() == n_edges());
794 for (int i = 0; i < boundary_ids.size(); i++)
795 {
796 edges[valid_to_all_edge(i)].boundary_id = boundary_ids[i];
797 }
798 }
799
800 void NCMesh2D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
801 {
802 body_ids_.resize(n_faces());
803 std::fill(body_ids_.begin(), body_ids_.end(), -1);
804
805 for (int e = 0; e < n_faces(); ++e)
806 {
807 const auto bary = face_barycenter(e);
808 body_ids_[e] = marker(e, element_vertices(e), bary);
809 elements[valid_to_all_elem(e)].body_id = body_ids_[e];
810 }
811 }
812
813 void NCMesh2D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
814 {
815 boundary_ids_.resize(n_edges());
816
817 for (int e = 0; e < n_edges(); ++e)
818 {
819 bool is_boundary = is_boundary_edge(e);
820 const auto p = edge_barycenter(e);
821
822 std::vector<int> vs = {edge_vertex(e, 0), edge_vertex(e, 1)};
823 std::sort(vs.begin(), vs.end());
824 boundary_ids_[e] = marker(e, vs, p, is_boundary);
825 edges[valid_to_all_edge(e)].boundary_id = boundary_ids_[e];
826 }
827 }
828
829 } // namespace mesh
830} // namespace polyfem
int V
Eigen::MatrixXd vec
Definition Assembler.cpp:72
int edge_id
int x
RowVectorNd face_barycenter(const int index) const override
face barycenter
Definition Mesh2D.cpp:69
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:656
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:660
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
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:494
int find_edge(Eigen::Vector2i v) const
Definition NCMesh2D.cpp:125
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition NCMesh2D.cpp:13
void traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector< follower_edge > &list) const
Definition NCMesh2D.cpp:522
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 NCMesh2D.cpp:800
Navigation::Index switch_face(Navigation::Index idx) const override
Definition NCMesh2D.cpp:667
std::vector< int > valid_to_all_edgeMap
Definition NCMesh2D.hpp:397
RowVectorNd point(const int global_index) const override
point coordinates
Definition NCMesh2D.hpp:245
Navigation::Index switch_vertex(Navigation::Index idx) const override
Definition NCMesh2D.cpp:620
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
Definition NCMesh2D.cpp:517
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition NCMesh2D.cpp:71
int face_edge(const int f_id, const int le_id) const
Definition NCMesh2D.hpp:195
void refine_elements(const std::vector< int > &ids)
Definition NCMesh2D.cpp:244
void prepare_mesh() override
method used to finalize the mesh.
Definition NCMesh2D.hpp:288
std::vector< ncVert > vertices
Definition NCMesh2D.hpp:389
std::vector< int > all_to_valid_elemMap
Definition NCMesh2D.hpp:395
static double element_weight_to_edge_weight(const int l, const Eigen::Vector2d &pos)
Definition NCMesh2D.cpp:427
void refine_element(int id_full)
Definition NCMesh2D.cpp:178
void set_boundary_ids(const std::vector< int > &boundary_ids) override
Set the boundary selection from a vector.
Definition NCMesh2D.cpp:791
int n_faces() const override
number of faces
Definition NCMesh2D.hpp:168
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition NCMesh2D.cpp:712
int all_to_valid_vertex(const int id) const
Definition NCMesh2D.hpp:318
std::vector< ncElem > elements
Definition NCMesh2D.hpp:388
void set_body_ids(const std::vector< int > &body_ids) override
Set the volume sections.
Definition NCMesh2D.cpp:782
Navigation::Index switch_edge(Navigation::Index idx) const override
Definition NCMesh2D.cpp:644
int get_edge(Eigen::Vector2i v)
Definition NCMesh2D.cpp:135
int add_element(Eigen::Vector3i v, int parent=-1)
Definition NCMesh2D.cpp:148
void coarsen_element(int id_full)
Definition NCMesh2D.cpp:254
int get_vertex(Eigen::Vector2i v)
Definition NCMesh2D.cpp:111
Navigation::Index get_index_from_face(int f, int lv=0) const override
Definition NCMesh2D.cpp:607
double edge_length(const int gid) const override
edge length
Definition NCMesh2D.cpp:695
int valid_to_all_edge(const int id) const
Definition NCMesh2D.hpp:336
int n_vertices() const override
number of vertices
Definition NCMesh2D.hpp:177
bool is_boundary_edge(const int edge_global_id) const override
is edge boundary
Definition NCMesh2D.hpp:219
void refine(const int n_refinement, const double t) override
refine the mesh
Definition NCMesh2D.cpp:23
std::pair< RowVectorNd, int > face_node(const Navigation::Index &index, const int n_new_nodes, const int i, const int j) const override
Definition NCMesh2D.cpp:86
void normalize() override
normalize the mesh
Definition NCMesh2D.cpp:683
std::vector< int > refineHistory
Definition NCMesh2D.hpp:399
int valid_to_all_elem(const int id) const
Definition NCMesh2D.hpp:349
std::vector< int > valid_to_all_elemMap
Definition NCMesh2D.hpp:395
int edge_vertex(const int e_id, const int lv_id) const override
id of the edge vertex
Definition NCMesh2D.hpp:192
std::vector< int > all_to_valid_edgeMap
Definition NCMesh2D.hpp:397
std::vector< int > all_to_valid_vertexMap
Definition NCMesh2D.hpp:396
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > midpointMap
Definition NCMesh2D.hpp:392
std::vector< ncBoundary > edges
Definition NCMesh2D.hpp:390
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
Definition NCMesh2D.cpp:39
void compute_elements_tag() override
compute element types, see ElementType
Definition NCMesh2D.cpp:703
RowVectorNd edge_barycenter(const int index) const override
edge barycenter
Definition NCMesh2D.cpp:717
int all_to_valid_edge(const int id) const
Definition NCMesh2D.hpp:330
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > edgeMap
Definition NCMesh2D.hpp:393
int find_vertex(Eigen::Vector2i v) const
Definition NCMesh2D.cpp:101
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition NCMesh2D.cpp:590
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition NCMesh2D.hpp:187
int valid_to_all_vertex(const int id) const
Definition NCMesh2D.hpp:323
void build_element_vertex_adjacency()
Definition NCMesh2D.cpp:375
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 NCMesh2D.cpp:813
int n_edges() const override
number of edges
Definition NCMesh2D.hpp:169
std::vector< int > valid_to_all_vertexMap
Definition NCMesh2D.hpp:396
bool load(const std::string &path) override
loads a mesh from the path
Definition NCMesh2D.cpp:548
std::pair< RowVectorNd, int > edge_node(const Navigation::Index &index, const int n_new_nodes, const int i) const override
Definition NCMesh2D.cpp:77
void update_elements_tag() override
Update elements types.
Definition NCMesh2D.cpp:707
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition NCMesh2D.cpp:492
int find(const Eigen::VectorXi &vec, int x)
Definition NCMesh2D.cpp:289
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.
double line_weight(Eigen::Matrix< double, 2, 2 > &e, Eigen::VectorXd &v)
Definition NCMesh2D.cpp:364
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13