PolyFEM
Loading...
Searching...
No Matches
Mesh.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
8
9#include <Eigen/Dense>
10#include <geogram/mesh/mesh.h>
11
12#include <memory>
13
14namespace polyfem
15{
16 namespace mesh
17 {
38
40 class Mesh
41 {
42 protected:
45 {
46 public:
47 int v1, v2;
48 std::vector<int> nodes_ids;
49 Eigen::MatrixXd nodes;
50 };
51
54 {
55 public:
56 int v1, v2, v3;
57 std::vector<int> nodes_ids;
58 Eigen::MatrixXd nodes;
59 };
60
63 {
64 public:
65 int v1, v2, v3, v4;
66 std::vector<int> nodes_ids;
67 Eigen::MatrixXd nodes;
68 };
69
70 public:
77 static std::unique_ptr<Mesh> create(const std::string &path, const bool non_conforming = false);
78
85 static std::unique_ptr<Mesh> create(GEO::Mesh &M, const bool non_conforming = false);
86
94 static std::unique_ptr<Mesh> create(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &cells, const bool non_conforming = false);
95
102 static std::unique_ptr<Mesh> create(const int dim, const bool non_conforming = false);
103
106 virtual std::unique_ptr<Mesh> copy() const = 0;
107
108 protected:
112 Mesh() = default;
113
114 public:
117 virtual ~Mesh() = default;
121 Mesh(Mesh &&) = default;
126 Mesh &operator=(Mesh &&) = default;
130 Mesh(const Mesh &) = default;
135 Mesh &operator=(const Mesh &) = default;
136
142 virtual void refine(const int n_refinement, const double t) = 0;
143
148 virtual bool is_volume() const = 0;
153 int dimension() const { return (is_volume() ? 3 : 2); }
158 virtual bool is_conforming() const = 0;
163 int n_elements() const { return (is_volume() ? n_cells() : n_faces()); }
168 int n_boundary_elements() const { return (is_volume() ? n_faces() : n_edges()); }
169
174 virtual int n_cells() const = 0;
179 virtual int n_faces() const = 0;
183 virtual int n_edges() const = 0;
187 virtual int n_vertices() const = 0;
188
193 virtual int n_face_vertices(const int f_id) const = 0;
198 virtual int n_cell_vertices(const int c_id) const = 0;
204 virtual int edge_vertex(const int e_id, const int lv_id) const = 0;
210 virtual int face_vertex(const int f_id, const int lv_id) const = 0;
216 virtual int cell_vertex(const int f_id, const int lv_id) const = 0;
222 int element_vertex(const int el_id, const int lv_id) const
223 {
224 return (is_volume() ? cell_vertex(el_id, lv_id) : face_vertex(el_id, lv_id));
225 }
226
231 std::vector<int> element_vertices(const int el_id) const
232 {
233 std::vector<int> res;
234 for (int i = 0; i < n_cell_vertices(el_id); ++i)
235 res.push_back(element_vertex(el_id, i));
236 std::sort(res.begin(), res.end());
237 return res;
238 }
239
240 int boundary_element_vertex(const int primitive_id, const int lv_id) const
241 {
242 return (is_volume() ? face_vertex(primitive_id, lv_id) : edge_vertex(primitive_id, lv_id));
243 }
244
249 virtual bool is_boundary_vertex(const int vertex_global_id) const = 0;
254 virtual bool is_boundary_edge(const int edge_global_id) const = 0;
259 virtual bool is_boundary_face(const int face_global_id) const = 0;
264 virtual bool is_boundary_element(const int element_global_id) const = 0;
265
266 virtual bool save(const std::string &path) const = 0;
267
268 private:
274 virtual bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) = 0;
275
276 public:
281 virtual void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes) = 0;
285 inline const Eigen::MatrixXi &orders() const { return orders_; }
289 inline bool is_rational() const { return is_rational_; }
293 inline void set_is_rational(const bool in_is_rational) { is_rational_ = in_is_rational; }
294
297 virtual void normalize() = 0;
298
300 virtual void compute_elements_tag() = 0;
302 virtual void update_elements_tag() { assert(false); }
303
308 virtual double edge_length(const int gid) const
309 {
310 assert(false);
311 return 0;
312 }
317 virtual double quad_area(const int gid) const
318 {
319 assert(false);
320 return 0;
321 }
326 virtual double tri_area(const int gid) const
327 {
328 assert(false);
329 return 0;
330 }
331
336 virtual RowVectorNd point(const int global_index) const = 0;
341 virtual void set_point(const int global_index, const RowVectorNd &p) = 0;
342
347 virtual RowVectorNd edge_barycenter(const int e) const = 0;
352 virtual RowVectorNd face_barycenter(const int f) const = 0;
357 virtual RowVectorNd cell_barycenter(const int c) const = 0;
358
362 void edge_barycenters(Eigen::MatrixXd &barycenters) const;
366 void face_barycenters(Eigen::MatrixXd &barycenters) const;
370 void cell_barycenters(Eigen::MatrixXd &barycenters) const;
374 virtual void compute_element_barycenters(Eigen::MatrixXd &barycenters) const = 0;
375
379 virtual void elements_boxes(std::vector<std::array<Eigen::Vector3d, 2>> &boxes) const = 0;
385 virtual void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const = 0;
386
391 virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const = 0;
392
397 bool is_spline_compatible(const int el_id) const;
402 bool is_cube(const int el_id) const;
407 bool is_polytope(const int el_id) const;
412 bool is_simplex(const int el_id) const;
417 bool is_prism(const int el_id) const;
418
423 bool is_pyramid(const int el_id) const;
424
428 const std::vector<ElementType> &elements_tag() const { return elements_tag_; }
433 void set_tag(const int el, const ElementType type) { elements_tag_[el] = type; }
434
438 void compute_node_ids(const std::function<int(const size_t, const RowVectorNd &, bool)> &marker);
439
443 virtual void load_boundary_ids(const std::string &path);
444
448 virtual void compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker) = 0;
449
453 virtual void compute_body_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &)> &marker) = 0;
454
458 virtual void set_boundary_ids(const std::vector<int> &boundary_ids)
459 {
460 assert(boundary_ids.size() == n_boundary_elements());
461 boundary_ids_ = boundary_ids;
462 }
466 virtual void set_body_ids(const std::vector<int> &body_ids)
467 {
468 assert(body_ids.size() == n_elements());
469 body_ids_ = body_ids;
470 }
471
476 virtual int get_default_boundary_id(const int primitive) const
477 {
478 if (is_volume() ? is_boundary_face(primitive) : is_boundary_edge(primitive))
479 return std::numeric_limits<int>::max(); // default for no selected boundary
480 else
481 return -1; // default for no boundary
482 }
483
488 virtual int get_boundary_id(const int primitive) const
489 {
490 return has_boundary_ids() ? (boundary_ids_.at(primitive)) : get_default_boundary_id(primitive);
491 }
492
497 virtual int get_node_id(const int node_id) const
498 {
499 if (has_node_ids())
500 return node_ids_.at(node_id);
501 else
502 return -1; // default for no boundary
503 }
504
508 void update_nodes(const Eigen::VectorXi &in_node_to_node);
509
514 virtual int get_body_id(const int primitive) const
515 {
516 if (has_body_ids())
517 return body_ids_.at(primitive);
518 else
519 return 0;
520 }
523 virtual const std::vector<int> &get_body_ids() const
524 {
525 return body_ids_;
526 }
530 bool has_node_ids() const { return !node_ids_.empty(); }
534 bool has_boundary_ids() const { return !boundary_ids_.empty(); }
538 virtual bool has_body_ids() const { return !body_ids_.empty(); }
539
544 virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const = 0;
550 virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1, const std::vector<bool> &valid_elements) const = 0;
551
552 // /// @brief generate a triangular representation of every face
553 // ///
554 // /// @param[out] tris triangles connectivity
555 // /// @param[out] pts triangles vertices
556 // /// @param[out] ranges connection to original faces
557 // virtual void triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const = 0;
558
563 const std::vector<double> &cell_weights(const int cell_index) const { return cell_weights_[cell_index]; }
567 void set_cell_weights(const std::vector<std::vector<double>> &in_cell_weights) { cell_weights_ = in_cell_weights; }
568
571 virtual void prepare_mesh() {};
572
576 bool has_poly() const
577 {
578 for (int i = 0; i < n_elements(); ++i)
579 {
580 if (is_polytope(i))
581 return true;
582 }
583
584 return false;
585 }
586
590 bool has_prism() const
591 {
592 for (int i = 0; i < n_elements(); ++i)
593 {
594 if (is_prism(i))
595 return true;
596 }
597
598 return false;
599 }
600
604 bool has_pyramids() const
605 {
606 for (int i = 0; i < n_elements(); ++i)
607 {
608 if (is_pyramid(i))
609 return true;
610 }
611
612 return false;
613 }
614
618 bool is_simplicial() const
619 {
620 for (int i = 0; i < n_elements(); ++i)
621 {
622 if (!is_simplex(i))
623 return false;
624 }
625
626 return true;
627 }
628
632 inline bool is_linear() const { return orders_.size() == 0 || orders_.maxCoeff() == 1; }
633
637 std::vector<std::pair<int, int>> edges() const;
641 std::vector<std::vector<int>> faces() const;
642
646 std::unordered_map<std::pair<int, int>, size_t, polyfem::utils::HashPair> edges_to_ids() const;
650 std::unordered_map<std::vector<int>, size_t, polyfem::utils::HashVector> faces_to_ids() const;
651
655 inline const Eigen::VectorXi &in_ordered_vertices() const { return in_ordered_vertices_; }
659 inline const Eigen::MatrixXi &in_ordered_edges() const { return in_ordered_edges_; }
663 inline const Eigen::MatrixXi &in_ordered_faces() const { return in_ordered_faces_; }
664
668 virtual void append(const Mesh &mesh);
669
673 void append(const std::unique_ptr<Mesh> &mesh)
674 {
675 if (mesh != nullptr)
676 append(*mesh);
677 }
678
682 void apply_affine_transformation(const MatrixNd &A, const VectorNd &b);
683
684 protected:
689 virtual bool load(const std::string &path) = 0;
694 virtual bool load(const GEO::Mesh &M) = 0;
695
697 std::vector<ElementType> elements_tag_;
699 std::vector<int> node_ids_;
701 std::vector<int> boundary_ids_;
703 std::vector<int> body_ids_;
705 Eigen::MatrixXi orders_;
707 bool is_rational_ = false;
708
710 std::vector<EdgeNodes> edge_nodes_;
712 std::vector<FaceNodes> face_nodes_;
714 std::vector<CellNodes> cell_nodes_;
716 std::vector<std::vector<double>> cell_weights_;
717
719 Eigen::VectorXi in_ordered_vertices_;
721 Eigen::MatrixXi in_ordered_edges_;
723 Eigen::MatrixXi in_ordered_faces_;
724 };
725 } // namespace mesh
726} // namespace polyfem
int V
Class to store the high-order cells nodes.
Definition Mesh.hpp:63
Eigen::MatrixXd nodes
Definition Mesh.hpp:67
std::vector< int > nodes_ids
Definition Mesh.hpp:66
Class to store the high-order edge nodes.
Definition Mesh.hpp:45
Eigen::MatrixXd nodes
Definition Mesh.hpp:49
std::vector< int > nodes_ids
Definition Mesh.hpp:48
Class to store the high-order face nodes.
Definition Mesh.hpp:54
std::vector< int > nodes_ids
Definition Mesh.hpp:57
Eigen::MatrixXd nodes
Definition Mesh.hpp:58
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:41
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:163
virtual int n_vertices() const =0
number of vertices
virtual int get_body_id(const int primitive) const
Get the volume selection of an element (cell in 3d, face in 2d)
Definition Mesh.hpp:514
virtual void compute_element_barycenters(Eigen::MatrixXd &barycenters) const =0
utility for 2d/3d.
virtual double edge_length(const int gid) const
edge length
Definition Mesh.hpp:308
void set_cell_weights(const std::vector< std::vector< double > > &in_cell_weights)
Set the cell weights for rational polynomial meshes.
Definition Mesh.hpp:567
const Eigen::MatrixXi & in_ordered_edges() const
Order of the input edges.
Definition Mesh.hpp:659
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:365
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:705
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:469
bool is_rational() const
check if curved mesh has rational polynomials elements
Definition Mesh.hpp:289
Mesh()=default
Construct a new Mesh object.
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:697
virtual RowVectorNd edge_barycenter(const int e) const =0
edge barycenter
virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const =0
Get all the edges.
int element_vertex(const int el_id, const int lv_id) const
id of the vertex of a element
Definition Mesh.hpp:222
bool is_rational_
stores if the mesh is rational
Definition Mesh.hpp:707
bool has_boundary_ids() const
checks if surface selections are available
Definition Mesh.hpp:534
virtual void compute_boundary_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &, bool)> &marker)=0
computes boundary selections based on a function
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
Definition Mesh.hpp:326
void cell_barycenters(Eigen::MatrixXd &barycenters) const
all cells barycenters
Definition Mesh.cpp:320
virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1, const std::vector< bool > &valid_elements) const =0
Get all the edges according to valid_elements selection.
bool is_simplicial() const
checks if the mesh is simplicial
Definition Mesh.hpp:618
virtual bool is_conforming() const =0
if the mesh is conforming
virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const =0
computes the bbox of the mesh
virtual void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const =0
constructs barycentric coodiantes for a point p.
virtual RowVectorNd face_barycenter(const int f) const =0
face barycenter
virtual void update_elements_tag()
Update elements types.
Definition Mesh.hpp:302
Mesh & operator=(const Mesh &)=default
virtual void compute_elements_tag()=0
compute element types, see ElementType
virtual bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)=0
build a mesh from matrices
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:352
bool has_node_ids() const
checks if points selections are available
Definition Mesh.hpp:530
virtual void set_body_ids(const std::vector< int > &body_ids)
Set the volume sections.
Definition Mesh.hpp:466
virtual RowVectorNd point(const int global_index) const =0
point coordinates
const Eigen::MatrixXi & orders() const
order of each element
Definition Mesh.hpp:285
virtual bool is_boundary_face(const int face_global_id) const =0
is face boundary
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:723
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:387
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:488
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 simplex
Definition Mesh.cpp:422
void face_barycenters(Eigen::MatrixXd &barycenters) const
all face barycenters
Definition Mesh.cpp:311
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
Definition Mesh.hpp:317
virtual bool has_body_ids() const
checks if volumes selections are available
Definition Mesh.hpp:538
bool has_prism() const
checks if the mesh has prisms
Definition Mesh.hpp:590
virtual void normalize()=0
normalize the mesh
bool is_spline_compatible(const int el_id) const
checks if element is spline compatible
Definition Mesh.cpp:332
virtual void prepare_mesh()
method used to finalize the mesh.
Definition Mesh.hpp:571
virtual void load_boundary_ids(const std::string &path)
loads the boundary selections for a file
Definition Mesh.cpp:399
void set_is_rational(const bool in_is_rational)
Set the is rational object.
Definition Mesh.hpp:293
Mesh & operator=(Mesh &&)=default
Copy constructor.
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:701
bool is_prism(const int el_id) const
checks if element is a prism
Definition Mesh.cpp:427
bool is_linear() const
check if the mesh is linear
Definition Mesh.hpp:632
void apply_affine_transformation(const MatrixNd &A, const VectorNd &b)
Apply an affine transformation to the vertex positions .
Definition Mesh.cpp:642
std::vector< int > node_ids_
list of node labels
Definition Mesh.hpp:699
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:485
bool has_pyramids() const
checks if the mesh has pyramids
Definition Mesh.hpp:604
virtual void refine(const int n_refinement, const double t)=0
refine the mesh
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:714
std::vector< std::vector< double > > cell_weights_
weights associates to cells for rational polynomail meshes
Definition Mesh.hpp:716
const std::vector< double > & cell_weights(const int cell_index) const
weights for rational polynomial meshes
Definition Mesh.hpp:563
virtual bool load(const GEO::Mesh &M)=0
loads a mesh from a geo mesh
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:437
void update_nodes(const Eigen::VectorXi &in_node_to_node)
Update the node ids to reorder them.
Definition Mesh.cpp:371
virtual int edge_vertex(const int e_id, const int lv_id) const =0
id of the edge vertex
virtual std::unique_ptr< Mesh > copy() const =0
Create a copy of the mesh.
void append(const std::unique_ptr< Mesh > &mesh)
appends a new mesh to the end of this, utility that takes pointer, calls other one
Definition Mesh.hpp:673
virtual bool is_boundary_edge(const int edge_global_id) const =0
is edge boundary
bool has_poly() const
checks if the mesh has polytopes
Definition Mesh.hpp:576
std::vector< int > element_vertices(const int el_id) const
list of vids of an element
Definition Mesh.hpp:231
const Eigen::VectorXi & in_ordered_vertices() const
Order of the input vertices.
Definition Mesh.hpp:655
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:703
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:476
int dimension() const
utily for dimension
Definition Mesh.hpp:153
virtual int n_cells() const =0
number of cells
virtual const std::vector< int > & get_body_ids() const
Get the volume selection of all elements (cells in 3d, faces in 2d)
Definition Mesh.hpp:523
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:712
virtual bool load(const std::string &path)=0
loads a mesh from the path
virtual int n_faces() const =0
number of faces
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:428
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:710
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:453
virtual void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes)=0
attach high order nodes
int n_boundary_elements() const
utitlity to return the number of boundary elements, faces or edges in 3d and 2d
Definition Mesh.hpp:168
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:721
void set_tag(const int el, const ElementType type)
changes the element type
Definition Mesh.hpp:433
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:504
bool is_pyramid(const int el_id) const
checks if element is a pyramid
Definition Mesh.cpp:432
Mesh(Mesh &&)=default
Construct a new Mesh object.
virtual bool save(const std::string &path) const =0
virtual RowVectorNd cell_barycenter(const int c) const =0
cell barycenter
const Eigen::MatrixXi & in_ordered_faces() const
Order of the input edges.
Definition Mesh.hpp:663
virtual int n_edges() const =0
number of edges
virtual void set_boundary_ids(const std::vector< int > &boundary_ids)
Set the boundary selection from a vector.
Definition Mesh.hpp:458
virtual void compute_body_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &)> &marker)=0
computes boundary selections based on a function
void edge_barycenters(Eigen::MatrixXd &barycenters) const
all edges barycenters
Definition Mesh.cpp:302
virtual int cell_vertex(const int f_id, const int lv_id) const =0
id of the vertex of a cell
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
virtual void elements_boxes(std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const =0
constructs a box around every element (3d cell, 2d face)
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:719
virtual bool is_boundary_element(const int element_global_id) const =0
is cell boundary
virtual int n_cell_vertices(const int c_id) const =0
number of vertices of a cell
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:497
virtual ~Mesh()=default
Destroy the Mesh object.
virtual int face_vertex(const int f_id, const int lv_id) const =0
id of the face vertex
Mesh(const Mesh &)=default
Construct a new Mesh object.
int boundary_element_vertex(const int primitive_id, const int lv_id) const
Definition Mesh.hpp:240
ElementType
Type of Element, check [Poly-Spline Finite Element Method] for a complete description.
Definition Mesh.hpp:23
@ 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.
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor, 3, 3 > MatrixNd
Definition Types.hpp:14
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13