Candidates

Candidates

class Candidates;

A class for storing and managing collision candidates.

Public Functions

Candidates() = default;
void build(const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const double inflation_radius = 0,
    
BroadPhasebroad_phase = nullptr);

Initialize the set of discrete collision detection candidates.

Parameters:
const CollisionMesh &mesh

The surface of the collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Surface vertex positions (rowwise).

const double inflation_radius = 0

Amount to inflate the bounding boxes.

BroadPhase *broad_phase = nullptr

Broad phase method to use.

void build(const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
    
const double inflation_radius = 0,
    
BroadPhasebroad_phase = nullptr);

Initialize the set of continuous collision detection candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The surface of the collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices_t0

Surface vertex starting positions (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> vertices_t1

Surface vertex ending positions (rowwise).

const double inflation_radius = 0

Amount to inflate the bounding boxes.

BroadPhase *broad_phase = nullptr

Broad phase method to use.

size_t size() const;

Get the number of collision candidates.

Returns:

The number of collision candidates.

bool empty() const;

Check if there are no collision candidates.

Returns:

True if there are no collision candidates, false otherwise.

void clear();

Clear all collision candidates.

CollisionStenciloperator[](size_t i);

Get a collision stencil by index.

Parameters:
size_t i

The index of the collision stencil.

Returns:

A reference to the collision stencil.

const CollisionStenciloperator[](size_t i) const;

Get a collision stencil by index.

Parameters:
size_t i

The index of the collision stencil.

Returns:

A const reference to the collision stencil.

bool is_step_collision_free(const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
    
const double min_distance = 0.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const;

Determine if the step is collision free from the set of candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices_t0

Surface vertex starting positions (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> vertices_t1

Surface vertex ending positions (rowwise).

const double min_distance = 0.0

The minimum distance allowable between any two elements.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

The narrow phase CCD algorithm to use.

Returns:

True if any collisions occur.

double compute_collision_free_stepsize(const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1,
    
const double min_distance = 0.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const;

Computes a maximal step size that is collision free using the set of collision candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices_t0

Surface vertex starting positions (rowwise). Assumed to be intersection free.

Eigen::ConstRef<Eigen::MatrixXd> vertices_t1

Surface vertex ending positions (rowwise).

const double min_distance = 0.0

The minimum distance allowable between any two elements.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

The narrow phase CCD algorithm to use.

Returns:

A step-size \(\in [0, 1]\) that is collision free. A value of 1.0 if a full step and 0.0 is no step.

double compute_noncandidate_conservative_stepsize(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> displacements,
    
const double dhat) const;

Computes a conservative bound on the largest-feasible step size for surface primitives not in collision.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> displacements

Surface vertex displacements (rowwise).

const double dhat

Barrier activation distance.

Returns:

A step-size \(\in [0, 1]\) that is collision free for non-candidate elements.

double compute_cfl_stepsize(const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices_t1const double dhat,
    
const double min_distance = 0.0,
    
BroadPhasebroad_phase = nullptr,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const;

Computes a CFL-inspired CCD maximum step step size.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices_t0

Surface vertex starting positions (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> vertices_t1

Surface vertex ending positions (rowwise).

const double dhat

Barrier activation distance.

const double min_distance = 0.0

The minimum distance allowable between any two elements.

BroadPhase *broad_phase = nullptr

The broad phase algorithm to use.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

The narrow phase CCD algorithm to use.

Returns:

A step-size \(\in [0, 1]\) that is collision free.

Eigen::VectorXd compute_per_vertex_safe_distances(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const double inflation_radius,
    
const double min_distance = 0.0) const;

Compute the maximum distance every vertex can move (independently) without colliding with any other element.

Note

Cap the value at the inflation radius used to build the candidates.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

const double inflation_radius

The inflation radius used to build the candidates.

const double min_distance = 0.0

The minimum distance allowable between any two elements.

Returns:

A vector of minimum distances, one for each vertex.

std::vector<VertexVertexCandidate> edge_vertex_to_vertex_vertex(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const std::function<bool(double)>& is_active
    
= default_is_active
) const;

Convert edge-vertex candidates to vertex-vertex candidates.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

const std::function<bool(double)> &is_active = default_is_active

(Optional) Function to determine if a candidate is active.

Returns:

Vertex-vertex candidates derived from edge-vertex candidates.

std::vector<VertexVertexCandidate> face_vertex_to_vertex_vertex(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const std::function<bool(double)>& is_active
    
= default_is_active
) const;

Convert face-vertex candidates to vertex-vertex candidates.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

const std::function<bool(double)> &is_active = default_is_active

(Optional) Function to determine if a candidate is active.

Returns:

Vertex-vertex candidates derived from face-vertex candidates.

std::vector<EdgeVertexCandidate> face_vertex_to_edge_vertex(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const std::function<bool(double)>& is_active
    
= default_is_active
) const;

Convert face-vertex candidates to edge-vertex candidates.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

const std::function<bool(double)> &is_active = default_is_active

(Optional) Function to determine if a candidate is active.

Returns:

Edge-vertex candidates derived from face-vertex candidates.

std::vector<EdgeVertexCandidate> edge_edge_to_edge_vertex(
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
const std::function<bool(double)>& is_active
    
= default_is_active
) const;

Convert edge-edge candidates to edge-vertex candidates.

Parameters:
const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

const std::function<bool(double)> &is_active = default_is_active

(Optional) Function to determine if a candidate is active.

Returns:

Edge-vertex candidates derived from edge-edge candidates.

bool save_obj(const std::stringfilename,
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Save the collision candidates to an OBJ file.

Parameters:
const std::string &filename

The name of the file to save the candidates to.

Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertex positions (rowwise).

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edge indices (rowwise).

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh face indices (rowwise).

Returns:

True if the file was saved successfully.

Public Members

std::vector<VertexVertexCandidate> vv_candidates;
std::vector<EdgeVertexCandidate> ev_candidates;
std::vector<EdgeEdgeCandidate> ee_candidates;
std::vector<FaceVertexCandidate> fv_candidates;
std::vector<PlaneVertexCandidate> pv_candidates;

Private Static Functions

static inline bool default_is_active(double candidate);

Collision Stencil

class CollisionStencil;

Inheritance diagram for ipc::CollisionStencil:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "1" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil" fillcolor="#BFBFBF"] "2" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate"] "3" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision"] "4" [label="ipc::EdgeEdgeTangentialCollision" tooltip="ipc::EdgeEdgeTangentialCollision"] "5" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate"] "6" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision"] "7" [label="ipc::EdgeVertexTangentialCollision" tooltip="ipc::EdgeVertexTangentialCollision"] "8" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate"] "9" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision"] "10" [label="ipc::FaceVertexTangentialCollision" tooltip="ipc::FaceVertexTangentialCollision"] "11" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "14" [label="ipc::PlaneVertexCandidate" tooltip="ipc::PlaneVertexCandidate"] "12" [label="ipc::PlaneVertexNormalCollision" tooltip="ipc::PlaneVertexNormalCollision"] "15" [label="ipc::PlaneVertexTangentialCollision" tooltip="ipc::PlaneVertexTangentialCollision"] "16" [label="ipc::TangentialCollision" tooltip="ipc::TangentialCollision"] "18" [label="ipc::VertexVertexCandidate" tooltip="ipc::VertexVertexCandidate"] "13" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision"] "17" [label="ipc::VertexVertexTangentialCollision" tooltip="ipc::VertexVertexTangentialCollision"] "2" -> "1" [dir=forward tooltip="public-inheritance"] "3" -> "2" [dir=forward tooltip="public-inheritance"] "3" -> "11" [dir=forward tooltip="public-inheritance"] "4" -> "2" [dir=forward tooltip="public-inheritance"] "4" -> "16" [dir=forward tooltip="public-inheritance"] "5" -> "1" [dir=forward tooltip="public-inheritance"] "6" -> "5" [dir=forward tooltip="public-inheritance"] "6" -> "11" [dir=forward tooltip="public-inheritance"] "7" -> "5" [dir=forward tooltip="public-inheritance"] "7" -> "16" [dir=forward tooltip="public-inheritance"] "8" -> "1" [dir=forward tooltip="public-inheritance"] "9" -> "8" [dir=forward tooltip="public-inheritance"] "9" -> "11" [dir=forward tooltip="public-inheritance"] "10" -> "8" [dir=forward tooltip="public-inheritance"] "10" -> "16" [dir=forward tooltip="public-inheritance"] "11" -> "1" [dir=forward tooltip="public-inheritance"] "14" -> "1" [dir=forward tooltip="public-inheritance"] "12" -> "11" [dir=forward tooltip="public-inheritance"] "12" -> "14" [dir=forward tooltip="public-inheritance"] "15" -> "14" [dir=forward tooltip="public-inheritance"] "15" -> "16" [dir=forward tooltip="public-inheritance"] "16" -> "1" [dir=forward tooltip="public-inheritance"] "18" -> "1" [dir=forward tooltip="public-inheritance"] "13" -> "11" [dir=forward tooltip="public-inheritance"] "13" -> "18" [dir=forward tooltip="public-inheritance"] "17" -> "16" [dir=forward tooltip="public-inheritance"] "17" -> "18" [dir=forward tooltip="public-inheritance"] }

A stencil representing a collision between at most four vertices.

Subclassed by ipc::EdgeEdgeCandidate, ipc::EdgeVertexCandidate, ipc::FaceVertexCandidate, ipc::NormalCollision, ipc::PlaneVertexCandidate, ipc::TangentialCollision, ipc::VertexVertexCandidate

Public Functions

virtual ~CollisionStencil() = default;
template <typename CollisionStencilType>
inline CollisionStencilTypeas();

Get this as a child collision type.

Template Parameters:
typename CollisionStencilType

The child collision type to get this as.

Throws:

std::bad_cast – if this is not of type CollisionStencilType.

Returns:

This as the child collision type.

template <typename CollisionStencilType>
inline const CollisionStencilTypeas() const;

Get this as a child collision type.

Template Parameters:
typename CollisionStencilType

The child collision type to get this as.

Throws:

std::bad_cast – if this is not of type CollisionStencilType.

Returns:

This as the child collision type.

virtual int num_vertices() const = 0;

Get the number of vertices in the collision stencil.

inline int dim(const int ndof) const;

Get the dimension of the collision stencil.

Parameters:
const int ndof

Number of degrees of freedom in the stencil.

Returns:

The dimension of the collision stencil.

virtual std::array<index_t, STENCIL_SIZE> vertex_ids(
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const
   
 = 0;

Get the vertex IDs of the collision stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

The vertex IDs of the collision stencil. Elements i > num_vertices() are -1.

inline std::array<VectorMax3d, STENCIL_SIZE> vertices(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Get the vertex attributes of the collision stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Vertex attributes

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

The vertex positions of the collision stencil. Elements i > num_vertices() are NaN.

inline VectorMax12d dof(Eigen::ConstRef<Eigen::MatrixXd> X,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Select this stencil’s DOF from the full matrix of DOF.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> X

Full matrix of DOF (rowwise).

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

This stencil’s DOF.

inline double compute_distance(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

inline VectorMax4d compute_coefficients(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the coefficients of the stencil s.t.

\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Coefficients of the stencil.

inline VectorMax3d compute_distance_vector(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance vector using the mesh vertices.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices.

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges.

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces.

Returns:

The distance vector (dim-dimensional).

inline VectorMax3d compute_normal(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces,
    
const bool flip_if_negative = true,
    
double* sign = nullptr) const;

Compute the normal of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

const bool flip_if_negative = true

If true, flip the normal if the point is on the negative side.

double *sign = nullptr

If not nullptr, set to the sign of the normal before any flipping.

Returns:

Normal of the stencil.

inline MatrixMax<double, 3, 12> compute_normal_jacobian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces,
    
const bool flip_if_negative = true) const;

Compute the Jacobian of the normal of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

const bool flip_if_negative = true

If true, flip the normal if the point is on the negative side.

Returns:

Jacobian of the normal of the stencil.

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

VectorMax3d compute_distance_vector(
    
Eigen::ConstRef<VectorMax12d> positions) const;

Compute the distance vector of the stencil: t = ∑ cᵢ xᵢ.

The distance vector is the vector between the closest points on the collision primitives. Its squared norm equals the squared distance.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

The distance vector (dim-dimensional, i.e., 2D or 3D).

VectorMax3d compute_distance_vector(
    
Eigen::ConstRef<VectorMax12d> positions,
    
VectorMax4dcoeffs) const;

Compute the distance vector and the coefficients together.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

VectorMax4d &coeffs

[out] The computed coefficients cᵢ.

Returns:

The distance vector.

MatrixMax<double, 12, 3> compute_distance_vector_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const;

Compute the Jacobian of the distance vector w.r.t.

positions: ∂t/∂x = [c₀I c₁I … cₙI]ᵀ.

Since ∂t/∂x has a very simple structure (block-diagonal with scalar coefficients times identity), many operations can be done without forming this matrix. This method is provided for completeness and verification.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

The Jacobian ∂t/∂x ∈ ℝ^{ndof × dim}.

VectorMax3d compute_normal(Eigen::ConstRef<VectorMax12d> positions,
    
bool flip_if_negative = truedouble* sign = nullptr) const;

Compute the normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

bool flip_if_negative = true

If true, flip the normal if the point is on the negative side.

double *sign = nullptr

If not nullptr, set to the sign of the normal before any flipping.

Returns:

Normal of the stencil.

MatrixMax<double, 3, 12> compute_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions,
    
bool flip_if_negative = true) const;

Compute the Jacobian of the normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

bool flip_if_negative = true

If true, flip the normal if the point is on the negative side.

Returns:

Jacobian of the normal of the stencil.

virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1double& toi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const
   
 = 0;

Perform narrow-phase CCD on the candidate.

Parameters:
Eigen::ConstRef<VectorMax12d> vertices_t0

[in] Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

[in] Stencil vertices at the end of the time step.

double &toi

[out] Computed time of impact (normalized).

const double min_distance = 0.0

[in] Minimum separation distance between primitives.

const double tmax = 1.0

[in] Maximum time (normalized) to look for collisions.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

[in] The narrow phase CCD algorithm to use.

Returns:

If the candidate had a collision over the time interval.

std::ostreamwrite_ccd_query(std::ostreamout,
    
Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1) const;

Write the CCD query to a stream.

Parameters:
std::ostream &out

Stream to write to.

Eigen::ConstRef<VectorMax12d> vertices_t0

Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

Stencil vertices at the end of the time step.

Returns:

The stream.

Public Static Functions

static VectorMax12d diag_distance_vector_outer(
    
Eigen::ConstRef<VectorMax4d> coeffsconst int d);

Compute diag((∂t/∂x)(∂t/∂x)ᵀ) efficiently (Eq.

11 of the paper).

Result is [c₀², c₀², c₀², c₁², c₁², c₁², …] (each cᵢ² repeated dim times). This is used for computing diag(∂²b/∂x²) without the full Hessian.

Parameters:
Eigen::ConstRef<VectorMax4d> coeffs

The coefficients cᵢ (from compute_coefficients).

const int d

The spatial dimension (2 or 3).

Returns:

The diagonal of (∂t/∂x)(∂t/∂x)ᵀ as a vector of size ndof.

static VectorMax12d diag_distance_vector_t_outer(
    
Eigen::ConstRef<VectorMax4d> coeffs,
    
Eigen::ConstRef<VectorMax3d> distance_vector);

Compute diag((∂t/∂x · t)(∂t/∂x · t)ᵀ) efficiently (Eq.

12).

Result is element-wise square of [c₀tᵀ, c₁tᵀ, …, cₙtᵀ]. This is used for computing diag(∂²b/∂x²) without the full Hessian.

Parameters:
Eigen::ConstRef<VectorMax4d> coeffs

The coefficients cᵢ.

Eigen::ConstRef<VectorMax3d> distance_vector

The distance vector t.

Returns:

The diagonal of (∂t/∂x·t)(∂t/∂x·t)ᵀ as a vector of size ndof.

static VectorMax3d contract_distance_vector_jacobian(
    
Eigen::ConstRef<VectorMax4d> coeffs,
    
Eigen::ConstRef<VectorMax12d> pconst int d);

Compute pᵀ(∂t/∂x) efficiently as ∑ cᵢ pᵢ (Eqs.

13-14).

Given p = [p₀, p₁, …, pₙ]ᵀ where pᵢ ∈ ℝ^dim, this computes pᵀ(∂t/∂x) = ∑ cᵢ pᵢ which is a dim-dimensional vector.

Parameters:
Eigen::ConstRef<VectorMax4d> coeffs

The coefficients cᵢ.

Eigen::ConstRef<VectorMax12d> p

A vector of size ndof (the direction for the quadratic form).

const int d

The spatial dimension (2 or 3).

Returns:

pᵀ(∂t/∂x) as a dim-dimensional vector.

Public Static Attributes

static constexpr int STENCIL_SIZE = 4;

The maximum number of vertices in a collision stencil.

Protected Functions

virtual VectorMax3d compute_unnormalized_normal(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Unnormalized normal of the stencil.

virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the Jacobian of the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Jacobian of the unnormalized normal of the stencil.

Vertex-Vertex Candidate

class VertexVertexCandidate : public virtual ipc::CollisionStencil;

Inheritance diagram for ipc::VertexVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::VertexVertexCandidate" tooltip="ipc::VertexVertexCandidate" fillcolor="#BFBFBF"] "3" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision"] "4" [label="ipc::VertexVertexTangentialCollision" tooltip="ipc::VertexVertexTangentialCollision"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::VertexVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::VertexVertexCandidate" tooltip="ipc::VertexVertexCandidate" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

A candidate for vertex-vertex collision detection.

Subclassed by ipc::VertexVertexNormalCollision, ipc::VertexVertexTangentialCollision

Public Functions

VertexVertexCandidate(index_t vertex0_idindex_t vertex1_id);
inline virtual int num_vertices() const override;

Get the number of vertices in the collision stencil.

inline virtual std::array<index_t, 4> vertex_ids(
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;

Get the indices of the vertices.

Parameters:
Eigen::ConstRef<Eigen::MatrixXi> edges

edge matrix of mesh

Eigen::ConstRef<Eigen::MatrixXi> faces

face matrix of mesh

Returns:

List of vertex indices

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1double& toi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const override;

Perform narrow-phase CCD on the candidate.

Parameters:
Eigen::ConstRef<VectorMax12d> vertices_t0

[in] Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

[in] Stencil vertices at the end of the time step.

double &toi

[out] Computed time of impact (normalized).

const double min_distance = 0.0

[in] Minimum separation distance between primitives.

const double tmax = 1.0

[in] Maximum time (normalized) to look for collisions.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

[in] The narrow phase CCD algorithm to use.

Returns:

If the candidate had a collision over the time interval.

bool operator==(const VertexVertexCandidateother) const;
bool operator!=(const VertexVertexCandidateother) const;
bool operator<(const VertexVertexCandidateother) const;

Compare EdgeVertexCandidates for sorting.

inline VectorMax4d compute_coefficients(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the coefficients of the stencil s.t.

\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Coefficients of the stencil.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

inline double compute_distance(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance of the stencil.

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

Public Members

index_t vertex0_id;

ID of the first vertex.

index_t vertex1_id;

ID of the second vertex.

Protected Functions

virtual VectorMax3d compute_unnormalized_normal(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Unnormalized normal of the stencil.

virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the Jacobian of the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Jacobian of the unnormalized normal of the stencil.

Friends

template <typename H>
inline friend H AbslHashValue(H hconst VertexVertexCandidatevv);

Edge-Vertex Candidate

class EdgeVertexCandidate : public virtual ipc::CollisionStencil;

Inheritance diagram for ipc::EdgeVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate" fillcolor="#BFBFBF"] "3" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision"] "4" [label="ipc::EdgeVertexTangentialCollision" tooltip="ipc::EdgeVertexTangentialCollision"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::EdgeVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

A candidate for edge-vertex collision detection.

Subclassed by ipc::EdgeVertexNormalCollision, ipc::EdgeVertexTangentialCollision

Public Functions

EdgeVertexCandidate(index_t edge_idindex_t vertex_id);
inline virtual int num_vertices() const override;

Get the number of vertices in the collision stencil.

inline virtual std::array<index_t, 4> vertex_ids(
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;

Get the vertex IDs for the edge-vertex pair.

Parameters:
Eigen::ConstRef<Eigen::MatrixXi> edges

The edge connectivity matrix

Eigen::ConstRef<Eigen::MatrixXi> faces

The face connectivity matrix

Returns:

An array of vertex IDs in the order: [vi, e0i, e1i, -1]

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1double& toi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const override;

Perform narrow-phase CCD on the candidate.

Parameters:
Eigen::ConstRef<VectorMax12d> vertices_t0

[in] Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

[in] Stencil vertices at the end of the time step.

double &toi

[out] Computed time of impact (normalized).

const double min_distance = 0.0

[in] Minimum separation distance between primitives.

const double tmax = 1.0

[in] Maximum time (normalized) to look for collisions.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

[in] The narrow phase CCD algorithm to use.

Returns:

If the candidate had a collision over the time interval.

inline virtual PointEdgeDistanceType known_dtype() const;
bool operator==(const EdgeVertexCandidateother) const;
bool operator!=(const EdgeVertexCandidateother) const;
bool operator<(const EdgeVertexCandidateother) const;

Compare EdgeVertexCandidates for sorting.

inline VectorMax4d compute_coefficients(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the coefficients of the stencil s.t.

\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Coefficients of the stencil.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

inline double compute_distance(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance of the stencil.

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

Public Members

index_t edge_id;

ID of the edge.

index_t vertex_id;

ID of the vertex.

Protected Functions

virtual VectorMax3d compute_unnormalized_normal(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Unnormalized normal of the stencil.

virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the Jacobian of the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Jacobian of the unnormalized normal of the stencil.

Friends

template <typename H>
inline friend H AbslHashValue(H hconst EdgeVertexCandidateev);

Edge-Edge Candidate

class EdgeEdgeCandidate : public virtual ipc::CollisionStencil;

Inheritance diagram for ipc::EdgeEdgeCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate" fillcolor="#BFBFBF"] "3" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision"] "4" [label="ipc::EdgeEdgeTangentialCollision" tooltip="ipc::EdgeEdgeTangentialCollision"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::EdgeEdgeCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

A candidate for edge-edge collision detection.

Subclassed by ipc::EdgeEdgeNormalCollision, ipc::EdgeEdgeTangentialCollision

Public Functions

EdgeEdgeCandidate(index_t edge0_idindex_t edge1_id);
inline virtual int num_vertices() const override;

Get the number of vertices in the collision stencil.

inline virtual std::array<index_t, 4> vertex_ids(
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;

Get the vertex IDs for the edge-edge pair.

Parameters:
Eigen::ConstRef<Eigen::MatrixXi> edges

The edge connectivity matrix

Eigen::ConstRef<Eigen::MatrixXi> faces

The face connectivity matrix

Returns:

An array of vertex IDs in the order: [ea0i, ea1i, eb0i, eb1i]

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1double& toi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const override;

Perform narrow-phase CCD on the candidate.

Parameters:
Eigen::ConstRef<VectorMax12d> vertices_t0

[in] Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

[in] Stencil vertices at the end of the time step.

double &toi

[out] Computed time of impact (normalized).

const double min_distance = 0.0

[in] Minimum separation distance between primitives.

const double tmax = 1.0

[in] Maximum time (normalized) to look for collisions.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

[in] The narrow phase CCD algorithm to use.

Returns:

If the candidate had a collision over the time interval.

inline virtual EdgeEdgeDistanceType known_dtype() const;
bool operator==(const EdgeEdgeCandidateother) const;
bool operator!=(const EdgeEdgeCandidateother) const;
bool operator<(const EdgeEdgeCandidateother) const;

Compare EdgeEdgeCandidates for sorting.

inline VectorMax4d compute_coefficients(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the coefficients of the stencil s.t.

\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Coefficients of the stencil.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

inline double compute_distance(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance of the stencil.

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

Public Members

index_t edge0_id;

ID of the first edge.

index_t edge1_id;

ID of the second edge.

Protected Functions

virtual VectorMax3d compute_unnormalized_normal(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Unnormalized normal of the stencil.

virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the Jacobian of the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Jacobian of the unnormalized normal of the stencil.

Friends

template <typename H>
inline friend H AbslHashValue(H hconst EdgeEdgeCandidateee);

Edge-Face Candidate

class EdgeFaceCandidate;

Candidate for intersection between edge and face.

Not included in Candidates because it is not a collision candidate.

Public Functions

EdgeFaceCandidate(index_t edge_idindex_t face_id);
bool operator==(const EdgeFaceCandidateother) const;
bool operator!=(const EdgeFaceCandidateother) const;
bool operator<(const EdgeFaceCandidateother) const;

Compare EdgeFaceCandidate for sorting.

Public Members

index_t edge_id;

ID of the edge.

index_t face_id;

ID of the face.

Friends

template <typename H>
inline friend H AbslHashValue(H hconst EdgeFaceCandidatefv);

Face-Vertex Candidate

class FaceVertexCandidate : public virtual ipc::CollisionStencil;

Inheritance diagram for ipc::FaceVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate" fillcolor="#BFBFBF"] "3" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision"] "4" [label="ipc::FaceVertexTangentialCollision" tooltip="ipc::FaceVertexTangentialCollision"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::FaceVertexCandidate:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

A candidate for face-vertex collision detection.

Subclassed by ipc::FaceVertexNormalCollision, ipc::FaceVertexTangentialCollision

Public Functions

FaceVertexCandidate(index_t face_idindex_t vertex_id);
inline virtual int num_vertices() const override;

Get the number of vertices in the collision stencil.

inline virtual std::array<index_t, 4> vertex_ids(
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;

Get the vertex IDs for the face-vertex pair.

Parameters:
Eigen::ConstRef<Eigen::MatrixXi> edges

The edge connectivity matrix

Eigen::ConstRef<Eigen::MatrixXi> faces

The face connectivity matrix

Returns:

An array of vertex IDs in the order: [vi, f0i, f1i, f2i]

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
    
Eigen::ConstRef<VectorMax12d> vertices_t1double& toi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_phase_ccd
    
= DEFAULT_NARROW_PHASE_CCD
) const override;

Perform narrow-phase CCD on the candidate.

Parameters:
Eigen::ConstRef<VectorMax12d> vertices_t0

[in] Stencil vertices at the start of the time step.

Eigen::ConstRef<VectorMax12d> vertices_t1

[in] Stencil vertices at the end of the time step.

double &toi

[out] Computed time of impact (normalized).

const double min_distance = 0.0

[in] Minimum separation distance between primitives.

const double tmax = 1.0

[in] Maximum time (normalized) to look for collisions.

const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD

[in] The narrow phase CCD algorithm to use.

Returns:

If the candidate had a collision over the time interval.

inline virtual PointTriangleDistanceType known_dtype() const;
bool operator==(const FaceVertexCandidateother) const;
bool operator!=(const FaceVertexCandidateother) const;
bool operator<(const FaceVertexCandidateother) const;

Compare FaceVertexCandidate for sorting.

inline VectorMax4d compute_coefficients(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the coefficients of the stencil s.t.

\(d(x) = \|\sum c_i \mathbf{x}_i\|^2\).

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Coefficients of the stencil.

virtual VectorMax4d compute_coefficients(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

inline double compute_distance(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance of the stencil.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance of the stencil.

virtual double compute_distance(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax12d compute_distance_gradient(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<Eigen::MatrixXd> vertices,
    
Eigen::ConstRef<Eigen::MatrixXi> edges,
    
Eigen::ConstRef<Eigen::MatrixXi> faces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
Eigen::ConstRef<Eigen::MatrixXd> vertices

Collision mesh vertices

Eigen::ConstRef<Eigen::MatrixXi> edges

Collision mesh edges

Eigen::ConstRef<Eigen::MatrixXi> faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
Eigen::ConstRef<VectorMax12d> positions) const
   
 = 0;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

Public Members

index_t face_id;

ID of the face.

index_t vertex_id;

ID of the vertex.

Protected Functions

virtual VectorMax3d compute_unnormalized_normal(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Unnormalized normal of the stencil.

virtual MatrixMax<double, 3, 12>
compute_unnormalized_normal_jacobian(
    
Eigen::ConstRef<VectorMax12d> positions) const override;

Compute the Jacobian of the unnormalized normal of the stencil.

Parameters:
Eigen::ConstRef<VectorMax12d> positions

Stencil’s vertex positions.

Returns:

Jacobian of the unnormalized normal of the stencil.

Friends

template <typename H>
inline friend H AbslHashValue(H hconst FaceVertexCandidatefv);