PolyFEM
Loading...
Searching...
No Matches
NCMesh3D.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <set>
5#include <tuple>
6
7namespace polyfem
8{
9 namespace mesh
10 {
11 class NCMesh3D : public Mesh3D
12 {
13 public:
15 {
16 int id;
17 double p1, p2;
18
19 follower_edge(int id_, double p1_, double p2_)
20 {
21 id = id_;
22 p1 = p1_;
23 p2 = p2_;
24 }
25 };
26
28 {
29 int id;
30 Eigen::Vector2d p1, p2, p3;
31
32 follower_face(int id_, Eigen::Vector2d p1_, Eigen::Vector2d p2_, Eigen::Vector2d p3_)
33 {
34 id = id_;
35 p1 = p1_;
36 p2 = p2_;
37 p3 = p3_;
38 }
39 };
40
41 struct ncVert
42 {
43 ncVert(const Eigen::VectorXd pos_) : pos(pos_) {};
44 ~ncVert() {};
45
46 int n_elem() const
47 {
48 return elem_list.size();
49 }
50
51 void add_element(const int e)
52 {
53 elem_list.insert(e);
54 };
55
56 void remove_element(const int e)
57 {
58 int num = elem_list.erase(e);
59 assert(num == 1);
60 }
61
62 std::set<int> elem_list; // elements that contain this edge/face
63
64 Eigen::VectorXd pos;
65 bool isboundary = false;
66
67 int edge = -1; // only used if the vertex is on the interior of an edge
68 int face = -1; // only used if the vertex is on the interior of an face
69
70 double weight = -1.; // only 2d, the local position of this vertex on the edge
71 };
72
74 {
75 ncBoundary(const Eigen::VectorXi vert) : vertices(vert)
76 {
77 weights.setConstant(-1);
78 };
80
81 int n_elem() const
82 {
83 return elem_list.size();
84 }
85
86 void add_element(const int e)
87 {
88 elem_list.insert(e);
89 };
90
91 void remove_element(const int e)
92 {
93 int num = elem_list.erase(e);
94 assert(num == 1);
95 }
96
97 int get_element() const
98 {
99 assert(n_elem() > 0);
100 auto it = elem_list.cbegin();
101 return *it;
102 }
103
104 int find_opposite_element(int e) const
105 {
106 assert(n_elem() == 2);
107 bool exist = false;
108 int oppo = -1;
109 for (int elem : elem_list)
110 if (elem == e)
111 exist = true;
112 else
113 oppo = elem;
114
115 assert(oppo >= 0);
116 assert(exist);
117 return oppo;
118 }
119
120 std::set<int> elem_list; // elements that contain this edge/face
121 Eigen::VectorXi vertices;
122
123 bool isboundary = false; // valid only after calling mark_boundary()
124 int boundary_id = 0;
125
126 int leader = -1; // if this edge/face lies on a larger edge/face
127 std::vector<int> followers; // followers of this edge/face
128
129 int leader_face = -1; // if this edge lies in the interior of a face
130
131 std::vector<int> global_ids; // only used for building basis
132
133 // the following only used if it's an edge
134 Eigen::Vector2d weights; // position of this edge on its leader edge
135 };
136
137 struct ncElem
138 {
139 ncElem(const int dim_, const Eigen::VectorXi vertices_, const int level_, const int parent_) : dim(dim_), level(level_), parent(parent_), geom_vertices(vertices_)
140 {
142
143 assert(geom_vertices.size() == dim + 1);
144 edges.setConstant(3 * (dim - 1), 1, -1);
145 faces.setConstant(4 * (dim - 2), 1, -1);
146 children.setConstant(std::round(pow(2, dim)), 1, -1);
147 };
149
150 bool is_valid() const
151 {
152 return (!is_ghost) && (!is_refined);
153 };
154
155 bool is_not_valid() const
156 {
157 return is_ghost || is_refined;
158 };
159
160 int dim;
161 int level; // level of refinement
163 Eigen::VectorXi geom_vertices;
164
165 Eigen::VectorXi vertices; // geom_vertices with a different order that used in polyfem
166
167 Eigen::VectorXi edges;
168 Eigen::VectorXi faces;
169 Eigen::VectorXi children;
170
172
173 bool is_refined = false;
174 bool is_ghost = false;
175 };
176
177 NCMesh3D() = default;
178 virtual ~NCMesh3D() = default;
179 NCMesh3D(NCMesh3D &&) = default;
181 NCMesh3D(const NCMesh3D &) = default;
182 NCMesh3D &operator=(const NCMesh3D &) = default;
183
184 void refine(const int n_refinement, const double t) override;
185
186 bool is_conforming() const override { return false; }
187
188 int n_cells() const override { return n_elements; }
189 int n_faces() const override
190 {
191 int n = 0;
192 for (const auto &face : faces)
193 if (face.n_elem())
194 n++;
195 return n;
196 }
197 int n_edges() const override
198 {
199 int n = 0;
200 for (const auto &edge : edges)
201 if (edge.n_elem())
202 n++;
203 return n;
204 }
205 int n_vertices() const override
206 {
207 int n = 0;
208 for (const auto &vert : vertices)
209 if (vert.n_elem())
210 n++;
211 return n;
212 }
213
214 int n_face_vertices(const int f_id) const override { return 3; }
215 int n_face_cells(const int f_id) const { return faces[valid_to_all_face(f_id)].n_elem(); }
216 int n_edge_cells(const int e_id) const { return edges[valid_to_all_edge(e_id)].n_elem(); }
217 int n_cell_vertices(const int c_id) const override { return 4; }
218 int n_cell_edges(const int c_id) const override { return 6; }
219 int n_cell_faces(const int c_id) const override { return 4; }
220 int cell_vertex(const int f_id, const int lv_id) const override { return all_to_valid_vertex(elements[valid_to_all_elem(f_id)].vertices(lv_id)); }
221 int cell_face(const int c_id, const int lf_id) const override { return all_to_valid_face(elements[valid_to_all_elem(c_id)].faces(lf_id)); }
222 int cell_edge(const int c_id, const int le_id) const override { return all_to_valid_edge(elements[valid_to_all_elem(c_id)].edges(le_id)); }
223 int face_vertex(const int f_id, const int lv_id) const override { return all_to_valid_vertex(faces[valid_to_all_face(f_id)].vertices(lv_id)); }
224 int face_edge(const int f_id, const int le_id) const;
225 int edge_vertex(const int e_id, const int lv_id) const override { return all_to_valid_vertex(edges[valid_to_all_edge(e_id)].vertices(lv_id)); }
226
227 inline int cell_ref_level(const int c_id) const { return elements[valid_to_all_elem(c_id)].level; }
228
229 bool is_boundary_vertex(const int vertex_global_id) const override { return vertices[valid_to_all_vertex(vertex_global_id)].isboundary; }
230 bool is_boundary_edge(const int edge_global_id) const override { return edges[valid_to_all_edge(edge_global_id)].isboundary; }
231 bool is_boundary_face(const int face_global_id) const override { return faces[valid_to_all_face(face_global_id)].isboundary; }
232 bool is_boundary_element(const int element_global_id) const override;
233
234 bool save(const std::string &path) const override
235 {
236 // TODO
237 return false;
238 }
239
240 bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override;
241
242 void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes) override;
243
244 void normalize() override;
245
246 void compute_elements_tag() override;
247 void update_elements_tag() override;
248
249 double edge_length(const int gid) const override;
250
251 RowVectorNd point(const int global_index) const override;
252 void set_point(const int global_index, const RowVectorNd &p) override;
253 RowVectorNd edge_barycenter(const int e) const override;
254 RowVectorNd face_barycenter(const int f) const override;
255 RowVectorNd cell_barycenter(const int c) const override;
256
257 void bounding_box(RowVectorNd &min, RowVectorNd &max) const override;
258
259 void compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker) override;
260
261 void compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker) override;
262 void set_boundary_ids(const std::vector<int> &boundary_ids) override;
263 void set_body_ids(const std::vector<int> &body_ids) override;
264
265 int get_boundary_id(const int primitive) const override { return faces[valid_to_all_face(primitive)].boundary_id; };
266 int get_body_id(const int primitive) const override { return elements[valid_to_all_elem(primitive)].body_id; };
267
268 // void triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const override;
269
270 RowVectorNd kernel(const int cell_id) const override;
271 // navigation wrapper
272 Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const override;
273 Navigation3D::Index get_index_from_element(int hi) const override;
274
275 Navigation3D::Index get_index_from_element_edge(int hi, int v0, int v1) const override;
276 Navigation3D::Index get_index_from_element_face(int hi, int v0, int v1, int v2) const override;
277
278 std::vector<uint32_t> vertex_neighs(const int v_gid) const override;
279 std::vector<uint32_t> edge_neighs(const int e_gid) const override;
280
281 int leader_edge_of_vertex(const int v_id) const
282 {
283 assert(adj_prepared);
284 return (vertices[valid_to_all_vertex(v_id)].edge < 0) ? -1 : all_to_valid_edge(vertices[valid_to_all_vertex(v_id)].edge);
285 }
286 int leader_edge_of_edge(const int e_id) const
287 {
288 assert(adj_prepared);
289 return (edges[valid_to_all_edge(e_id)].leader < 0) ? -1 : all_to_valid_edge(edges[valid_to_all_edge(e_id)].leader);
290 }
291
292 int leader_face_of_vertex(const int v_id) const
293 {
294 assert(adj_prepared);
295 int id_full = vertices[valid_to_all_vertex(v_id)].face;
296 return (id_full < 0) ? -1 : all_to_valid_face(id_full);
297 }
298 int leader_face_of_edge(const int e_id) const
299 {
300 assert(adj_prepared);
301 int id_full = edges[valid_to_all_edge(e_id)].leader_face;
302 return (id_full < 0) ? -1 : all_to_valid_face(id_full);
303 }
304 int leader_face_of_face(const int f_id) const
305 {
306 assert(adj_prepared);
307 int id_full = faces[valid_to_all_face(f_id)].leader;
308 return (id_full < 0) ? -1 : all_to_valid_face(id_full);
309 }
310
311 // number of follower edges of a leader edge
312 int n_follower_edges(const int e_id) const
313 {
314 assert(adj_prepared);
315 return edges[valid_to_all_edge(e_id)].followers.size();
316 }
317 int n_follower_faces(const int e_id) const
318 {
319 assert(adj_prepared);
320 return faces[valid_to_all_face(e_id)].followers.size();
321 }
322
323 // Navigation in a surface mesh
328
329 // Iterate in a mesh
332
333 void get_vertex_elements_neighs(const int v_id, std::vector<int> &ids) const override;
334 void get_edge_elements_neighs(const int e_id, std::vector<int> &ids) const override;
335 void get_face_elements_neighs(const int f_id, std::vector<int> &ids) const;
336
337 void refine_element(int id_full);
338 void refine_elements(const std::vector<int> &ids);
339
340 void coarsen_element(int id_full);
341
342 void mark_boundary();
343
354
355 void build_index_mapping();
356
357 std::array<int, 4> get_ordered_vertices_from_tet(const int element_index) const override;
358
359 void append(const Mesh &mesh) override;
360
361 std::unique_ptr<Mesh> copy() const override;
362
363 private:
365 {
366 long operator()(const Eigen::Vector2i &a) const
367 {
368 return (long)((long)984120265 * a[0] + (long)125965121 * a[1]);
369 }
370 };
371
373 {
374 long operator()(const Eigen::Vector3i &a) const
375 {
376 return (long)((long)984120265 * a[0] + (long)125965121 * a[1] + (long)495698413 * a[2]);
377 }
378 };
379
380 protected:
381 bool load(const std::string &path) override;
382 bool load(const GEO::Mesh &M) override;
383
384 // index map from vertices to valid ones, and its inverse
385 inline int all_to_valid_vertex(const int id) const
386 {
387 assert(index_prepared);
388 return all_to_valid_vertexMap[id];
389 };
390 inline int valid_to_all_vertex(const int id) const
391 {
392 assert(index_prepared);
393 return valid_to_all_vertexMap[id];
394 };
395
396 // index map from edges to valid ones, and its inverse
397 inline int all_to_valid_edge(const int id) const
398 {
399 assert(index_prepared);
400 return all_to_valid_edgeMap[id];
401 };
402 inline int valid_to_all_edge(const int id) const
403 {
404 assert(index_prepared);
405 return valid_to_all_edgeMap[id];
406 };
407
408 // index map from faces to valid ones, and its inverse
409 inline int all_to_valid_face(const int id) const
410 {
411 assert(index_prepared);
412 return all_to_valid_faceMap[id];
413 };
414 inline int valid_to_all_face(const int id) const
415 {
416 assert(index_prepared);
417 return valid_to_all_faceMap[id];
418 };
419
420 // index map from elements to valid ones, and its inverse
421 inline int all_to_valid_elem(const int id) const
422 {
423 assert(index_prepared);
424 return all_to_valid_elemMap[id];
425 };
426 inline int valid_to_all_elem(const int id) const
427 {
428 assert(index_prepared);
429 return valid_to_all_elemMap[id];
430 };
431
432 // find the mid-point of edge v[0]v[1], return -1 if not exists
433 int find_vertex(Eigen::Vector2i v) const;
434 int find_vertex(const int v1, const int v2) const { return find_vertex(Eigen::Vector2i(v1, v2)); };
435
436 // find the mid-point of edge v[0]v[1], create one if not exists
437 int get_vertex(Eigen::Vector2i v);
438
439 // find the edge v[0]v[1], return -1 if not exists
440 int find_edge(Eigen::Vector2i v) const;
441 int find_edge(const int v1, const int v2) const { return find_edge(Eigen::Vector2i(v1, v2)); };
442
443 // find the edge v[0]v[1], create one if not exists
444 int get_edge(Eigen::Vector2i v);
445 int get_edge(const int v1, const int v2) { return get_edge(Eigen::Vector2i(v1, v2)); };
446
447 // find the face, return -1 if not exists
448 int find_face(Eigen::Vector3i v) const;
449 int find_face(const int v1, const int v2, const int v3) const { return find_face(Eigen::Vector3i(v1, v2, v3)); };
450 // find the face, create one if not exists
451 int get_face(Eigen::Vector3i v);
452 int get_face(const int v1, const int v2, const int v3) { return get_face(Eigen::Vector3i(v1, v2, v3)); };
453
454 void traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector<follower_edge> &list) const;
456
457 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;
459
461
462 int add_element(Eigen::Vector4i v, int parent = -1);
463
464 int n_elements = 0;
465
466 bool index_prepared = false;
467 bool adj_prepared = false;
468
469 std::vector<ncElem> elements;
470 std::vector<ncVert> vertices;
471 std::vector<ncBoundary> edges;
472 std::vector<ncBoundary> faces;
473
474 std::unordered_map<Eigen::Vector2i, int, ArrayHasher2D> midpointMap;
475 std::unordered_map<Eigen::Vector2i, int, ArrayHasher2D> edgeMap;
476 std::unordered_map<Eigen::Vector3i, int, ArrayHasher3D> faceMap;
477
482
483 std::vector<int> refineHistory;
484
485 // elementAdj(i, j) = 1 iff element i touches element j
486 Eigen::SparseMatrix<bool, Eigen::RowMajor> elementAdj;
487 };
488 } // namespace mesh
489} // namespace polyfem
int V
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int all_to_valid_vertex(const int id) const
Definition NCMesh3D.hpp:385
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 leader_face_of_vertex(const int v_id) const
Definition NCMesh3D.hpp:292
int find_edge(const int v1, const int v2) const
Definition NCMesh3D.hpp:441
int cell_ref_level(const int c_id) const
Definition NCMesh3D.hpp:227
int n_faces() const override
number of faces
Definition NCMesh3D.hpp:189
int n_cell_vertices(const int c_id) const override
number of vertices of a cell
Definition NCMesh3D.hpp:217
NCMesh3D & operator=(const NCMesh3D &)=default
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
int n_edge_cells(const int e_id) const
Definition NCMesh3D.hpp:216
int get_boundary_id(const int primitive) const override
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition NCMesh3D.hpp:265
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
int n_follower_faces(const int e_id) const
Definition NCMesh3D.hpp:317
bool save(const std::string &path) const override
Definition NCMesh3D.hpp:234
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 n_follower_edges(const int e_id) const
Definition NCMesh3D.hpp:312
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
int leader_edge_of_vertex(const int v_id) const
Definition NCMesh3D.hpp:281
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
int find_face(const int v1, const int v2, const int v3) const
Definition NCMesh3D.hpp:449
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
int find_vertex(const int v1, const int v2) const
Definition NCMesh3D.hpp:434
virtual ~NCMesh3D()=default
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
int leader_face_of_edge(const int e_id) const
Definition NCMesh3D.hpp:298
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
NCMesh3D & operator=(NCMesh3D &&)=default
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
int get_face(const int v1, const int v2, const int v3)
Definition NCMesh3D.hpp:452
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
int leader_edge_of_edge(const int e_id) const
Definition NCMesh3D.hpp:286
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
int get_body_id(const int primitive) const override
Get the volume selection of an element (cell in 3d, face in 2d)
Definition NCMesh3D.hpp:266
bool is_boundary_edge(const int edge_global_id) const override
is edge boundary
Definition NCMesh3D.hpp:230
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
NCMesh3D(const NCMesh3D &)=default
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 get_edge(const int v1, const int v2)
Definition NCMesh3D.hpp:445
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
int cell_edge(const int c_id, const int le_id) const override
Definition NCMesh3D.hpp:222
Navigation3D::Index switch_face(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:476
int leader_face_of_face(const int f_id) const
Definition NCMesh3D.hpp:304
Eigen::SparseMatrix< bool, Eigen::RowMajor > elementAdj
Definition NCMesh3D.hpp:486
bool is_conforming() const override
if the mesh is conforming
Definition NCMesh3D.hpp:186
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
int n_face_cells(const int f_id) const
Definition NCMesh3D.hpp:215
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
NCMesh3D(NCMesh3D &&)=default
std::vector< ncBoundary > faces
Definition NCMesh3D.hpp:472
Navigation3D::Index switch_edge(Navigation3D::Index idx) const override
Definition NCMesh3D.cpp:468
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
long operator()(const Eigen::Vector2i &a) const
Definition NCMesh3D.hpp:366
long operator()(const Eigen::Vector3i &a) const
Definition NCMesh3D.hpp:374
follower_edge(int id_, double p1_, double p2_)
Definition NCMesh3D.hpp:19
follower_face(int id_, Eigen::Vector2d p1_, Eigen::Vector2d p2_, Eigen::Vector2d p3_)
Definition NCMesh3D.hpp:32
ncBoundary(const Eigen::VectorXi vert)
Definition NCMesh3D.hpp:75
ncElem(const int dim_, const Eigen::VectorXi vertices_, const int level_, const int parent_)
Definition NCMesh3D.hpp:139
void remove_element(const int e)
Definition NCMesh3D.hpp:56
void add_element(const int e)
Definition NCMesh3D.hpp:51
ncVert(const Eigen::VectorXd pos_)
Definition NCMesh3D.hpp:43