PolyFEM
Loading...
Searching...
No Matches
Mesh.cpp
Go to the documentation of this file.
1
7
11
14
15#include <geogram/mesh/mesh_io.h>
16#include <geogram/mesh/mesh_geometry.h>
17
18#include <Eigen/Geometry>
19
20#include <igl/boundary_facets.h>
21#include <igl/oriented_facets.h>
22#include <igl/edges.h>
23
24#include <filesystem>
25#include <unordered_set>
26
28namespace polyfem::mesh
29{
30 using namespace polyfem::io;
31 using namespace polyfem::utils;
32
33 namespace
34 {
35 std::vector<int> sort_face(const Eigen::RowVectorXi f)
36 {
37 std::vector<int> sorted_face(f.data(), f.data() + f.size());
38 std::sort(sorted_face.begin(), sorted_face.end());
39 return sorted_face;
40 }
41
42 // Constructs a list of unique faces represented in a given mesh (V,T)
43 //
44 // Inputs:
45 // T: #T × 4 matrix of indices of tet corners
46 // Outputs:
47 // F: #F × 3 list of faces in no particular order
48 template <typename DerivedT, typename DerivedF>
49 void get_faces(
50 const Eigen::MatrixBase<DerivedT> &T,
51 Eigen::PlainObjectBase<DerivedF> &F)
52 {
53 assert(T.cols() == 4);
54 assert(T.rows() >= 1);
55
56 Eigen::MatrixXi BF, OF;
57 igl::boundary_facets(T, BF);
58 igl::oriented_facets(T, OF); // boundary facets + duplicated interior faces
59 assert((OF.rows() + BF.rows()) % 2 == 0);
60 const int num_faces = (OF.rows() + BF.rows()) / 2;
61 F.resize(num_faces, 3);
62 F.topRows(BF.rows()) = BF;
63 std::unordered_set<std::vector<int>, HashVector> processed_faces;
64 for (int fi = 0; fi < BF.rows(); fi++)
65 {
66 processed_faces.insert(sort_face(BF.row(fi)));
67 }
68
69 for (int fi = 0; fi < OF.rows(); fi++)
70 {
71 std::vector<int> sorted_face = sort_face(OF.row(fi));
72 const auto iter = processed_faces.find(sorted_face);
73 if (iter == processed_faces.end())
74 {
75 F.row(processed_faces.size()) = OF.row(fi);
76 processed_faces.insert(sorted_face);
77 }
78 }
79
80 assert(F.rows() == processed_faces.size());
81 }
82 } // namespace
83
84 std::unique_ptr<Mesh> Mesh::create(const int dim, const bool non_conforming)
85 {
86 assert(dim == 2 || dim == 3);
87 if (dim == 2 && non_conforming)
88 return std::make_unique<NCMesh2D>();
89 else if (dim == 2 && !non_conforming)
90 return std::make_unique<CMesh2D>();
91 else if (dim == 3 && non_conforming)
92 return std::make_unique<NCMesh3D>();
93 else if (dim == 3 && !non_conforming)
94 return std::make_unique<CMesh3D>();
95 throw std::runtime_error("Invalid dimension");
96 }
97
98 std::unique_ptr<Mesh> Mesh::create(GEO::Mesh &meshin, const bool non_conforming)
99 {
100 if (is_planar(meshin))
101 {
102 generate_edges(meshin);
103 std::unique_ptr<Mesh> mesh = create(2, non_conforming);
104 if (mesh->load(meshin))
105 {
106 mesh->in_ordered_vertices_ = Eigen::VectorXi::LinSpaced(meshin.vertices.nb(), 0, meshin.vertices.nb() - 1);
107 assert(mesh->in_ordered_vertices_[0] == 0);
108 assert(mesh->in_ordered_vertices_[1] == 1);
109 assert(mesh->in_ordered_vertices_[2] == 2);
110 assert(mesh->in_ordered_vertices_[mesh->in_ordered_vertices_.size() - 1] == meshin.vertices.nb() - 1);
111
112 mesh->in_ordered_edges_.resize(meshin.edges.nb(), 2);
113
114 for (int e = 0; e < (int)meshin.edges.nb(); ++e)
115 {
116 for (int lv = 0; lv < 2; ++lv)
117 {
118 mesh->in_ordered_edges_(e, lv) = meshin.edges.vertex(e, lv);
119 }
120 assert(mesh->in_ordered_edges_(e, 0) != mesh->in_ordered_edges_(e, 1));
121 }
122 assert(mesh->in_ordered_edges_.size() > 0);
123
124 mesh->in_ordered_faces_.resize(0, 0);
125
126 return mesh;
127 }
128 }
129 else
130 {
131 std::unique_ptr<Mesh> mesh = create(3, non_conforming);
132 meshin.cells.connect();
133 if (mesh->load(meshin))
134 {
135 mesh->in_ordered_vertices_ = Eigen::VectorXi::LinSpaced(meshin.vertices.nb(), 0, meshin.vertices.nb() - 1);
136 assert(mesh->in_ordered_vertices_[0] == 0);
137 assert(mesh->in_ordered_vertices_[1] == 1);
138 assert(mesh->in_ordered_vertices_[2] == 2);
139 assert(mesh->in_ordered_vertices_[mesh->in_ordered_vertices_.size() - 1] == meshin.vertices.nb() - 1);
140
141 mesh->in_ordered_edges_.resize(meshin.edges.nb(), 2);
142
143 for (int e = 0; e < (int)meshin.edges.nb(); ++e)
144 {
145 for (int lv = 0; lv < 2; ++lv)
146 {
147 mesh->in_ordered_edges_(e, lv) = meshin.edges.vertex(e, lv);
148 }
149 }
150 assert(mesh->in_ordered_edges_.size() > 0);
151
152 mesh->in_ordered_faces_.resize(meshin.facets.nb(), meshin.facets.nb_vertices(0));
153
154 for (int f = 0; f < (int)meshin.edges.nb(); ++f)
155 {
156 assert(mesh->in_ordered_faces_.cols() == meshin.facets.nb_vertices(f));
157
158 for (int lv = 0; lv < mesh->in_ordered_faces_.cols(); ++lv)
159 {
160 mesh->in_ordered_faces_(f, lv) = meshin.facets.vertex(f, lv);
161 }
162 }
163 assert(mesh->in_ordered_faces_.size() > 0);
164
165 return mesh;
166 }
167 }
168
169 logger().error("Failed to load mesh");
170 return nullptr;
171 }
172
173 std::unique_ptr<Mesh> Mesh::create(const std::string &path, const bool non_conforming)
174 {
175 if (!std::filesystem::exists(path))
176 {
177 logger().error(path.empty() ? "No mesh provided!" : "Mesh file does not exist: {}", path);
178 return nullptr;
179 }
180
181 std::string lowername = path;
182 std::transform(lowername.begin(), lowername.end(), lowername.begin(), ::tolower);
183
184 if (StringUtils::endswith(lowername, ".hybrid"))
185 {
186 std::unique_ptr<Mesh> mesh = create(3, non_conforming);
187 if (mesh->load(path))
188 {
189 // TODO add in_ordered_vertices_, in_ordered_edges_, in_ordered_faces_
190 return mesh;
191 }
192 }
193 else if (StringUtils::endswith(lowername, ".msh"))
194 {
195 Eigen::MatrixXd vertices;
196 Eigen::MatrixXi cells;
197 std::vector<std::vector<int>> elements;
198 std::vector<std::vector<double>> weights;
199 std::vector<int> body_ids;
200
201 if (!MshReader::load(path, vertices, cells, elements, weights, body_ids))
202 {
203 logger().error("Failed to load MSH mesh: {}", path);
204 return nullptr;
205 }
206
207 const int dim = vertices.cols();
208 std::unique_ptr<Mesh> mesh = create(vertices, cells, non_conforming);
209
210 // Only tris and tets
211 if ((dim == 2 && cells.cols() == 3) || (dim == 3 && cells.cols() == 4))
212 {
213 mesh->attach_higher_order_nodes(vertices, elements);
214 mesh->set_cell_weights(weights);
215 // TODO: not clear?
216 }
217
218 for (const auto &w : weights)
219 {
220 if (!w.empty())
221 {
222 mesh->set_is_rational(true);
223 break;
224 }
225 }
226
227 mesh->set_body_ids(body_ids);
228
229 return mesh;
230 }
231 else
232 {
233 GEO::Mesh tmp;
234 if (GEO::mesh_load(path, tmp))
235 {
236 return create(tmp, non_conforming);
237 }
238 }
239 logger().error("Failed to load mesh: {}", path);
240 return nullptr;
241 }
242
243 std::unique_ptr<Mesh> Mesh::create(
244 const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &cells, const bool non_conforming)
245 {
246 const int dim = vertices.cols();
247
248 std::unique_ptr<Mesh> mesh = create(dim, non_conforming);
249
250 mesh->build_from_matrices(vertices, cells);
251
252 mesh->in_ordered_vertices_ = Eigen::VectorXi::LinSpaced(vertices.rows(), 0, vertices.rows() - 1);
253 assert(mesh->in_ordered_vertices_[0] == 0);
254 assert(mesh->in_ordered_vertices_[1] == 1);
255 assert(mesh->in_ordered_vertices_[2] == 2);
256 assert(mesh->in_ordered_vertices_[mesh->in_ordered_vertices_.size() - 1] == vertices.rows() - 1);
257
258 if (dim == 2)
259 {
260 std::unordered_set<std::pair<int, int>, HashPair> edges;
261 for (int f = 0; f < cells.rows(); ++f)
262 {
263 for (int lv = 0; lv < cells.cols(); ++lv)
264 {
265 const int v0 = cells(f, lv);
266 const int v1 = cells(f, (lv + 1) % cells.cols());
267 edges.emplace(std::pair<int, int>(std::min(v0, v1), std::max(v0, v1)));
268 }
269 }
270 mesh->in_ordered_edges_.resize(edges.size(), 2);
271 int index = 0;
272 for (auto it = edges.begin(); it != edges.end(); ++it)
273 {
274 mesh->in_ordered_edges_(index, 0) = it->first;
275 mesh->in_ordered_edges_(index, 1) = it->second;
276 ++index;
277 }
278
279 assert(mesh->in_ordered_edges_.size() > 0);
280
281 mesh->in_ordered_faces_.resize(0, 0);
282 }
283 else
284 {
285 if (cells.cols() == 4)
286 {
287 get_faces(cells, mesh->in_ordered_faces_);
288 igl::edges(mesh->in_ordered_faces_, mesh->in_ordered_edges_);
289 }
290 // else TODO
291 }
292
293 return mesh;
294 }
295
297
298 void Mesh::edge_barycenters(Eigen::MatrixXd &barycenters) const
299 {
300 barycenters.resize(n_edges(), dimension());
301 for (int e = 0; e < n_edges(); ++e)
302 {
303 barycenters.row(e) = edge_barycenter(e);
304 }
305 }
306
307 void Mesh::face_barycenters(Eigen::MatrixXd &barycenters) const
308 {
309 barycenters.resize(n_faces(), dimension());
310 for (int f = 0; f < n_faces(); ++f)
311 {
312 barycenters.row(f) = face_barycenter(f);
313 }
314 }
315
316 void Mesh::cell_barycenters(Eigen::MatrixXd &barycenters) const
317 {
318 barycenters.resize(n_cells(), dimension());
319 for (int c = 0; c < n_cells(); ++c)
320 {
321 barycenters.row(c) = cell_barycenter(c);
322 }
323 }
324
326
327 // Queries on the tags
328 bool Mesh::is_spline_compatible(const int el_id) const
329 {
330 if (is_volume())
331 {
334 // || elements_tag_[el_id] == ElementType::SIMPLE_SINGULAR_INTERIOR_CUBE
335 // || elements_tag_[el_id] == ElementType::SIMPLE_SINGULAR_BOUNDARY_CUBE;
336 }
337 else
338 {
341 // || elements_tag_[el_id] == ElementType::INTERFACE_CUBE
342 // || elements_tag_[el_id] == ElementType::SIMPLE_SINGULAR_INTERIOR_CUBE;
343 }
344 }
345
346 // -----------------------------------------------------------------------------
347
358
359 // -----------------------------------------------------------------------------
360
361 bool Mesh::is_polytope(const int el_id) const
362 {
365 }
366
367 void Mesh::update_nodes(const Eigen::VectorXi &in_node_to_node)
368 {
369 if (in_node_to_node.size() <= 0 || node_ids_.empty())
370 {
371 node_ids_.clear();
372 return;
373 }
374
375 const auto tmp = node_ids_;
376
377 for (int n = 0; n < n_vertices(); ++n)
378 {
379 node_ids_[in_node_to_node[n]] = tmp[n];
380 }
381 }
382
383 void Mesh::compute_node_ids(const std::function<int(const size_t, const RowVectorNd &, bool)> &marker)
384 {
385 node_ids_.resize(n_vertices());
386
387 for (int n = 0; n < n_vertices(); ++n)
388 {
389 bool is_boundary = is_boundary_vertex(n);
390 const auto p = point(n);
391 node_ids_[n] = marker(n, p, is_boundary);
392 }
393 }
394
395 void Mesh::load_boundary_ids(const std::string &path)
396 {
398
399 std::ifstream file(path);
400
401 std::string line;
402 int bindex = 0;
403 while (std::getline(file, line))
404 {
405 std::istringstream iss(line);
406 int v;
407 iss >> v;
408 boundary_ids_[bindex] = v;
409
410 ++bindex;
411 }
412
413 assert(boundary_ids_.size() == size_t(bindex));
414
415 file.close();
416 }
417
418 bool Mesh::is_simplex(const int el_id) const
419 {
420 return elements_tag_[el_id] == ElementType::SIMPLEX;
421 }
422
423 std::vector<std::pair<int, int>> Mesh::edges() const
424 {
425 std::vector<std::pair<int, int>> res;
426 res.reserve(n_edges());
427
428 for (int e_id = 0; e_id < n_edges(); ++e_id)
429 {
430 const int e0 = edge_vertex(e_id, 0);
431 const int e1 = edge_vertex(e_id, 1);
432
433 res.emplace_back(std::min(e0, e1), std::max(e0, e1));
434 }
435
436 return res;
437 }
438
439 std::vector<std::vector<int>> Mesh::faces() const
440 {
441 std::vector<std::vector<int>> res(n_faces());
442
443 for (int f_id = 0; f_id < n_faces(); ++f_id)
444 {
445 auto &tmp = res[f_id];
446 for (int lv_id = 0; lv_id < n_face_vertices(f_id); ++lv_id)
447 tmp.push_back(face_vertex(f_id, lv_id));
448
449 std::sort(tmp.begin(), tmp.end());
450 }
451
452 return res;
453 }
454
455 std::unordered_map<std::pair<int, int>, size_t, HashPair> Mesh::edges_to_ids() const
456 {
457 std::unordered_map<std::pair<int, int>, size_t, HashPair> res;
458 res.reserve(n_edges());
459
460 for (int e_id = 0; e_id < n_edges(); ++e_id)
461 {
462 const int e0 = edge_vertex(e_id, 0);
463 const int e1 = edge_vertex(e_id, 1);
464
465 res[std::pair<int, int>(std::min(e0, e1), std::max(e0, e1))] = e_id;
466 }
467
468 return res;
469 }
470
471 std::unordered_map<std::vector<int>, size_t, HashVector> Mesh::faces_to_ids() const
472 {
473 std::unordered_map<std::vector<int>, size_t, HashVector> res;
474 res.reserve(n_faces());
475
476 for (int f_id = 0; f_id < n_faces(); ++f_id)
477 {
478 std::vector<int> f;
479 f.reserve(n_face_vertices(f_id));
480 for (int lv_id = 0; lv_id < n_face_vertices(f_id); ++lv_id)
481 f.push_back(face_vertex(f_id, lv_id));
482 std::sort(f.begin(), f.end());
483
484 res[f] = f_id;
485 }
486
487 return res;
488 }
489
490 void Mesh::append(const Mesh &mesh)
491 {
492 const int n_vertices = this->n_vertices();
493
494 elements_tag_.insert(elements_tag_.end(), mesh.elements_tag_.begin(), mesh.elements_tag_.end());
495
496 // --------------------------------------------------------------------
497
498 // Initialize node_ids_ if it is not initialized yet.
499 if (!has_node_ids() && mesh.has_node_ids())
500 {
501 node_ids_.resize(n_vertices);
502 for (int i = 0; i < node_ids_.size(); ++i)
503 node_ids_[i] = get_node_id(i); // results in default if node_ids_ is empty
504 }
505
506 if (mesh.has_node_ids())
507 {
508 node_ids_.insert(node_ids_.end(), mesh.node_ids_.begin(), mesh.node_ids_.end());
509 }
510 else if (has_node_ids()) // && !mesh.has_node_ids()
511 {
512 node_ids_.resize(n_vertices + mesh.n_vertices());
513 for (int i = 0; i < mesh.n_vertices(); ++i)
514 node_ids_[n_vertices + i] = mesh.get_node_id(i); // results in default if node_ids_ is empty
515 }
516
517 assert(node_ids_.empty() || node_ids_.size() == n_vertices + mesh.n_vertices());
518
519 // --------------------------------------------------------------------
520
521 // Initialize boundary_ids_ if it is not initialized yet.
522 if (!has_boundary_ids() && mesh.has_boundary_ids())
523 {
525 for (int i = 0; i < boundary_ids_.size(); ++i)
527 }
528
529 if (mesh.has_boundary_ids())
530 {
531 boundary_ids_.insert(boundary_ids_.end(), mesh.boundary_ids_.begin(), mesh.boundary_ids_.end());
532 }
533 else if (has_boundary_ids()) // && !mesh.has_boundary_ids()
534 {
536 for (int i = 0; i < mesh.n_boundary_elements(); ++i)
537 boundary_ids_[n_boundary_elements() + i] = mesh.get_boundary_id(i); // results in default if mesh.boundary_ids_ is empty
538 }
539
540 // --------------------------------------------------------------------
541
542 // Initialize body_ids_ if it is not initialized yet.
543 if (!has_body_ids() && mesh.has_body_ids())
544 body_ids_ = std::vector<int>(n_elements(), 0); // 0 is the default body_id
545
546 if (mesh.has_body_ids())
547 body_ids_.insert(body_ids_.end(), mesh.body_ids_.begin(), mesh.body_ids_.end());
548 else if (has_body_ids()) // && !mesh.has_body_ids()
549 body_ids_.resize(n_elements() + mesh.n_elements(), 0); // 0 is the default body_id
550
551 // --------------------------------------------------------------------
552
553 if (orders_.size() == 0)
554 orders_.setOnes(n_elements(), 1);
555 Eigen::MatrixXi mesh_orders = mesh.orders_;
556 if (mesh_orders.size() == 0)
557 mesh_orders.setOnes(mesh.n_elements(), 1);
558 assert(orders_.cols() == mesh_orders.cols());
559 orders_.conservativeResize(orders_.rows() + mesh_orders.rows(), orders_.cols());
560 orders_.bottomRows(mesh_orders.rows()) = mesh_orders;
561
563
564 // --------------------------------------------------------------------
565 for (const auto &n : mesh.edge_nodes_)
566 {
567 auto tmp = n;
568 tmp.v1 += n_vertices;
569 tmp.v2 += n_vertices;
570 edge_nodes_.push_back(tmp);
571 }
572 for (const auto &n : mesh.face_nodes_)
573 {
574 auto tmp = n;
575 tmp.v1 += n_vertices;
576 tmp.v2 += n_vertices;
577 tmp.v3 += n_vertices;
578 face_nodes_.push_back(tmp);
579 }
580 for (const auto &n : mesh.cell_nodes_)
581 {
582 auto tmp = n;
583 tmp.v1 += n_vertices;
584 tmp.v2 += n_vertices;
585 tmp.v3 += n_vertices;
586 tmp.v4 += n_vertices;
587 cell_nodes_.push_back(tmp);
588 }
589 cell_weights_.insert(cell_weights_.end(), mesh.cell_weights_.begin(), mesh.cell_weights_.end());
590 // --------------------------------------------------------------------
591
592 assert(in_ordered_vertices_.cols() == mesh.in_ordered_vertices_.cols());
593 in_ordered_vertices_.conservativeResize(in_ordered_vertices_.rows() + mesh.in_ordered_vertices_.rows(), in_ordered_vertices_.cols());
594 in_ordered_vertices_.bottomRows(mesh.in_ordered_vertices_.rows()) = mesh.in_ordered_vertices_.array() + n_vertices;
595
596 if (in_ordered_edges_.size() == 0 || mesh.in_ordered_edges_.size() == 0)
597 in_ordered_edges_.resize(0, 0);
598 else
599 {
600 assert(in_ordered_edges_.cols() == mesh.in_ordered_edges_.cols());
602 }
603
604 if (in_ordered_faces_.size() == 0 || mesh.in_ordered_faces_.size() == 0)
605 in_ordered_faces_.resize(0, 0);
606 else
607 {
608 assert(in_ordered_faces_.cols() == mesh.in_ordered_faces_.cols());
610 }
611 }
612
613 namespace
614 {
615 template <typename T>
616 void transform_high_order_nodes(std::vector<T> &nodes, const MatrixNd &A, const VectorNd &b)
617 {
618 for (T &n : nodes)
619 {
620 if (n.nodes.size())
621 {
622 n.nodes = (n.nodes * A.transpose()).rowwise() + b.transpose();
623 }
624 }
625 }
626 } // namespace
627
629 {
630 for (int i = 0; i < n_vertices(); ++i)
631 {
632 VectorNd p = point(i).transpose();
633 p = A * p + b;
634 set_point(i, p.transpose());
635 }
636
637 transform_high_order_nodes(edge_nodes_, A, b);
638 transform_high_order_nodes(face_nodes_, A, b);
639 transform_high_order_nodes(cell_nodes_, A, b);
640 }
641} // namespace polyfem::mesh
static bool load(const std::string &path, Eigen::MatrixXd &vertices, Eigen::MatrixXi &cells, std::vector< std::vector< int > > &elements, std::vector< std::vector< double > > &weights, std::vector< int > &body_ids)
Definition MshReader.cpp:29
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:158
virtual int n_vertices() const =0
number of vertices
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:361
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:666
std::unordered_map< std::pair< int, int >, size_t, polyfem::utils::HashPair > edges_to_ids() const
map from edge (pair of v id) to the id of the edge
Definition Mesh.cpp:455
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:658
virtual RowVectorNd edge_barycenter(const int e) const =0
edge barycenter
bool is_rational_
stores if the mesh is rational
Definition Mesh.hpp:668
bool has_boundary_ids() const
checks if surface selections are available
Definition Mesh.hpp:524
void cell_barycenters(Eigen::MatrixXd &barycenters) const
all cells barycenters
Definition Mesh.cpp:316
virtual RowVectorNd face_barycenter(const int f) const =0
face barycenter
virtual void set_point(const int global_index, const RowVectorNd &p)=0
Set the point.
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:348
bool has_node_ids() const
checks if points selections are available
Definition Mesh.hpp:520
virtual RowVectorNd point(const int global_index) const =0
point coordinates
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:684
void compute_node_ids(const std::function< int(const size_t, const RowVectorNd &, bool)> &marker)
computes boundary selections based on a function
Definition Mesh.cpp:383
virtual int get_boundary_id(const int primitive) const
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:478
virtual bool is_boundary_vertex(const int vertex_global_id) const =0
is vertex boundary
bool is_simplex(const int el_id) const
checks if element is simples compatible
Definition Mesh.cpp:418
void face_barycenters(Eigen::MatrixXd &barycenters) const
all face barycenters
Definition Mesh.cpp:307
virtual bool has_body_ids() const
checks if volumes selections are available
Definition Mesh.hpp:528
bool is_spline_compatible(const int el_id) const
checks if element is spline compatible
Definition Mesh.cpp:328
virtual void load_boundary_ids(const std::string &path)
loads the boundary selections for a file
Definition Mesh.cpp:395
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:662
void apply_affine_transformation(const MatrixNd &A, const VectorNd &b)
Apply an affine transformation to the vertex positions .
Definition Mesh.cpp:628
std::vector< int > node_ids_
list of node labels
Definition Mesh.hpp:660
std::unordered_map< std::vector< int >, size_t, polyfem::utils::HashVector > faces_to_ids() const
map from face (tuple of v id) to the id of the face
Definition Mesh.cpp:471
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:675
std::vector< std::vector< double > > cell_weights_
weights associates to cells for rational polynomail meshes
Definition Mesh.hpp:677
virtual bool is_volume() const =0
checks if mesh is volume
std::vector< std::pair< int, int > > edges() const
list of sorted edges.
Definition Mesh.cpp:423
void update_nodes(const Eigen::VectorXi &in_node_to_node)
Update the node ids to reorder them.
Definition Mesh.cpp:367
virtual int edge_vertex(const int e_id, const int lv_id) const =0
id of the edge vertex
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:664
static std::unique_ptr< Mesh > create(const std::string &path, const bool non_conforming=false)
factory to build the proper mesh
Definition Mesh.cpp:173
virtual int get_default_boundary_id(const int primitive) const
Get the default boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:466
int dimension() const
utily for dimension
Definition Mesh.hpp:148
virtual int n_cells() const =0
number of cells
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:673
virtual int n_faces() const =0
number of faces
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:671
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:439
int n_boundary_elements() const
utitlity to return the number of boundary elements, faces or edges in 3d and 2d
Definition Mesh.hpp:163
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:682
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:490
virtual RowVectorNd cell_barycenter(const int c) const =0
cell barycenter
virtual int n_edges() const =0
number of edges
void edge_barycenters(Eigen::MatrixXd &barycenters) const
all edges barycenters
Definition Mesh.cpp:298
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:680
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:487
virtual int face_vertex(const int f_id, const int lv_id) const =0
id of the face vertex
bool is_planar(const GEO::Mesh &M, const double tol=1e-5)
Determine if the given mesh is planar (2D or tiny z-range).
Definition MeshUtils.cpp:31
@ REGULAR_INTERIOR_CUBE
Triangle/tet element.
@ REGULAR_BOUNDARY_CUBE
Quad/Hex incident to more than 1 singular vertices (should not happen in 2D)
@ MULTI_SINGULAR_BOUNDARY_CUBE
Quad incident to exactly 1 singular vertex (in 2D); hex incident to exactly 1 singular interior edge,...
@ MULTI_SINGULAR_INTERIOR_CUBE
Quad/hex incident to exactly 1 singular vertex (in 2D) or edge (in 3D)
@ SIMPLE_SINGULAR_INTERIOR_CUBE
Regular quad/hex inside a 3^n patch.
@ INTERFACE_CUBE
Boundary hex that is not regular nor SimpleSingularBoundaryCube.
@ INTERIOR_POLYTOPE
Quad/hex that is at the interface with a polytope (if a cube has both external boundary and and inter...
@ SIMPLE_SINGULAR_BOUNDARY_CUBE
Boundary quad/hex, where all boundary vertices/edges are incident to at most 2 quads/hexes.
@ BOUNDARY_POLYTOPE
Interior polytope.
void generate_edges(GEO::Mesh &M)
assing edges to M
bool endswith(const std::string &str, const std::string &suffix)
void append_rows(DstMat &dst, const SrcMat &src)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:9
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor, 3, 3 > MatrixNd
Definition Types.hpp:12
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:11