PolyFEM
Loading...
Searching...
No Matches
CMesh3D.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
5
8#include <geogram/mesh/mesh.h>
9#include <iostream>
10#include <Eigen/Dense>
11#include <vector>
12#include <array>
13
14namespace polyfem
15{
16 namespace mesh
17 {
18 class CMesh3D : public Mesh3D
19 {
20 public:
21 CMesh3D() = default;
22 virtual ~CMesh3D() = default;
23 CMesh3D(CMesh3D &&) = default;
24 CMesh3D &operator=(CMesh3D &&) = default;
25 CMesh3D(const CMesh3D &) = default;
26 CMesh3D &operator=(const CMesh3D &) = default;
27
28 bool is_conforming() const override { return true; }
29
30 void refine(const int n_refinement, const double t) override;
31
32 int n_cells() const override { return int(mesh_.elements.size()); }
33 int n_faces() const override { return int(mesh_.faces.size()); }
34 int n_edges() const override { return int(mesh_.edges.size()); }
35 int n_vertices() const override { return int(mesh_.points.cols()); }
36
37 inline int n_face_vertices(const int f_id) const override { return mesh_.faces[f_id].vs.size(); }
38 inline int n_cell_vertices(const int c_id) const override { return mesh_.elements[c_id].vs.size(); }
39 inline int n_cell_edges(const int c_id) const override { return mesh_.elements[c_id].es.size(); }
40 inline int n_cell_faces(const int c_id) const override { return mesh_.elements[c_id].fs.size(); }
41 inline int cell_vertex(const int c_id, const int lv_id) const override { return mesh_.elements[c_id].vs[lv_id]; }
42 inline int cell_face(const int c_id, const int lf_id) const override { return mesh_.elements[c_id].fs[lf_id]; }
43 inline int cell_edge(const int c_id, const int le_id) const override { return mesh_.elements[c_id].es[le_id]; }
44 inline int face_vertex(const int f_id, const int lv_id) const override { return mesh_.faces[f_id].vs[lv_id]; }
45 inline int edge_vertex(const int e_id, const int lv_id) const override { return mesh_.edges[e_id].vs[lv_id]; }
46
47 void elements_boxes(std::vector<std::array<Eigen::Vector3d, 2>> &boxes) const override;
48 void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const override;
49
50 bool is_boundary_vertex(const int vertex_global_id) const override { return mesh_.vertices[vertex_global_id].boundary; }
51 bool is_boundary_edge(const int edge_global_id) const override { return mesh_.edges[edge_global_id].boundary; }
52 bool is_boundary_face(const int face_global_id) const override { return mesh_.faces[face_global_id].boundary; }
53 bool is_boundary_element(const int element_global_id) const override;
54
55 bool save(const std::string &path) const override;
56
57 bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override;
58
59 void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes) override;
60
61 void normalize() override;
62
63 double quad_area(const int gid) const override;
64
65 void compute_elements_tag() override;
66
67 RowVectorNd kernel(const int cell_id) const override;
68 RowVectorNd point(const int vertex_id) const override;
69 void set_point(const int global_index, const RowVectorNd &p) override;
70 RowVectorNd edge_barycenter(const int e) const override;
71 RowVectorNd face_barycenter(const int f) const override;
72 RowVectorNd cell_barycenter(const int c) const override;
73
74 void bounding_box(RowVectorNd &min, RowVectorNd &max) const override;
75
76 // navigation wrapper
77 Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const override { return Navigation3D::get_index_from_element_face(mesh_, hi, lf, lv); }
79
80 Navigation3D::Index get_index_from_element_edge(int hi, int v0, int v1) const override { return Navigation3D::get_index_from_element_edge(mesh_, hi, v0, v1); }
81 Navigation3D::Index get_index_from_element_face(int hi, int v0, int v1, int v2) const override { return Navigation3D::get_index_from_element_tri(mesh_, hi, v0, v1, v2); }
82
83 inline std::vector<uint32_t> vertex_neighs(const int v_gid) const override { return mesh_.vertices[v_gid].neighbor_hs; }
84 inline std::vector<uint32_t> edge_neighs(const int e_gid) const override { return mesh_.edges[e_gid].neighbor_hs; }
85
86 // Navigation in a surface mesh
91
92 // Iterate in a mesh
95
96 void get_vertex_elements_neighs(const int v_id, std::vector<int> &ids) const override
97 {
98 ids.clear();
99 ids.insert(ids.begin(), mesh_.vertices[v_id].neighbor_hs.begin(), mesh_.vertices[v_id].neighbor_hs.end());
100 }
101 void get_edge_elements_neighs(const int e_id, std::vector<int> &ids) const override
102 {
103 ids.clear();
104 ids.insert(ids.begin(), mesh_.edges[e_id].neighbor_hs.begin(), mesh_.edges[e_id].neighbor_hs.end());
105 }
106
107 void compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker) override;
108
109 void compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker) override;
110
111 // void triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const override;
112
113 // used for sweeping 2D mesh
115 {
116 std::cerr << "never user this function" << std::endl;
117 return mesh_;
118 }
119 static void geomesh_2_mesh_storage(const GEO::Mesh &gm, Mesh3DStorage &m);
120
121 void append(const Mesh &mesh) override;
122
123 std::unique_ptr<Mesh> copy() const override;
124
125 protected:
126 bool load(const std::string &path) override;
127 bool load(const GEO::Mesh &M) override;
128
129 private:
131 };
132 } // namespace mesh
133} // namespace polyfem
int V
void get_edge_elements_neighs(const int e_id, std::vector< int > &ids) const override
Definition CMesh3D.hpp:101
RowVectorNd cell_barycenter(const int c) const override
cell barycenter
Definition CMesh3D.cpp:1508
RowVectorNd edge_barycenter(const int e) const override
edge barycenter
Definition CMesh3D.cpp:1487
Navigation3D::Index get_index_from_element_face(int hi, int v0, int v1, int v2) const override
Definition CMesh3D.hpp:81
int n_cell_faces(const int c_id) const override
Definition CMesh3D.hpp:40
Navigation3D::Index get_index_from_element(int hi) const override
Definition CMesh3D.hpp:78
bool save(const std::string &path) const override
Definition CMesh3D.cpp:416
int cell_face(const int c_id, const int lf_id) const override
Definition CMesh3D.hpp:42
static void geomesh_2_mesh_storage(const GEO::Mesh &gm, Mesh3DStorage &m)
Definition CMesh3D.cpp:1522
RowVectorNd point(const int vertex_id) const override
point coordinates
Definition CMesh3D.cpp:1251
int cell_edge(const int c_id, const int le_id) const override
Definition CMesh3D.hpp:43
void normalize() override
normalize the mesh
Definition CMesh3D.cpp:1022
bool is_conforming() const override
if the mesh is conforming
Definition CMesh3D.hpp:28
int edge_vertex(const int e_id, const int lv_id) const override
id of the edge vertex
Definition CMesh3D.hpp:45
int n_cell_edges(const int c_id) const override
Definition CMesh3D.hpp:39
Navigation3D::Index switch_element(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:90
void refine(const int n_refinement, const double t) override
refine the mesh
Definition CMesh3D.cpp:19
int n_edges() const override
number of edges
Definition CMesh3D.hpp:34
int face_vertex(const int f_id, const int lv_id) const override
id of the face vertex
Definition CMesh3D.hpp:44
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 CMesh3D.cpp:1210
CMesh3D(CMesh3D &&)=default
RowVectorNd face_barycenter(const int f) const override
face barycenter
Definition CMesh3D.cpp:1494
CMesh3D & operator=(CMesh3D &&)=default
int n_vertices() const override
number of vertices
Definition CMesh3D.hpp:35
Mesh3DStorage mesh_
Definition CMesh3D.hpp:130
bool is_boundary_face(const int face_global_id) const override
is face boundary
Definition CMesh3D.hpp:52
CMesh3D & operator=(const CMesh3D &)=default
bool is_boundary_edge(const int edge_global_id) const override
is edge boundary
Definition CMesh3D.hpp:51
double quad_area(const int gid) const override
area of a quad face of an hex mesh
Definition CMesh3D.cpp:1466
std::vector< uint32_t > vertex_neighs(const int v_gid) const override
Definition CMesh3D.hpp:83
void get_vertex_elements_neighs(const int v_id, std::vector< int > &ids) const override
Definition CMesh3D.hpp:96
Navigation3D::Index next_around_face(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:94
bool is_boundary_vertex(const int vertex_global_id) const override
is vertex boundary
Definition CMesh3D.hpp:50
Mesh3DStorage & mesh_storge()
Definition CMesh3D.hpp:114
RowVectorNd kernel(const int cell_id) const override
Definition CMesh3D.cpp:1257
Navigation3D::Index next_around_edge(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:93
void elements_boxes(std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const override
constructs a box around every element (3d cell, 2d face)
Definition CMesh3D.cpp:1563
void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes) override
attach high order nodes
Definition CMesh3D.cpp:512
void set_point(const int global_index, const RowVectorNd &p) override
Set the point.
Definition CMesh3D.cpp:1240
void bounding_box(RowVectorNd &min, RowVectorNd &max) const override
computes the bbox of the mesh
Definition CMesh3D.cpp:1015
virtual ~CMesh3D()=default
Navigation3D::Index get_index_from_element_edge(int hi, int v0, int v1) const override
Definition CMesh3D.hpp:80
void append(const Mesh &mesh) override
appends a new mesh to the end of this
Definition CMesh3D.cpp:1599
int n_cells() const override
number of cells
Definition CMesh3D.hpp:32
int n_faces() const override
number of faces
Definition CMesh3D.hpp:33
Navigation3D::Index switch_vertex(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:87
bool load(const std::string &path) override
loads a mesh from the path
Definition CMesh3D.cpp:80
std::vector< uint32_t > edge_neighs(const int e_gid) const override
Definition CMesh3D.hpp:84
Navigation3D::Index get_index_from_element(int hi, int lf, int lv) const override
Definition CMesh3D.hpp:77
std::unique_ptr< Mesh > copy() const override
Create a copy of the mesh.
Definition CMesh3D.cpp:1611
void compute_elements_tag() override
compute element types, see ElementType
Definition CMesh3D.cpp:1264
Navigation3D::Index switch_edge(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:88
int n_cell_vertices(const int c_id) const override
number of vertices of a cell
Definition CMesh3D.hpp:38
Navigation3D::Index switch_face(Navigation3D::Index idx) const override
Definition CMesh3D.hpp:89
void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const override
constructs barycentric coodiantes for a point p.
Definition CMesh3D.cpp:1585
bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) override
build a mesh from matrices
Definition CMesh3D.cpp:469
bool is_boundary_element(const int element_global_id) const override
is cell boundary
Definition CMesh3D.cpp:1189
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition CMesh3D.hpp:37
int cell_vertex(const int c_id, const int lv_id) const override
id of the vertex of a cell
Definition CMesh3D.hpp:41
CMesh3D(const CMesh3D &)=default
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 CMesh3D.cpp:1228
std::vector< Element > elements
std::vector< Vertex > vertices
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
Index next_around_3Dedge(const Mesh3DStorage &M, Index idx)
Index switch_element(const Mesh3DStorage &M, Index idx)
Index get_index_from_element_tri(const Mesh3DStorage &M, int hi, int v0, int v1, int v2)
Index switch_vertex(const Mesh3DStorage &M, Index idx)
Index get_index_from_element_face(const Mesh3DStorage &M, int hi)
Index get_index_from_element_edge(const Mesh3DStorage &M, int hi, int v0, int v1)
Index switch_edge(const Mesh3DStorage &M, Index idx)
Index next_around_2Dface(const Mesh3DStorage &M, Index idx)
Index switch_face(const Mesh3DStorage &M, Index idx)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13