Collision Mesh¶
-
class CollisionMesh¶
A class for encapsolating the transformation/selections needed to go from a volumetric FE mesh to a surface collision mesh.
Public Functions¶
-
inline CollisionMesh()¶
Construct a new Collision Mesh object.
Collision Mesh objects are immutable, so use the other constructors.
-
CollisionMesh(const Eigen::MatrixXd &vertices_at_rest, const Eigen::MatrixXi &edges, const Eigen::MatrixXi &faces, const Eigen::SparseMatrix<double> &displacement_map = Eigen::SparseMatrix<double>())¶
Construct a new Collision Mesh object directly from the collision mesh vertices.
- Parameters:¶
- const Eigen::MatrixXd &vertices_at_rest¶
The vertices of the collision mesh at rest.
- const Eigen::MatrixXi &edges¶
The edges of the collision mesh.
- const Eigen::MatrixXi &faces¶
The faces of the collision mesh.
- const Eigen::SparseMatrix<double> &displacement_map = Eigen::SparseMatrix<double>()¶
The displacement mapping from displacements on the full mesh to the collision mesh.
-
CollisionMesh(const std::vector<bool> &include_vertex, const Eigen::MatrixXd &full_vertices_at_rest, const Eigen::MatrixXi &edges, const Eigen::MatrixXi &faces, const Eigen::SparseMatrix<double> &displacement_map = Eigen::SparseMatrix<double>())¶
Construct a new Collision Mesh object from a full mesh vertices.
- Parameters:¶
- const std::vector<bool> &include_vertex¶
Vector of bools indicating whether each vertex should be included in the collision mesh.
- const Eigen::MatrixXd &full_vertices_at_rest¶
The vertices of the full mesh at rest.
- const Eigen::MatrixXi &edges¶
The edges of the collision mesh indexed into the full mesh vertices.
- const Eigen::MatrixXi &faces¶
The faces of the collision mesh indexed into the full mesh vertices.
- const Eigen::SparseMatrix<double> &displacement_map = Eigen::SparseMatrix<double>()¶
The displacement mapping from displacements on the full mesh to the collision mesh.
-
void init_adjacencies()¶
Initialize vertex-vertex and edge-vertex adjacencies.
-
void init_area_jacobians()¶
Initialize vertex and edge areas.
-
inline ~CollisionMesh()¶
Destroy the Collision Mesh object.
-
inline size_t num_vertices() const¶
Get the number of vertices in the collision mesh.
-
inline size_t num_edges() const¶
Get the number of edges in the collision mesh.
-
inline size_t num_faces() const¶
Get the number of faces in the collision mesh.
-
inline size_t dim() const¶
Get the dimension of the mesh.
-
inline size_t ndof() const¶
Get the number of degrees of freedom in the collision mesh.
-
inline size_t full_num_vertices() const¶
Get the number of vertices in the full mesh.
-
inline size_t full_ndof() const¶
Get the number of degrees of freedom in the full mesh.
-
inline const Eigen::MatrixXd &vertices_at_rest() const¶
Get the vertices of the collision mesh at rest.
-
inline const Eigen::MatrixXi &edges() const¶
Get the edges of the collision mesh.
-
inline const Eigen::MatrixXi &faces() const¶
Get the faces of the collision mesh.
-
inline const Eigen::MatrixXi &faces_to_edges() const¶
Get the mapping from faces to edges of the collision mesh.
-
Eigen::MatrixXd vertices(const Eigen::MatrixXd &full_vertices) const¶
Compute the vertex positions from the positions of the full mesh.
-
Eigen::MatrixXd displace_vertices(const Eigen::MatrixXd &full_displacements) const¶
Compute the vertex positions from vertex displacements on the full mesh.
-
Eigen::MatrixXd map_displacements(const Eigen::MatrixXd &full_displacements) const¶
Map vertex displacements on the full mesh to vertex displacements on the collision mesh.
-
inline size_t to_full_vertex_id(const size_t id) const¶
Map a vertex ID to the corresponding vertex ID in the full mesh.
-
Eigen::VectorXd to_full_dof(const Eigen::VectorXd &x) const¶
Map a vector quantity on the collision mesh to the full mesh.
This is useful for mapping gradients from the collision mesh to the full mesh (i.e., applies the chain-rule).
- Parameters:¶
- Returns:¶
Vector quantity on the full mesh with size equal to full_ndof().
-
Eigen::SparseMatrix<double> to_full_dof(const Eigen::SparseMatrix<double> &X) const¶
Map a matrix quantity on the collision mesh to the full mesh.
This is useful for mapping Hessians from the collision mesh to the full mesh (i.e., applies the chain-rule).
- Parameters:¶
- Returns:¶
Matrix quantity on the full mesh with size equal to full_ndof() × full_ndof().
-
inline const std::vector<unordered_set<int>> &vertex_vertex_adjacencies() const¶
Get the vertex-vertex adjacency matrix.
-
inline const std::vector<unordered_set<int>> &edge_vertex_adjacencies() const¶
Get the edge-vertex adjacency matrix.
-
inline bool are_adjacencies_initialized() const¶
Determine if the adjacencies have been initialized by calling init_adjacencies().
-
inline bool is_vertex_on_boundary(const int vi) const¶
Is a vertex on the boundary of the collision mesh?
-
inline const Eigen::VectorXd &vertex_areas() const¶
Get the barycentric area of the vertices.
-
inline const Eigen::SparseVector<double> &vertex_area_gradient(const size_t vi) const¶
Get the gradient of the barycentric area of a vertex wrt the rest positions of all points.
-
inline const Eigen::VectorXd &edge_areas() const¶
Get the barycentric area of the edges.
-
inline const Eigen::SparseVector<double> &edge_area_gradient(const size_t ei) const¶
Get the gradient of the barycentric area of an edge wrt the rest positions of all points.
-
inline bool are_area_jacobians_initialized() const¶
Determine if the area Jacobians have been initialized by calling init_area_jacobians().
Public Members¶
-
std::function<bool(size_t, size_t)> can_collide = default_can_collide¶
A function that takes two vertex IDs (row numbers in V) and returns true if the vertices (and faces or edges containing the vertices) can collide.
By default all primitives can collide with all other primitives.
Public Static Functions¶
-
static inline CollisionMesh build_from_full_mesh(const Eigen::MatrixXd &full_vertices_at_rest, const Eigen::MatrixXi &edges, const Eigen::MatrixXi &faces)¶
Helper function that automatically builds include_vertex using construct_is_on_surface.
-
static std::vector<bool> construct_is_on_surface(const int num_vertices, const Eigen::MatrixXi &edges)¶
Construct a vector of bools indicating whether each vertex is on the surface.
-
inline CollisionMesh()¶