PolyFEM
Loading...
Searching...
No Matches
NCMesh3D.cpp
Go to the documentation of this file.
4
6
7#include <igl/writeMESH.h>
8
9#include <geogram/mesh/mesh_io.h>
10#include <fstream>
11
12using namespace polyfem::utils;
13
14namespace polyfem
15{
16 namespace mesh
17 {
18 int NCMesh3D::face_edge(const int f_id, const int le_id) const
19 {
20 const int v0 = faces[valid_to_all_face(f_id)].vertices(le_id);
21 const int v1 = faces[valid_to_all_face(f_id)].vertices((le_id + 1) % 3);
22 const int e_id = find_edge(v0, v1);
23 return all_to_valid_edge(e_id);
24 }
25
26 void NCMesh3D::refine(const int n_refinement, const double t)
27 {
28 if (n_refinement <= 0)
29 return;
30 std::vector<bool> refine_mask(elements.size(), false);
31 for (int i = 0; i < elements.size(); i++)
32 if (elements[i].is_valid())
33 refine_mask[i] = true;
34
35 for (int i = 0; i < refine_mask.size(); i++)
36 if (refine_mask[i])
38
39 refine(n_refinement - 1, t);
40 }
41
42 bool NCMesh3D::is_boundary_element(const int element_global_id) const
43 {
44 assert(index_prepared);
45 for (int lv = 0; lv < n_cell_edges(element_global_id); lv++)
46 if (is_boundary_vertex(cell_vertex(element_global_id, lv)))
47 return true;
48
49 return false;
50 }
51
52 bool NCMesh3D::build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
53 {
54 n_elements = 0;
55 elements.clear();
56 vertices.clear();
57 edges.clear();
58 midpointMap.clear();
59 edgeMap.clear();
60 faceMap.clear();
61 refineHistory.clear();
62 index_prepared = false;
63 adj_prepared = false;
64
65 vertices.reserve(V.rows());
66 for (int i = 0; i < V.rows(); i++)
67 {
68 vertices.emplace_back(V.row(i));
69 }
70 for (int i = 0; i < F.rows(); i++)
71 {
72 add_element(F.row(i), -1);
73 }
74
76
77 return true;
78 }
79
80 void NCMesh3D::attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes)
81 {
82 for (int f = 0; f < n_cells(); ++f)
83 if (nodes[f].size() != 4)
84 throw std::runtime_error("NCMesh doesn't support high order mesh!");
85 }
86
88 {
89 polyfem::RowVectorNd min, max;
90 bounding_box(min, max);
91
92 auto extent = max - min;
93 double scale = extent.maxCoeff();
94
95 for (auto &v : vertices)
96 v.pos = (v.pos - min.transpose()) / scale;
97 }
98
107 double NCMesh3D::edge_length(const int gid) const
108 {
109 const int v1 = edge_vertex(gid, 0);
110 const int v2 = edge_vertex(gid, 1);
111
112 return (point(v1) - point(v2)).norm();
113 }
114
115 RowVectorNd NCMesh3D::point(const int global_index) const
116 {
117 return vertices[valid_to_all_vertex(global_index)].pos.transpose();
118 }
119
120 void NCMesh3D::set_point(const int global_index, const RowVectorNd &p)
121 {
122 vertices[valid_to_all_vertex(global_index)].pos = p;
123 }
124
126 {
127 const int v1 = edge_vertex(e, 0);
128 const int v2 = edge_vertex(e, 1);
129
130 return 0.5 * (point(v1) + point(v2));
131 }
133 {
134 const int v1 = face_vertex(f, 0);
135 const int v2 = face_vertex(f, 1);
136 const int v3 = face_vertex(f, 2);
137
138 return (point(v1) + point(v2) + point(v3)) / 3.;
139 }
141 {
142 const int v1 = face_vertex(c, 0);
143 const int v2 = face_vertex(c, 1);
144 const int v3 = face_vertex(c, 2);
145
146 return (point(v1) + point(v2) + point(v3)) / 3.;
147 }
148
150 {
151 min = vertices[0].pos;
152 max = vertices[0].pos;
153
154 for (const auto &v : vertices)
155 {
156 for (int d = 0; d < 3; d++)
157 {
158 if (v.pos[d] > max[d])
159 max[d] = v.pos[d];
160 if (v.pos[d] < min[d])
161 min[d] = v.pos[d];
162 }
163 }
164 }
165
166 void NCMesh3D::compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker)
167 {
168 boundary_ids_.resize(n_faces());
169 std::fill(boundary_ids_.begin(), boundary_ids_.end(), -1);
170
171 for (int f = 0; f < n_faces(); ++f)
172 {
173 const bool is_boundary = is_boundary_face(f);
174 std::vector<int> vs(n_face_vertices(f));
175 const auto p = face_barycenter(f);
176
177 for (int vid = 0; vid < vs.size(); ++vid)
178 vs[vid] = face_vertex(f, vid);
179
180 std::sort(vs.begin(), vs.end());
181 boundary_ids_[f] = marker(f, vs, p, is_boundary);
182
183 faces[valid_to_all_face(f)].boundary_id = boundary_ids_[f];
184 }
185 }
186
187 void NCMesh3D::compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker)
188 {
189 body_ids_.resize(n_cells());
190 std::fill(body_ids_.begin(), body_ids_.end(), -1);
191
192 for (int e = 0; e < n_cells(); ++e)
193 {
194 const auto bary = cell_barycenter(e);
195 body_ids_[e] = marker(e, element_vertices(e), bary);
196 elements[valid_to_all_elem(e)].body_id = body_ids_[e];
197 }
198 }
199 void NCMesh3D::set_boundary_ids(const std::vector<int> &boundary_ids)
200 {
201 assert(boundary_ids.size() == n_faces());
202 for (int i = 0; i < boundary_ids.size(); i++)
203 {
204 faces[valid_to_all_face(i)].boundary_id = boundary_ids[i];
205 }
206 }
207 void NCMesh3D::set_body_ids(const std::vector<int> &body_ids)
208 {
209 assert(body_ids.size() == n_cells());
210 for (int i = 0; i < body_ids.size(); i++)
211 {
212 elements[valid_to_all_elem(i)].body_id = body_ids[i];
213 }
214 }
215
216 // void NCMesh3D::triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const
217 // {
218 // ranges.clear();
219
220 // std::vector<Eigen::MatrixXi> local_tris(n_cells());
221 // std::vector<Eigen::MatrixXd> local_pts(n_cells());
222 // Eigen::MatrixXi tets;
223
224 // int total_tris = 0;
225 // int total_pts = 0;
226
227 // ranges.push_back(0);
228
229 // Eigen::MatrixXd face_barys;
230 // face_barycenters(face_barys);
231
232 // Eigen::MatrixXd cell_barys;
233 // cell_barycenters(cell_barys);
234
235 // for (std::size_t e = 0; e < n_cells(); ++e)
236 // {
237 // const int n_vertices = n_cell_vertices(e);
238 // const int n_faces = n_cell_faces(e);
239
240 // Eigen::MatrixXd local_pt(n_vertices + n_faces, 3);
241
242 // std::map<int, int> global_to_local;
243
244 // for (int i = 0; i < n_vertices; ++i)
245 // {
246 // const int global_index = cell_vertex(e, i);
247 // local_pt.row(i) = point(global_index);
248 // global_to_local[global_index] = i;
249 // }
250
251 // int n_local_faces = 0;
252 // for (int i = 0; i < n_faces; ++i)
253 // {
254 // const int f_id = cell_face(e, i);
255 // n_local_faces += n_face_vertices(f_id);
256
257 // local_pt.row(n_vertices + i) = face_barys.row(f_id);
258 // }
259
260 // Eigen::MatrixXi local_faces(n_local_faces, 3);
261
262 // int face_index = 0;
263 // for (int i = 0; i < n_faces; ++i)
264 // {
265 // const int f_id = cell_face(e, i);
266
267 // const Eigen::RowVector3d e0 = (point(face_vertex(f_id, 0)) - local_pt.row(n_vertices + i));
268 // const Eigen::RowVector3d e1 = (point(face_vertex(f_id, 1)) - local_pt.row(n_vertices + i));
269 // const Eigen::RowVector3d normal = e0.cross(e1);
270 // // const Eigen::RowVector3d check_dir = (node_from_element(e)-p);
271 // const Eigen::RowVector3d check_dir = (cell_barys.row(e) - point(face_vertex(f_id, 1)));
272
273 // const bool reverse_order = normal.dot(check_dir) > 0;
274
275 // for (int j = 0; j < n_face_vertices(f_id); ++j)
276 // {
277 // const int jp = (j + 1) % n_face_vertices(f_id);
278 // if (reverse_order)
279 // {
280 // local_faces(face_index, 0) = global_to_local[face_vertex(f_id, jp)];
281 // local_faces(face_index, 1) = global_to_local[face_vertex(f_id, j)];
282 // }
283 // else
284 // {
285 // local_faces(face_index, 0) = global_to_local[face_vertex(f_id, j)];
286 // local_faces(face_index, 1) = global_to_local[face_vertex(f_id, jp)];
287 // }
288 // local_faces(face_index, 2) = n_vertices + i;
289
290 // ++face_index;
291 // }
292 // }
293
294 // local_pts[e] = local_pt;
295 // local_tris[e] = local_faces;
296
297 // total_tris += local_tris[e].rows();
298 // total_pts += local_pts[e].rows();
299
300 // ranges.push_back(total_tris);
301
302 // assert(local_pts[e].rows() == local_pt.rows());
303 // }
304
305 // tris.resize(total_tris, 3);
306 // pts.resize(total_pts, 3);
307
308 // int tri_index = 0;
309 // int pts_index = 0;
310 // for (std::size_t i = 0; i < local_tris.size(); ++i)
311 // {
312 // tris.block(tri_index, 0, local_tris[i].rows(), local_tris[i].cols()) = local_tris[i].array() + pts_index;
313 // tri_index += local_tris[i].rows();
314
315 // pts.block(pts_index, 0, local_pts[i].rows(), local_pts[i].cols()) = local_pts[i];
316 // pts_index += local_pts[i].rows();
317 // }
318 // }
319
320 RowVectorNd NCMesh3D::kernel(const int cell_id) const
321 {
322 assert(false);
323 return RowVectorNd();
324 }
326 {
328
329 idx.element = hi;
330
331 idx.element_patch = lf;
332 idx.face = cell_face(hi, lf);
333
334 if (lv >= 4)
335 lv = lv % 4;
336 idx.face_corner = lv;
337 idx.vertex = face_vertex(idx.face, idx.face_corner);
338
339 idx.edge = face_edge(idx.face, idx.face_corner);
340
341 return idx;
342 }
344 {
346
347 idx.element = hi;
348 idx.element_patch = 0;
349 idx.face = cell_face(hi, 0);
350 idx.face_corner = 0;
351 idx.vertex = face_vertex(idx.face, idx.face_corner);
352 idx.edge = face_edge(idx.face, 0);
353
354 return idx;
355 }
356
358 {
360
361 idx.element = hi;
362 idx.vertex = v0;
364
365 for (int i = 0; i < 4; i++)
366 {
367 const int f_id = cell_face(idx.element, i);
368 for (int j = 0; j < 3; j++)
369 {
370 const int e_id = face_edge(f_id, j);
371 if (idx.edge == e_id)
372 {
373 idx.element_patch = i;
374 idx.face = f_id;
375
376 for (int d = 0; d < n_face_vertices(f_id); d++)
377 if (face_vertex(f_id, d) == idx.vertex)
378 idx.face_corner = d;
379
380 assert(switch_vertex(idx).vertex == v1);
381
382 return idx;
383 }
384 }
385 }
386
387 assert(false);
388 return idx;
389 }
390 Navigation3D::Index NCMesh3D::get_index_from_element_face(int hi, int v0, int v1, int v2) const
391 {
393
394 idx.element = hi;
395 idx.vertex = v0;
397
398 for (int d = 0; d < n_face_vertices(idx.face); d++)
399 if (face_vertex(idx.face, d) == idx.vertex)
400 idx.face_corner = d;
401
402 int v0_ = v0, v1_ = v1;
403
404 // let v0_ < v1_
405 if (v0_ > v1_)
406 std::swap(v0_, v1_);
407
408 for (int i = 0; i < 4; i++)
409 {
410 const int fid = cell_face(idx.element, i);
411 if (fid != idx.face)
412 continue;
413
414 idx.element_patch = i;
415
416 for (int j = 0; j < 3; j++)
417 {
418 const int eid = face_edge(fid, j);
419 const int ev0 = edge_vertex(eid, 0);
420 const int ev1 = edge_vertex(eid, 1);
421
422 if ((ev0 == v0_ && ev1 == v1_) || (ev0 == v1_ && ev1 == v0_))
423 {
424 idx.edge = eid;
425
426 assert(switch_vertex(idx).vertex == v1);
427 assert(switch_vertex(switch_edge(idx)).vertex == v2);
428
429 return idx;
430 }
431 }
432 }
433
434 assert(false);
435 return idx;
436 }
437
438 std::vector<uint32_t> NCMesh3D::vertex_neighs(const int v_gid) const
439 {
440 std::vector<uint32_t> hs;
441 for (auto h : vertices[valid_to_all_vertex(v_gid)].elem_list)
442 hs.push_back(all_to_valid_elem(h));
443
444 return hs;
445 }
446 std::vector<uint32_t> NCMesh3D::edge_neighs(const int e_gid) const
447 {
448 std::vector<uint32_t> hs;
449 for (auto h : edges[valid_to_all_edge(e_gid)].elem_list)
450 hs.push_back(all_to_valid_elem(h));
451
452 return hs;
453 }
454
456 {
457 if (idx.vertex == edge_vertex(idx.edge, 0))
458 idx.vertex = edge_vertex(idx.edge, 1);
459 else
460 idx.vertex = edge_vertex(idx.edge, 0);
461
462 for (int d = 0; d < n_face_vertices(idx.face); d++)
463 if (face_vertex(idx.face, d) == idx.vertex)
464 idx.face_corner = d;
465
466 return idx;
467 }
469 {
470 if (idx.edge == face_edge(idx.face, idx.face_corner))
471 idx.edge = face_edge(idx.face, (idx.face_corner + 2) % 3);
472 else
473 idx.edge = face_edge(idx.face, idx.face_corner);
474 return idx;
475 }
477 {
478 for (int i = 0; i < 4; i++)
479 {
480 const int fid = cell_face(idx.element, i);
481 if (fid == idx.face)
482 continue;
483 if (idx.edge == face_edge(fid, 0) || idx.edge == face_edge(fid, 1) || idx.edge == face_edge(fid, 2))
484 {
485 idx.face = fid;
486 idx.element_patch = i;
487
488 for (int d = 0; d < n_face_vertices(fid); d++)
489 if (face_vertex(fid, d) == idx.vertex)
490 idx.face_corner = d;
491
492 break;
493 }
494 }
495 return idx;
496 }
498 {
499 const auto &face = faces[valid_to_all_face(idx.face)];
500 if (face.n_elem() != 2)
501 {
502 idx.element = -1;
503 return idx;
504 }
505 if (all_to_valid_elem(face.get_element()) == idx.element)
506 idx.element = all_to_valid_elem(face.find_opposite_element(face.get_element()));
507 else
508 idx.element = all_to_valid_elem(face.get_element());
509
510 for (int f = 0; f < n_cell_faces(idx.element); f++)
511 if (cell_face(idx.element, f) == idx.face)
512 idx.element_patch = f;
513 return idx;
514 }
515
524
525 void NCMesh3D::get_vertex_elements_neighs(const int v_id, std::vector<int> &ids) const
526 {
527 ids.clear();
528 for (auto h : vertices[valid_to_all_vertex(v_id)].elem_list)
529 ids.push_back(all_to_valid_elem(h));
530 }
531 void NCMesh3D::get_edge_elements_neighs(const int e_id, std::vector<int> &ids) const
532 {
533 ids.clear();
534 for (auto h : edges[valid_to_all_edge(e_id)].elem_list)
535 ids.push_back(all_to_valid_elem(h));
536 }
537 void NCMesh3D::get_face_elements_neighs(const int f_id, std::vector<int> &ids) const
538 {
539 ids.clear();
540 for (auto h : faces[valid_to_all_face(f_id)].elem_list)
541 ids.push_back(all_to_valid_elem(h));
542 }
543
544 void NCMesh3D::refine_element(int id_full)
545 {
546 assert(elements[id_full].is_valid());
547 if (elements[id_full].is_not_valid())
548 throw std::runtime_error("Cannot refine an invalid element!");
549
550 const auto v = elements[id_full].vertices;
551 elements[id_full].is_refined = true;
552 n_elements--;
553
554 for (int f = 0; f < elements[id_full].faces.size(); f++)
555 faces[elements[id_full].faces(f)].remove_element(id_full);
556
557 for (int e = 0; e < elements[id_full].edges.size(); e++)
558 edges[elements[id_full].edges(e)].remove_element(id_full);
559
560 for (int i = 0; i < v.size(); i++)
561 vertices[v(i)].remove_element(id_full);
562
563 if (elements[id_full].children(0) >= 0)
564 {
565 for (int c = 0; c < elements[id_full].children.size(); c++)
566 {
567 const int child_id = elements[id_full].children(c);
568 auto &elem = elements[child_id];
569 elem.is_ghost = false;
570 n_elements++;
571
572 for (int f = 0; f < elem.faces.size(); f++)
573 faces[elem.faces(f)].add_element(child_id);
574
575 for (int e = 0; e < elem.edges.size(); e++)
576 edges[elem.edges(e)].add_element(child_id);
577
578 for (int v = 0; v < elem.vertices.size(); v++)
579 vertices[elem.vertices(v)].add_element(child_id);
580 }
581 }
582 else
583 {
584 // create mid-points if not exist
585 const int v1 = v[0];
586 const int v2 = v[1];
587 const int v3 = v[2];
588 const int v4 = v[3];
589 const int v5 = get_vertex(Eigen::Vector2i(v1, v2));
590 const int v8 = get_vertex(Eigen::Vector2i(v3, v2));
591 const int v6 = get_vertex(Eigen::Vector2i(v1, v3));
592 const int v7 = get_vertex(Eigen::Vector2i(v1, v4));
593 const int v9 = get_vertex(Eigen::Vector2i(v4, v2));
594 const int v10 = get_vertex(Eigen::Vector2i(v3, v4));
595
596 // inherite line singularity flag from parent edge
597 for (int i = 0; i < v.size(); i++)
598 for (int j = 0; j < i; j++)
599 {
600 int mid_id = find_vertex(v[i], v[j]);
601 int edge_id = find_edge(v[i], v[j]);
602 int edge1 = get_edge(v[i], mid_id);
603 int edge2 = get_edge(v[j], mid_id);
604 edges[edge1].boundary_id = edges[edge_id].boundary_id;
605 edges[edge2].boundary_id = edges[edge_id].boundary_id;
606 }
607
608 for (int i = 0; i < v.size(); i++)
609 for (int j = 0; j < i; j++)
610 for (int k = 0; k < j; k++)
611 {
612 int fid = find_face(v[i], v[j], v[k]);
613 int vij = find_vertex(v[i], v[j]);
614 int vjk = find_vertex(v[j], v[k]);
615 int vik = find_vertex(v[i], v[k]);
616 int facei = get_face(v[i], vij, vik);
617 int facej = get_face(v[j], vjk, vij);
618 int facek = get_face(v[k], vjk, vik);
619 int facem = get_face(vij, vjk, vik);
620 faces[facei].boundary_id = faces[fid].boundary_id;
621 faces[facej].boundary_id = faces[fid].boundary_id;
622 faces[facek].boundary_id = faces[fid].boundary_id;
623 faces[facem].boundary_id = faces[fid].boundary_id;
624 }
625
626 // create children
627 elements[id_full].children(0) = elements.size();
628 add_element(Eigen::Vector4i(v1, v5, v6, v7), id_full);
629 elements[id_full].children(1) = elements.size();
630 add_element(Eigen::Vector4i(v5, v2, v8, v9), id_full);
631 elements[id_full].children(2) = elements.size();
632 add_element(Eigen::Vector4i(v6, v8, v3, v10), id_full);
633 elements[id_full].children(3) = elements.size();
634 add_element(Eigen::Vector4i(v7, v9, v10, v4), id_full);
635 elements[id_full].children(4) = elements.size();
636 add_element(Eigen::Vector4i(v5, v6, v7, v9), id_full);
637 elements[id_full].children(5) = elements.size();
638 add_element(Eigen::Vector4i(v5, v9, v8, v6), id_full);
639 elements[id_full].children(6) = elements.size();
640 add_element(Eigen::Vector4i(v6, v7, v9, v10), id_full);
641 elements[id_full].children(7) = elements.size();
642 add_element(Eigen::Vector4i(v6, v10, v9, v8), id_full);
643 }
644
645 refineHistory.push_back(id_full);
646 }
647 void NCMesh3D::refine_elements(const std::vector<int> &ids)
648 {
649 std::vector<int> full_ids(ids.size());
650 for (int i = 0; i < ids.size(); i++)
651 full_ids[i] = valid_to_all_elem(ids[i]);
652
653 for (int i : full_ids)
655 }
656
658 {
659 const int parent_id = elements[id_full].parent;
660 auto &parent = elements[parent_id];
661
662 for (int i = 0; i < parent.children.size(); i++)
663 if (elements[parent.children(i)].is_not_valid())
664 throw std::runtime_error("Invalid siblings in coarsening!");
665
666 // remove elements
667 for (int c = 0; c < parent.children.size(); c++)
668 {
669 auto &elem = elements[parent.children(c)];
670 elem.is_ghost = true;
671 n_elements--;
672
673 for (int f = 0; f < elem.faces.size(); f++)
674 faces[elem.faces(f)].remove_element(parent.children(c));
675
676 for (int e = 0; e < elem.edges.size(); e++)
677 edges[elem.edges(e)].remove_element(parent.children(c));
678
679 for (int v = 0; v < elem.vertices.size(); v++)
680 vertices[elem.vertices(v)].remove_element(parent.children(c));
681 }
682
683 // add element
684 parent.is_refined = false;
685 n_elements++;
686
687 for (int f = 0; f < parent.faces.size(); f++)
688 faces[parent.faces(f)].add_element(parent_id);
689
690 for (int e = 0; e < parent.edges.size(); e++)
691 edges[parent.edges(e)].add_element(parent_id);
692
693 for (int v = 0; v < parent.vertices.size(); v++)
694 vertices[parent.vertices(v)].add_element(parent_id);
695
696 refineHistory.push_back(parent_id);
697 }
698
700 {
701 for (auto &face : faces)
702 {
703 if (face.n_elem() == 1)
704 face.isboundary = true;
705 else
706 face.isboundary = false;
707 }
708
709 for (auto &face : faces)
710 {
711 if (face.leader >= 0)
712 {
713 face.isboundary = false;
714 faces[face.leader].isboundary = false;
715 }
716 }
717
718 for (auto &vert : vertices)
719 vert.isboundary = false;
720
721 for (auto &edge : edges)
722 edge.isboundary = false;
723
724 for (auto &face : faces)
725 {
726 if (face.isboundary && face.n_elem())
727 {
728 for (int j = 0; j < 3; j++)
729 vertices[face.vertices(j)].isboundary = true;
730
731 for (int j = 0; j < 3; j++)
732 for (int i = 0; i < j; i++)
733 edges[find_edge(face.vertices(i), face.vertices(j))].isboundary = true;
734 }
735 }
736
737 // for (int v = 0; v < n_vertices(); v++)
738 // std::cout << v << ": " << point(v) << "\n";
739
740 // for (int f = 0; f < n_faces(); f++)
741 // {
742 // for (int v = 0; v < n_face_vertices(f); v++)
743 // std::cout << face_vertex(f, v) << ", ";
744 // std::cout << faces[valid_to_all_face(f)].isboundary << std::endl;
745 // }
746 }
747
749 {
750 all_to_valid_elemMap.assign(elements.size(), -1);
752
753 for (int i = 0, e = 0; i < elements.size(); i++)
754 {
755 if (elements[i].is_not_valid())
756 continue;
759 e++;
760 }
761
762 const int n_verts = n_vertices();
763
764 all_to_valid_vertexMap.assign(vertices.size(), -1);
765 valid_to_all_vertexMap.resize(n_verts);
766
767 for (int i = 0, j = 0; i < vertices.size(); i++)
768 {
769 if (vertices[i].n_elem() == 0)
770 continue;
773 j++;
774 }
775
776 all_to_valid_edgeMap.assign(edges.size(), -1);
778
779 for (int i = 0, j = 0; i < edges.size(); i++)
780 {
781 if (edges[i].n_elem() == 0)
782 continue;
785 j++;
786 }
787
788 all_to_valid_faceMap.assign(faces.size(), -1);
790
791 for (int i = 0, j = 0; i < faces.size(); i++)
792 {
793 if (faces[i].n_elem() == 0)
794 continue;
797 j++;
798 }
799 index_prepared = true;
800 }
801
802 void NCMesh3D::append(const Mesh &mesh)
803 {
804 assert(typeid(mesh) == typeid(NCMesh3D));
805 Mesh::append(mesh);
806
807 const NCMesh3D &mesh3d = dynamic_cast<const NCMesh3D &>(mesh);
808
809 const int n_v = n_vertices();
810 const int n_f = n_cells();
811
812 vertices.reserve(n_v + mesh3d.n_vertices());
813 for (int i = 0; i < mesh3d.n_vertices(); i++)
814 {
815 vertices.emplace_back(mesh3d.vertices[i].pos);
816 }
817 for (int i = 0; i < mesh3d.n_cells(); i++)
818 {
819 Eigen::Vector4i cell = mesh3d.elements[i].vertices;
820 cell = cell.array() + n_v;
821 add_element(cell, -1);
822 }
823
824 prepare_mesh();
825 }
826
827 std::unique_ptr<Mesh> NCMesh3D::copy() const
828 {
829 return std::make_unique<NCMesh3D>(*this);
830 }
831
832 bool NCMesh3D::load(const std::string &path)
833 {
834 if (!StringUtils::endswith(path, ".HYBRID"))
835 {
836 GEO::Mesh M;
837 GEO::mesh_load(path, M);
838 return load(M);
839 }
840 return false;
841 }
842 bool NCMesh3D::load(const GEO::Mesh &M)
843 {
844 assert(M.vertices.dimension() == 3);
845
846 Eigen::MatrixXd V(M.vertices.nb(), 2);
847 Eigen::MatrixXi C(M.cells.nb(), 4);
848
849 for (int v = 0; v < V.rows(); v++)
850 V.row(v) << M.vertices.point(v)[0], M.vertices.point(v)[1], M.vertices.point(v)[2];
851
852 for (int c = 0; c < C.rows(); c++)
853 {
854 if (M.cells.type(c) != GEO::MESH_TET)
855 throw std::runtime_error("NCMesh3D only supports tet mesh!");
856 for (int i = 0; i < C.cols(); i++)
857 C(c, i) = M.cells.vertex(c, i);
858 }
859
860 n_elements = 0;
861 vertices.reserve(V.rows());
862 for (int i = 0; i < V.rows(); i++)
863 {
864 vertices.emplace_back(V.row(i));
865 }
866 for (int i = 0; i < C.rows(); i++)
867 {
868 add_element(C.row(i), -1);
869 }
870
871 prepare_mesh();
872
873 return true;
874 }
875
876 int NCMesh3D::find_vertex(Eigen::Vector2i v) const
877 {
878 std::sort(v.data(), v.data() + v.size());
879 auto search = midpointMap.find(v);
880 if (search != midpointMap.end())
881 return search->second;
882 else
883 return -1;
884 }
885 int NCMesh3D::get_vertex(Eigen::Vector2i v)
886 {
887 std::sort(v.data(), v.data() + v.size());
888 int id = find_vertex(v);
889 if (id < 0)
890 {
891 Eigen::VectorXd v_mid = (vertices[v[0]].pos + vertices[v[1]].pos) / 2.;
892 id = vertices.size();
893 vertices.emplace_back(v_mid);
894 midpointMap.emplace(v, id);
895 }
896 return id;
897 }
898 int NCMesh3D::find_edge(Eigen::Vector2i v) const
899 {
900 std::sort(v.data(), v.data() + v.size());
901 auto search = edgeMap.find(v);
902 if (search != edgeMap.end())
903 return search->second;
904 else
905 return -1;
906 }
907 int NCMesh3D::get_edge(Eigen::Vector2i v)
908 {
909 std::sort(v.data(), v.data() + v.size());
910 int id = find_edge(v);
911 if (id < 0)
912 {
913 edges.emplace_back(v);
914 id = edges.size() - 1;
915 edgeMap.emplace(v, id);
916 }
917 return id;
918 }
919 int NCMesh3D::find_face(Eigen::Vector3i v) const
920 {
921 std::sort(v.data(), v.data() + v.size());
922 auto search = faceMap.find(v);
923 if (search != faceMap.end())
924 return search->second;
925 else
926 return -1;
927 }
928 int NCMesh3D::get_face(Eigen::Vector3i v)
929 {
930 std::sort(v.data(), v.data() + v.size());
931 int id = find_face(v);
932 if (id < 0)
933 {
934 faces.emplace_back(v);
935 id = faces.size() - 1;
936 faceMap.emplace(v, id);
937 }
938 return id;
939 }
940 void NCMesh3D::traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector<follower_edge> &list) const
941 {
942 int v_mid = find_vertex(v);
943 std::vector<follower_edge> list1, list2;
944 if (v_mid >= 0)
945 {
946 double p_mid = (p1 + p2) / 2;
947 traverse_edge(Eigen::Vector2i(v[0], v_mid), p1, p_mid, depth + 1, list1);
948 list.insert(
949 list.end(),
950 std::make_move_iterator(list1.begin()),
951 std::make_move_iterator(list1.end()));
952 traverse_edge(Eigen::Vector2i(v_mid, v[1]), p_mid, p2, depth + 1, list2);
953 list.insert(
954 list.end(),
955 std::make_move_iterator(list2.begin()),
956 std::make_move_iterator(list2.end()));
957 }
958 if (depth > 0)
959 {
960 int follower_id = find_edge(v);
961 if (follower_id >= 0 && edges[follower_id].n_elem() > 0)
962 list.emplace_back(follower_id, p1, p2);
963 }
964 }
966 {
967 for (auto &edge : edges)
968 {
969 edge.leader = -1;
970 edge.followers.clear();
971 edge.weights.setConstant(-1);
972 }
973
974 for (int e_id = 0; e_id < edges.size(); e_id++)
975 {
976 auto &edge = edges[e_id];
977 if (edge.n_elem() == 0)
978 continue;
979 std::vector<follower_edge> followers;
980 traverse_edge(edge.vertices, 0, 1, 0, followers);
981 for (auto &s : followers)
982 {
983 if (edges[s.id].leader >= 0 && std::abs(edges[s.id].weights(1) - edges[s.id].weights(0)) < std::abs(s.p2 - s.p1))
984 continue;
985 edge.followers.push_back(s.id);
986 edges[s.id].leader = e_id;
987 edges[s.id].weights << s.p1, s.p2;
988 }
989 }
990
991 // In 3d, it's possible for one edge to have both leader and follower edges, but we don't care this case.
992 for (auto &edge : edges)
993 if (edge.leader >= 0 && edge.followers.size())
994 edge.followers.clear();
995 }
996 void NCMesh3D::traverse_face(int v1, int v2, int v3, Eigen::Vector2d p1, Eigen::Vector2d p2, Eigen::Vector2d p3, int depth, std::vector<follower_face> &face_list, std::vector<int> &edge_list) const
997 {
998 int v12 = find_vertex(Eigen::Vector2i(v1, v2));
999 int v23 = find_vertex(Eigen::Vector2i(v3, v2));
1000 int v31 = find_vertex(Eigen::Vector2i(v1, v3));
1001 std::vector<follower_face> list1, list2, list3, list4;
1002 std::vector<int> list1_, list2_, list3_, list4_;
1003 if (depth > 0)
1004 {
1005 edge_list.push_back(find_edge(v1, v2));
1006 edge_list.push_back(find_edge(v3, v2));
1007 edge_list.push_back(find_edge(v1, v3));
1008 }
1009 if (v12 >= 0 && v23 >= 0 && v31 >= 0)
1010 {
1011 auto p12 = (p1 + p2) / 2, p23 = (p2 + p3) / 2, p31 = (p1 + p3) / 2;
1012 traverse_face(v1, v12, v31, p1, p12, p31, depth + 1, list1, list1_);
1013 traverse_face(v12, v2, v23, p12, p2, p23, depth + 1, list2, list2_);
1014 traverse_face(v31, v23, v3, p31, p23, p3, depth + 1, list3, list3_);
1015 traverse_face(v12, v23, v31, p12, p23, p31, depth + 1, list4, list4_);
1016
1017 face_list.insert(face_list.end(), std::make_move_iterator(list1.begin()), std::make_move_iterator(list1.end()));
1018 face_list.insert(face_list.end(), std::make_move_iterator(list2.begin()), std::make_move_iterator(list2.end()));
1019 face_list.insert(face_list.end(), std::make_move_iterator(list3.begin()), std::make_move_iterator(list3.end()));
1020 face_list.insert(face_list.end(), std::make_move_iterator(list4.begin()), std::make_move_iterator(list4.end()));
1021
1022 edge_list.insert(edge_list.end(), std::make_move_iterator(list1_.begin()), std::make_move_iterator(list1_.end()));
1023 edge_list.insert(edge_list.end(), std::make_move_iterator(list2_.begin()), std::make_move_iterator(list2_.end()));
1024 edge_list.insert(edge_list.end(), std::make_move_iterator(list3_.begin()), std::make_move_iterator(list3_.end()));
1025 edge_list.insert(edge_list.end(), std::make_move_iterator(list4_.begin()), std::make_move_iterator(list4_.end()));
1026 }
1027 if (depth > 0)
1028 {
1029 int follower_id = find_face(Eigen::Vector3i(v1, v2, v3));
1030 if (follower_id >= 0 && faces[follower_id].n_elem() > 0)
1031 face_list.emplace_back(follower_id, p1, p2, p3);
1032 }
1033 }
1035 {
1036 Eigen::Matrix<int, 4, 3> fv;
1037 fv.row(0) << 0, 1, 2;
1038 fv.row(1) << 0, 1, 3;
1039 fv.row(2) << 1, 2, 3;
1040 fv.row(3) << 2, 0, 3;
1041
1042 for (auto &face : faces)
1043 {
1044 face.leader = -1;
1045 face.followers.clear();
1046 }
1047
1048 for (auto &edge : edges)
1049 {
1050 edge.leader_face = -1;
1051 }
1052
1053 for (int f_id = 0; f_id < faces.size(); f_id++)
1054 {
1055 auto &face = faces[f_id];
1056 if (face.n_elem() == 0)
1057 continue;
1058 std::vector<follower_face> followers;
1059 std::vector<int> interior_edges;
1060 traverse_face(face.vertices(0), face.vertices(1), face.vertices(2), Eigen::Vector2d(0, 0), Eigen::Vector2d(1, 0), Eigen::Vector2d(0, 1), 0, followers, interior_edges); // order is important
1061 for (auto &s : followers)
1062 {
1063 faces[s.id].leader = f_id;
1064 face.followers.push_back(s.id);
1065 }
1066 for (int s : interior_edges)
1067 if (s >= 0 && edges[s].leader < 0 && edges[s].n_elem() > 0)
1068 edges[s].leader_face = f_id;
1069 }
1070 }
1072 {
1073 for (auto &vert : vertices)
1074 {
1075 vert.edge = -1;
1076 vert.face = -1;
1077 }
1078
1079 for (auto &small_edge : edges)
1080 {
1081 // invalid edges
1082 if (small_edge.n_elem() == 0)
1083 continue;
1084
1085 // not follower edges
1086 int large_edge = small_edge.leader;
1087 if (large_edge < 0)
1088 continue;
1089 assert(edges[large_edge].leader < 0);
1090
1091 // follower edges
1092 for (int j = 0; j < 2; j++)
1093 {
1094 const int v_id = small_edge.vertices(j);
1095 // hanging nodes
1096 if (v_id != edges[large_edge].vertices(0) && v_id != edges[large_edge].vertices(1))
1097 vertices[v_id].edge = large_edge;
1098 }
1099 }
1100
1101 for (auto &small_face : faces)
1102 {
1103 // invalid faces
1104 if (small_face.n_elem() == 0)
1105 continue;
1106
1107 // not follower faces
1108 int large_face = small_face.leader;
1109 if (large_face < 0)
1110 continue;
1111
1112 // follower faces
1113 for (int j = 0; j < 3; j++)
1114 {
1115 const int v_id = small_face.vertices(j);
1116 // hanging nodes
1117 if (v_id != faces[large_face].vertices(0) && v_id != faces[large_face].vertices(1) && v_id != faces[large_face].vertices(2))
1118 vertices[v_id].face = large_face;
1119 }
1120 }
1121 }
1122 std::array<int, 4> NCMesh3D::get_ordered_vertices_from_tet(const int element_index) const
1123 {
1124 return std::array<int, 4>({{cell_vertex(element_index, 0), cell_vertex(element_index, 1), cell_vertex(element_index, 2), cell_vertex(element_index, 3)}});
1125 }
1126 int NCMesh3D::add_element(Eigen::Vector4i v, int parent)
1127 {
1128 const int id = elements.size();
1129 const int level = (parent < 0) ? 0 : elements[parent].level + 1;
1130
1131 Eigen::Vector3d e1 = vertices[v[1]].pos - vertices[v[0]].pos;
1132 Eigen::Vector3d e2 = vertices[v[2]].pos - vertices[v[0]].pos;
1133 Eigen::Vector3d e3 = vertices[v[3]].pos - vertices[v[0]].pos;
1134 if ((e1.cross(e2)).dot(e3) < 0)
1135 std::swap(v[2], v[3]);
1136
1137 e2 = vertices[v[2]].pos - vertices[v[0]].pos;
1138 e3 = vertices[v[3]].pos - vertices[v[0]].pos;
1139 assert((e1.cross(e2)).dot(e3) > 0);
1140
1141 elements.emplace_back(3, v, level, parent);
1142
1143 if (parent >= 0)
1144 elements[id].body_id = elements[parent].body_id;
1145
1146 // add faces if not exist
1147 const int face012 = get_face(v[0], v[1], v[2]);
1148 const int face013 = get_face(v[0], v[1], v[3]);
1149 const int face123 = get_face(v[1], v[2], v[3]);
1150 const int face203 = get_face(v[2], v[0], v[3]);
1151
1152 faces[face012].add_element(id);
1153 faces[face013].add_element(id);
1154 faces[face123].add_element(id);
1155 faces[face203].add_element(id);
1156
1157 elements[id].faces << face012, face013, face123, face203;
1158
1159 // add edges if not exist
1160 const int edge01 = get_edge(v[0], v[1]);
1161 const int edge12 = get_edge(v[1], v[2]);
1162 const int edge20 = get_edge(v[2], v[0]);
1163 const int edge03 = get_edge(v[0], v[3]);
1164 const int edge13 = get_edge(v[1], v[3]);
1165 const int edge23 = get_edge(v[2], v[3]);
1166
1167 edges[edge01].add_element(id);
1168 edges[edge12].add_element(id);
1169 edges[edge20].add_element(id);
1170 edges[edge03].add_element(id);
1171 edges[edge13].add_element(id);
1172 edges[edge23].add_element(id);
1173
1174 elements[id].edges << edge01, edge12, edge20, edge03, edge13, edge23;
1175
1176 n_elements++;
1177
1178 for (int i = 0; i < v.size(); i++)
1179 vertices[v[i]].add_element(id);
1180
1181 return id;
1182 }
1183 } // namespace mesh
1184} // namespace polyfem
int V
int edge_id
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_vertex(Eigen::Vector2i v) const
Definition NCMesh3D.cpp:876
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition NCMesh3D.cpp:42
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 NCMesh3D.cpp:187
std::vector< int > all_to_valid_elemMap
Definition NCMesh3D.hpp:478
int n_faces() const override
number of faces
Definition NCMesh3D.hpp:189
int valid_to_all_edge(const int id) const
Definition NCMesh3D.hpp:402
Navigation3D::Index next_around_edge(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:516
int face_edge(const int f_id, const int le_id) const
Definition NCMesh3D.cpp:18
int n_cell_faces(const int c_id) const override
Definition NCMesh3D.hpp:219
void traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector< follower_edge > &list) const
Definition NCMesh3D.cpp:940
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > midpointMap
Definition NCMesh3D.hpp:474
int n_cells() const override
number of cells
Definition NCMesh3D.hpp:188
void build_element_vertex_adjacency()
void prepare_mesh() override
method used to finalize the mesh.
Definition NCMesh3D.hpp:344
std::vector< ncBoundary > edges
Definition NCMesh3D.hpp:471
int valid_to_all_face(const int id) const
Definition NCMesh3D.hpp:414
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition NCMesh3D.cpp:149
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition NCMesh3D.cpp:802
std::vector< int > refineHistory
Definition NCMesh3D.hpp:483
double edge_length(const int gid) const override
edge length
Definition NCMesh3D.cpp:107
int n_edges() const override
number of edges
Definition NCMesh3D.hpp:197
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 NCMesh3D.cpp:166
int get_edge(Eigen::Vector2i v)
Definition NCMesh3D.cpp:907
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition NCMesh3D.cpp:120
int add_element(Eigen::Vector4i v, int parent=-1)
int get_face(Eigen::Vector3i v)
Definition NCMesh3D.cpp:928
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition NCMesh3D.cpp:80
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition NCMesh3D.hpp:214
RowVectorNd cell_barycenter(const int c) const override
cell barycenter
Definition NCMesh3D.cpp:140
std::vector< uint32_t > edge_neighs(const int e_gid) const override
Definition NCMesh3D.cpp:446
int all_to_valid_elem(const int id) const
Definition NCMesh3D.hpp:421
int n_vertices() const override
number of vertices
Definition NCMesh3D.hpp:205
std::vector< ncVert > vertices
Definition NCMesh3D.hpp:470
Navigation3D::Index get_index_from_element_edge(int hi, int v0, int v1) const override
Definition NCMesh3D.cpp:357
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > edgeMap
Definition NCMesh3D.hpp:475
void traverse_face(int v1, int v2, int v3, Eigen::Vector2d p1, Eigen::Vector2d p2, Eigen::Vector2d p3, int depth, std::vector< follower_face > &face_list, std::vector< int > &edge_list) const
Definition NCMesh3D.cpp:996
void refine(const int n_refinement, const double t) override
refine the mesh
Definition NCMesh3D.cpp:26
void normalize() override
normalize the mesh
Definition NCMesh3D.cpp:87
std::vector< int > all_to_valid_vertexMap
Definition NCMesh3D.hpp:479
int all_to_valid_edge(const int id) const
Definition NCMesh3D.hpp:397
void get_face_elements_neighs(const int f_id, std::vector< int > &ids) const
Definition NCMesh3D.cpp:537
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
Definition NCMesh3D.cpp:827
int valid_to_all_elem(const int id) const
Definition NCMesh3D.hpp:426
void refine_element(int id_full)
Definition NCMesh3D.cpp:544
Navigation3D::Index switch_vertex(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:455
RowVectorNd kernel(const int cell_id) const override
Definition NCMesh3D.cpp:320
bool is_boundary_vertex(const int vertex_global_id) const override
is vertex boundary
Definition NCMesh3D.hpp:229
int cell_face(const int c_id, const int lf_id) const override
Definition NCMesh3D.hpp:221
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
Definition NCMesh3D.cpp:52
void compute_elements_tag() override
compute element types, see ElementType
Definition NCMesh3D.cpp:99
RowVectorNd point(const int global_index) const override
point coordinates
Definition NCMesh3D.cpp:115
std::vector< int > all_to_valid_edgeMap
Definition NCMesh3D.hpp:480
int find_edge(Eigen::Vector2i v) const
Definition NCMesh3D.cpp:898
Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const override
Definition NCMesh3D.cpp:325
bool is_boundary_face(const int face_global_id) const override
is face boundary
Definition NCMesh3D.hpp:231
int cell_vertex(const int f_id, const int lv_id) const override
id of the vertex of a cell
Definition NCMesh3D.hpp:220
std::vector< int > all_to_valid_faceMap
Definition NCMesh3D.hpp:481
void coarsen_element(int id_full)
Definition NCMesh3D.cpp:657
Navigation3D::Index next_around_face(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:520
std::vector< int > valid_to_all_elemMap
Definition NCMesh3D.hpp:478
std::vector< int > valid_to_all_edgeMap
Definition NCMesh3D.hpp:480
RowVectorNd face_barycenter(const int f) const override
face barycenter
Definition NCMesh3D.cpp:132
int n_cell_edges(const int c_id) const override
Definition NCMesh3D.hpp:218
Navigation3D::Index switch_element(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:497
void get_vertex_elements_neighs(const int v_id, std::vector< int > &ids) const override
Definition NCMesh3D.cpp:525
std::vector< ncElem > elements
Definition NCMesh3D.hpp:469
std::unordered_map< Eigen::Vector3i, int, ArrayHasher3D > faceMap
Definition NCMesh3D.hpp:476
RowVectorNd edge_barycenter(const int e) const override
edge barycenter
Definition NCMesh3D.cpp:125
std::vector< int > valid_to_all_faceMap
Definition NCMesh3D.hpp:481
std::vector< uint32_t > vertex_neighs(const int v_gid) const override
Definition NCMesh3D.cpp:438
std::vector< int > valid_to_all_vertexMap
Definition NCMesh3D.hpp:479
std::array< int, 4 > get_ordered_vertices_from_tet(const int element_index) const override
void update_elements_tag() override
Update elements types.
Definition NCMesh3D.cpp:103
int valid_to_all_vertex(const int id) const
Definition NCMesh3D.hpp:390
int edge_vertex(const int e_id, const int lv_id) const override
id of the edge vertex
Definition NCMesh3D.hpp:225
void refine_elements(const std::vector< int > &ids)
Definition NCMesh3D.cpp:647
int get_vertex(Eigen::Vector2i v)
Definition NCMesh3D.cpp:885
int find_face(Eigen::Vector3i v) const
Definition NCMesh3D.cpp:919
void get_edge_elements_neighs(const int e_id, std::vector< int > &ids) const override
Definition NCMesh3D.cpp:531
Navigation3D::Index switch_face(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:476
int face_vertex(const int f_id, const int lv_id) const override
id of the face vertex
Definition NCMesh3D.hpp:223
int all_to_valid_face(const int id) const
Definition NCMesh3D.hpp:409
bool load(const std::string &path) override
loads a mesh from the path
Definition NCMesh3D.cpp:832
Navigation3D::Index get_index_from_element_face(int hi, int v0, int v1, int v2) const override
Definition NCMesh3D.cpp:390
void set_body_ids(const std::vector< int > &body_ids) override
Set the volume sections.
Definition NCMesh3D.cpp:207
void set_boundary_ids(const std::vector< int > &boundary_ids) override
Set the boundary selection from a vector.
Definition NCMesh3D.cpp:199
std::vector< ncBoundary > faces
Definition NCMesh3D.hpp:472
Navigation3D::Index switch_edge(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:468
bool endswith(const std::string &str, const std::string &suffix)
std::tuple< bool, int, Tree > is_valid(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u, const double threshold)
Definition Jacobian.cpp:241
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13