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

CollisionMesh() = default

Construct a new Collision Mesh object.

Collision Mesh objects are immutable, so use the other constructors.

CollisionMesh(const Eigen::MatrixXd& rest_positions,
    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 &rest_positions

The vertices of the collision mesh at rest (#V × dim).

const Eigen::MatrixXi &edges

The edges of the collision mesh (#E × 2).

const Eigen::MatrixXi &faces

The faces of the collision mesh (#F × 3).

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_rest_positions
,
    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_rest_positions

The vertices of the full mesh at rest (#V × dim).

const Eigen::MatrixXi &edges

The edges of the collision mesh indexed into the full mesh vertices (#E × 2).

const Eigen::MatrixXi &faces

The faces of the collision mesh indexed into the full mesh vertices (#F × 3).

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.

~CollisionMesh() = default

Destroy the Collision Mesh object.

inline size_t num_vertices() const

Get the number of vertices in the collision mesh.

inline size_t num_codim_vertices() const

Get the number of codimensional vertices in the collision mesh.

inline size_t num_codim_edges() const

Get the number of codimensional edges 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& rest_positions() const

Get the vertices of the collision mesh at rest (#V × dim).

inline const Eigen::VectorXi& codim_vertices() const

Get the indices of codimensional vertices of the collision mesh (#CV x 1).

inline const Eigen::VectorXi& codim_edges() const

Get the indices of codimensional edges of the collision mesh (#CE x 1).

inline const Eigen::MatrixXi& edges() const

Get the edges of the collision mesh (#E × 2).

inline const Eigen::MatrixXi& faces() const

Get the faces of the collision mesh (#F × 3).

inline const Eigen::MatrixXi& faces_to_edges() const

Get the mapping from faces to edges of the collision mesh (#F × 3).

Eigen::MatrixXd vertices(
    const
 Eigen::MatrixXd& full_positions
)
 const

Compute the vertex positions from the positions of the full mesh.

Parameters:
const Eigen::MatrixXd &full_positions

The vertex positions of the full mesh (#FV × dim).

Returns:

The vertex positions of the collision mesh (#V × dim).

Eigen::MatrixXd displace_vertices(
    const
 Eigen::MatrixXd& full_displacements
)
 const

Compute the vertex positions from vertex displacements on the full mesh.

Parameters:
const Eigen::MatrixXd &full_displacements

The vertex displacements on the full mesh (#FV × dim).

Returns:

The vertex positions of the collision mesh (#V × dim).

Eigen::MatrixXd map_displacements(
    const
 Eigen::MatrixXd& full_displacements
)
 const

Map vertex displacements on the full mesh to vertex displacements on the collision mesh.

Parameters:
const Eigen::MatrixXd &full_displacements

The vertex displacements on the full mesh (#FV × dim).

Returns:

The vertex displacements on the collision mesh (#V × dim).

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.

Parameters:
const size_t id

Vertex ID in the collision mesh.

Returns:

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:
const Eigen::VectorXd &x

Vector quantity on the collision mesh with size equal to ndof().

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:
const Eigen::SparseMatrix<double> &X

Matrix quantity on the collision mesh with size equal to ndof() × ndof().

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>>&
vertex_edge_adjacencies
() const

Get the vertex-edge 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?

Parameters:
const int vi

Vertex ID.

Returns:

True if the vertex is on the boundary of the collision mesh.

inline double vertex_area(const size_t vi) const

Get the barycentric area of a vertex.

Parameters:
const size_t vi

Vertex ID.

Returns:

Barycentric area of vertex vi.

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.

Parameters:
const size_t vi

Vertex ID.

Returns:

Gradient of the barycentric area of vertex vi wrt the rest positions of all points.

inline double edge_area(const size_t ei) const

Get the barycentric area of an edge.

Parameters:
const size_t ei

Edge ID.

Returns:

Barycentric area of edge ei.

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.

Parameters:
const size_t ei

Edge ID.

Returns:

Gradient of the barycentric area of edge ei 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 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_rest_positions
,
    const
 Eigen::MatrixXi& edges
, const Eigen::MatrixXi& faces)

Helper function that automatically builds include_vertex using construct_is_on_surface.

Parameters:
const Eigen::MatrixXd &full_rest_positions

The full vertices at rest (#FV × dim).

const Eigen::MatrixXi &edges

The edge matrix of mesh (#E × 2).

const Eigen::MatrixXi &faces

The face matrix of mesh (#F × 3).

Returns:

Constructed CollisionMesh.

static std::vector<bool> construct_is_on_surface(
    const
 long num_vertices
, const Eigen::MatrixXi& edges,
    const
 Eigen::VectorXi& codim_vertices = Eigen::VectorXi()
)

Construct a vector of bools indicating whether each vertex is on the surface.

Parameters:
const long num_vertices

The number of vertices in the mesh.

const Eigen::MatrixXi &edges

The surface edges of the mesh (#E × 2).

const Eigen::VectorXi &codim_vertices = Eigen::VectorXi()

The indices of codimensional vertices (#CV x 1).

Returns:

A vector of bools indicating whether each vertex is on the surface.

static Eigen::MatrixXi construct_faces_to_edges(
    const
 Eigen::MatrixXi& faces
, const Eigen::MatrixXi& edges)

Construct a matrix that maps from the faces’ edges to rows in the edges matrix.

Parameters:
const Eigen::MatrixXi &faces

The face matrix of mesh (#F × 3).

const Eigen::MatrixXi &edges

The edge matrix of mesh (#E × 2).

Returns:

Matrix that maps from the faces’ edges to rows in the edges matrix.

Protected Functions

void init_codim_vertices()
void init_codim_edges()
void init_selection_matrices(const int dim)

Initialize the selection matrix from full vertices/DOF to collision vertices/DOF.

void init_areas()

Initialize vertex and edge areas.

Protected Attributes

Eigen::MatrixXd m_full_rest_positions

The full vertex positions at rest (#FV × dim).

Eigen::MatrixXd m_rest_positions

The vertex positions at rest (#V × dim).

Eigen::VectorXi m_codim_vertices

The indices of codimensional vertices (#CV x 1).

Eigen::VectorXi m_codim_edges

The indices of codimensional edges (#CE x 1).

Eigen::MatrixXi m_edges

Edges as rows of indicies into vertices (#E × 2).

Eigen::MatrixXi m_faces

Triangular faces as rows of indicies into vertices (#F × 3).

Eigen::MatrixXi m_faces_to_edges

Map from faces edges to rows of edges (#F × 3).

Eigen::VectorXi m_full_vertex_to_vertex

Map from full vertices to collision vertices.

Note

Negative values indicate full vertex is dropped.

Eigen::VectorXi m_vertex_to_full_vertex

Map from collision vertices to full vertices.

Eigen::SparseMatrix<double> m_select_vertices

Selection matrix S ∈ ℝ^{collision×full} for vertices.

Eigen::SparseMatrix<double> m_select_dof

Selection matrix S ∈ ℝ^{collision×full} for DOF.

Eigen::SparseMatrix<double> m_displacement_map

Mapping from full displacements to collision displacements.

Note

this is premultiplied by m_select_vertices

Eigen::SparseMatrix<double> m_displacement_dof_map

Mapping from full displacements to collision displacements.

Note

this is premultiplied by m_select_dof

std::vector<unordered_set<int>> m_vertex_vertex_adjacencies

Vertices adjacent to vertices.

std::vector<unordered_set<int>> m_vertex_edge_adjacencies

Edges adjacent to vertices.

std::vector<unordered_set<int>> m_edge_vertex_adjacencies

Vertices adjacent to edges.

std::vector<bool> m_is_vertex_on_boundary

Is vertex on the boundary of the triangle mesh in 3D or polyline in 2D?

Eigen::VectorXd m_vertex_areas

Vertex areas 2D: 1/2 sum of length of connected edges 3D: 1/3 sum of area of connected triangles.

Eigen::VectorXd m_edge_areas

Edge areas 3D: 1/3 sum of area of connected triangles.

std::vector<Eigen::SparseVector<double>> m_vertex_area_jacobian

The rows of the Jacobian of the vertex areas vector.

std::vector<Eigen::SparseVector<double>> m_edge_area_jacobian

The rows of the Jacobian of the edge areas vector.

Protected Static Functions

static Eigen::SparseMatrix<double> vertex_matrix_to_dof_matrix(
    const
 Eigen::SparseMatrix<double>& M_V
, int dim)

Convert a matrix meant for M_V * vertices to M_dof * x by duplicating the entries dim times.

Private Static Functions

static inline int default_can_collide(size_t, size_t)

By default all primitives can collide with all other primitives.


Last update: Mar 28, 2024