5#include <igl/writeOBJ.h>
25 if (n_refinement <= 0)
27 std::vector<bool> refine_mask(
elements.size(),
false);
28 for (
int i = 0; i <
elements.size(); i++)
30 refine_mask[i] =
true;
32 for (
int i = 0; i < refine_mask.size(); i++)
36 refine(n_refinement - 1, t);
42 mesh_.clear(
false,
false);
57 for (
int i = 0; i <
V.rows(); i++)
61 for (
int i = 0; i <
F.rows(); i++)
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!");
82 const double t = i / (n_new_nodes + 1.0);
84 return std::make_pair((1 - t) * v1 + t * v2, -1);
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;
98 return std::make_pair(b1 * v1 + b2 * v2 + b3 * v3, -1);
103 std::sort(v.data(), v.data() + v.size());
106 return search->second;
113 std::sort(v.data(), v.data() + v.size());
127 std::sort(v.data(), v.data() + v.size());
130 return search->second;
137 std::sort(v.data(), v.data() + v.size());
141 edges.emplace_back(v);
142 id =
edges.size() - 1;
151 const int level = (parent < 0) ? 0 :
elements[parent].level + 1;
152 elements.emplace_back(2, v, level, parent);
157 for (
int i = 0; i < v.size(); i++)
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]));
165 elements[id].edges << edge01, edge12, edge20;
167 edges[edge01].add_element(
id);
168 edges[edge12].add_element(
id);
169 edges[edge20].add_element(
id);
181 if (elem.is_not_valid())
182 throw std::runtime_error(
"Cannot refine an invalid element!");
184 const auto v = elem.vertices;
185 elem.is_refined =
true;
189 for (
int e = 0; e < 3; e++)
190 edges[elem.edges(e)].remove_element(id_full);
192 for (
int i = 0; i < v.size(); i++)
195 if (elem.children(0) >= 0)
197 for (
int c = 0; c < elem.children.size(); c++)
199 auto &child =
elements[elem.children(c)];
200 child.is_ghost =
false;
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++;
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]));
216 for (
int i = 0; i < v.size(); i++)
217 for (
int j = 0; j < i; j++)
229 add_element(Eigen::Vector3i(v[0], v01, v20), id_full);
231 add_element(Eigen::Vector3i(v[1], v12, v01), id_full);
233 add_element(Eigen::Vector3i(v[2], v20, v12), id_full);
235 add_element(Eigen::Vector3i(v12, v20, v01), id_full);
246 std::vector<int> full_ids(ids.size());
247 for (
int i = 0; i < ids.size(); i++)
250 for (
int i : full_ids)
256 const int parent_id =
elements[id_full].parent;
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!");
264 for (
int i = 0; i < parent.children.size(); i++)
266 auto &elem =
elements[parent.children(i)];
267 elem.is_ghost =
true;
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--;
276 parent.is_refined =
false;
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++;
291 for (
int i = 0; i <
vec.size(); i++)
301 for (
auto &edge :
edges)
304 edge.followers.clear();
305 edge.weights.setConstant(-1);
309 std::vector<follower_edge> followers;
310 for (
int e_id = 0; e_id <
elements.size(); e_id++)
312 const auto &element =
elements[e_id];
313 if (element.is_not_valid())
315 for (
int edge_local = 0; edge_local < 3; edge_local++)
317 v << element.vertices[edge_local], element.vertices[(edge_local + 1) % 3];
318 int edge_global = element.edges[edge_local];
319 assert(edge_global >= 0);
321 for (
auto &s : followers)
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;
334 for (
auto &edge :
edges)
336 if (edge.n_elem() == 1)
337 edge.isboundary =
true;
339 edge.isboundary =
false;
342 for (
auto &edge :
edges)
344 if (edge.leader >= 0)
346 edge.isboundary =
false;
347 edges[edge.leader].isboundary =
false;
352 vert.isboundary =
false;
354 for (
auto &edge :
edges)
356 if (edge.isboundary && edge.n_elem())
358 for (
int j = 0; j < 2; j++)
359 vertices[edge.vertices(j)].isboundary =
true;
364 double line_weight(Eigen::Matrix<double, 2, 2> &e, Eigen::VectorXd &v)
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)
383 Eigen::VectorXi vertexEdgeAdjacency;
384 vertexEdgeAdjacency.setConstant(
vertices.size(), 1, -1);
386 for (
int small_edge = 0; small_edge <
edges.size(); small_edge++)
388 if (
edges[small_edge].n_elem() == 0)
391 int large_edge =
edges[small_edge].leader;
395 int large_elem =
edges[large_edge].get_element();
396 for (
int j = 0; j < 2; j++)
398 int v_id =
edges[small_edge].vertices(j);
400 vertexEdgeAdjacency[v_id] = large_edge;
406 if (element.is_not_valid())
408 for (
int v_local = 0; v_local < 3; v_local++)
410 int v_global = element.vertices[v_local];
411 if (vertexEdgeAdjacency[v_global] < 0)
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];
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;
434 assert(fabs(pos(1)) < 1e-12);
438 assert(fabs(pos(0) + pos(1) - 1) < 1e-12);
442 assert(fabs(pos(0)) < 1e-12);
455 for (
int i = 0, e = 0; i <
elements.size(); i++)
469 for (
int i = 0, j = 0; i <
vertices.size(); i++)
481 for (
int i = 0, j = 0; i <
edges.size(); i++)
483 if (
edges[i].n_elem() == 0)
494 assert(
typeid(mesh) ==
typeid(
NCMesh2D));
507 for (
int i = 0; i < mesh2d.
n_faces(); i++)
509 Eigen::Vector3i face = mesh2d.
elements[i].vertices;
510 face = face.array() + n_v;
519 return std::make_unique<NCMesh2D>(*
this);
525 std::vector<follower_edge> list1, list2;
528 double p_mid = (p1 + p2) / 2;
529 traverse_edge(Eigen::Vector2i(v[0], v_mid), p1, p_mid, depth + 1, list1);
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);
537 std::make_move_iterator(list2.begin()),
538 std::make_move_iterator(list2.end()));
543 if (follower_id >= 0 &&
edges[follower_id].n_elem() > 0)
544 list.emplace_back(follower_id, p1, p2);
557 mesh_.clear(
false,
false);
561 Eigen::MatrixXd
V(mesh_.vertices.nb(), 2);
562 Eigen::MatrixXi
F(mesh_.facets.nb(), 3);
564 for (
int v = 0; v <
V.rows(); v++)
566 const double *ptr = mesh_.vertices.point_ptr(v);
567 V.row(v) << ptr[0], ptr[1];
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);
576 for (
int i = 0; i <
V.rows(); i++)
580 for (
int i = 0; i <
F.rows(); i++)
597 for (
int d = 0; d < 2; d++)
599 if (v.pos[d] > max[d])
601 if (v.pos[d] < min[d])
631 for (
int i = 0; i < edge.vertices.size(); i++)
632 if (edge.vertices(i) != v1)
634 v2 = edge.vertices(i);
653 for (
int i = 0; i < elem.edges.size(); i++)
655 const auto &edge =
edges[elem.edges(i)];
657 if (valid_edge_id != idx.
edge &&
find(edge.vertices, idx.
vertex) >= 0)
659 idx2.
edge = valid_edge_id;
675 if (edge.n_elem() == 2)
688 auto extent = max - min;
689 double scale = std::max(extent(0), extent(1));
692 v.pos = (v.pos - min.transpose()) / scale;
784 assert(body_ids.size() ==
n_faces());
785 for (
int i = 0; i < body_ids.size(); i++)
793 assert(boundary_ids.size() ==
n_edges());
794 for (
int i = 0; i < boundary_ids.size(); i++)
805 for (
int e = 0; e <
n_faces(); ++e)
817 for (
int e = 0; e <
n_edges(); ++e)
823 std::sort(vs.begin(), vs.end());
RowVectorNd face_barycenter(const int index) const override
face barycenter
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
std::vector< ElementType > elements_tag_
list of element types
std::vector< int > boundary_ids_
list of surface labels
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
std::vector< int > body_ids_
list of volume labels
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
int find_edge(Eigen::Vector2i v) const
bool is_boundary_element(const int element_global_id) const override
is cell boundary
void traverse_edge(Eigen::Vector2i v, double p1, double p2, int depth, std::vector< follower_edge > &list) const
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
Navigation::Index switch_face(Navigation::Index idx) const override
std::vector< int > valid_to_all_edgeMap
RowVectorNd point(const int global_index) const override
point coordinates
Navigation::Index switch_vertex(Navigation::Index idx) const override
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
int face_edge(const int f_id, const int le_id) const
void refine_elements(const std::vector< int > &ids)
void prepare_mesh() override
method used to finalize the mesh.
std::vector< ncVert > vertices
void build_edge_follower_chain()
std::vector< int > all_to_valid_elemMap
static double element_weight_to_edge_weight(const int l, const Eigen::Vector2d &pos)
void refine_element(int id_full)
void set_boundary_ids(const std::vector< int > &boundary_ids) override
Set the boundary selection from a vector.
int n_faces() const override
number of faces
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
int all_to_valid_vertex(const int id) const
std::vector< ncElem > elements
void set_body_ids(const std::vector< int > &body_ids) override
Set the volume sections.
Navigation::Index switch_edge(Navigation::Index idx) const override
int get_edge(Eigen::Vector2i v)
int add_element(Eigen::Vector3i v, int parent=-1)
void coarsen_element(int id_full)
int get_vertex(Eigen::Vector2i v)
Navigation::Index get_index_from_face(int f, int lv=0) const override
double edge_length(const int gid) const override
edge length
int valid_to_all_edge(const int id) const
int n_vertices() const override
number of vertices
bool is_boundary_edge(const int edge_global_id) const override
is edge boundary
void refine(const int n_refinement, const double t) override
refine the mesh
std::pair< RowVectorNd, int > face_node(const Navigation::Index &index, const int n_new_nodes, const int i, const int j) const override
void normalize() override
normalize the mesh
std::vector< int > refineHistory
int valid_to_all_elem(const int id) const
std::vector< int > valid_to_all_elemMap
int edge_vertex(const int e_id, const int lv_id) const override
id of the edge vertex
std::vector< int > all_to_valid_edgeMap
std::vector< int > all_to_valid_vertexMap
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > midpointMap
std::vector< ncBoundary > edges
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
void compute_elements_tag() override
compute element types, see ElementType
RowVectorNd edge_barycenter(const int index) const override
edge barycenter
int all_to_valid_edge(const int id) const
std::unordered_map< Eigen::Vector2i, int, ArrayHasher2D > edgeMap
int find_vertex(Eigen::Vector2i v) const
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
int n_face_vertices(const int f_id) const override
number of vertices of a face
int valid_to_all_vertex(const int id) const
void build_element_vertex_adjacency()
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
int n_edges() const override
number of edges
std::vector< int > valid_to_all_vertexMap
bool load(const std::string &path) override
loads a mesh from the path
std::pair< RowVectorNd, int > edge_node(const Navigation::Index &index, const int n_new_nodes, const int i) const override
void update_elements_tag() override
Update elements types.
void build_index_mapping()
void append(const Mesh &mesh) override
appends a new mesh to the end of this
int find(const Eigen::VectorXi &vec, int x)
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)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd