Candidates¶
Candidates¶
- class Candidates¶
Public Functions
- Candidates() = default¶
-
void build(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices,
const double inflation_radius = 0,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD)¶ Initialize the set of discrete collision detection candidates.
- Parameters:¶
- const CollisionMesh &mesh¶
The surface of the collision mesh.
- const Eigen::MatrixXd &vertices¶
Surface vertex positions (rowwise).
- const double inflation_radius = 0¶
Amount to inflate the bounding boxes.
- const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD¶
Broad phase method to use.
-
void build(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices_t0,
const Eigen::MatrixXd& vertices_t1,
const double inflation_radius = 0,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD)¶ Initialize the set of continuous collision detection candidates.
Note
Assumes the trajectory is linear.
- Parameters:¶
- const CollisionMesh &mesh¶
The surface of the collision mesh.
- const Eigen::MatrixXd &vertices_t0¶
Surface vertex starting positions (rowwise).
- const Eigen::MatrixXd &vertices_t1¶
Surface vertex ending positions (rowwise).
- const double inflation_radius = 0¶
Amount to inflate the bounding boxes.
- const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD¶
Broad phase method to use.
- size_t size() const¶
- bool empty() const¶
- void clear()¶
- ContinuousCollisionCandidate& operator[](size_t i)¶
- const ContinuousCollisionCandidate& operator[](size_t i) const¶
-
bool is_step_collision_free(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices_t0,
const Eigen::MatrixXd& vertices_t1,
const double min_distance = 0.0,
const NarrowPhaseCCD& narrow_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.
- const Eigen::MatrixXd &vertices_t0¶
Surface vertex starting positions (rowwise).
- const 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 CollisionMesh& mesh,
const Eigen::MatrixXd& vertices_t0,
const Eigen::MatrixXd& vertices_t1,
const double min_distance = 0.0,
const NarrowPhaseCCD& narrow_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.
- const Eigen::MatrixXd &vertices_t0¶
Surface vertex starting positions (rowwise). Assumed to be intersection free.
- const 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 CollisionMesh& mesh, const Eigen::MatrixXd& displacements,
const double dhat) const¶ Computes a conservative bound on the largest-feasible step size for surface primitives not in collision.
-
double compute_cfl_stepsize(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices_t0,
const Eigen::MatrixXd& vertices_t1, const double dhat,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD,
const double min_distance = 0.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const¶ Computes a CFL-inspired CCD maximum step step size.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices_t0¶
Surface vertex starting positions (rowwise).
- const 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.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
Public Members
- std::vector<VertexVertexCandidate> vv_candidates¶
- std::vector<EdgeVertexCandidate> ev_candidates¶
- std::vector<EdgeEdgeCandidate> ee_candidates¶
- std::vector<FaceVertexCandidate> fv_candidates¶
Collision Stencil¶
- class CollisionStencil¶
Inheritence diagram for ipc::CollisionStencil:
A stencil representing a collision between at most four vertices.
Subclassed by ipc::ContinuousCollisionCandidate, ipc::NormalCollision, ipc::TangentialCollision
Public Functions
- virtual ~CollisionStencil() = default¶
- virtual int num_vertices() const = 0¶
Get the number of vertices in the collision stencil.
-
virtual std::array<long, 4> vertex_ids(const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const = 0¶ Get the vertex IDs of the collision stencil.
- Parameters:¶
- Returns:¶
The vertex IDs of the collision stencil. Size is always 4, but elements i > num_vertices() are -1.
-
template <typename T>
inline std::array<VectorMax3<T>, 4> vertices(
const MatrixX<T>& vertices, const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const¶ Get the vertex attributes of the collision stencil.
- Template Parameters:¶
- typename T¶
Type of the attributes
- Parameters:¶
- Returns:¶
The vertex positions of the collision stencil. Size is always 4, but elements i > num_vertices() are NaN.
-
template <typename T>
inline VectorMax12<T> dof(const MatrixX<T>& X,
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const¶ Select this stencil’s DOF from the full matrix of DOF.
-
virtual double compute_distance(
const VectorMax12d& positions) const = 0¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
-
virtual VectorMax12d compute_distance_gradient(
const 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)
Continuous Collision Candidate¶
-
class ContinuousCollisionCandidate
: public virtual ipc::CollisionStencil¶ Inheritence diagram for ipc::ContinuousCollisionCandidate:
Collaboration diagram for ipc::ContinuousCollisionCandidate:
Virtual class for candidates that support CCD.
Subclassed by ipc::EdgeEdgeCandidate, ipc::EdgeVertexCandidate, ipc::FaceVertexCandidate, ipc::VertexVertexCandidate
Public Functions
- virtual ~ContinuousCollisionCandidate() = default¶
-
virtual bool ccd(const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const = 0¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- const VectorMax12d &vertices_t0¶
[in] Stencil vertices at the start of the time step.
- const 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::ostream& write_ccd_query(std::ostream& out,
const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1) const¶ Write the CCD query to a stream.
Vertex-Vertex Candidate¶
-
class VertexVertexCandidate
: public ipc::ContinuousCollisionCandidate¶ Inheritence diagram for ipc::VertexVertexCandidate:
Collaboration diagram for ipc::VertexVertexCandidate:
Subclassed by ipc::VertexVertexNormalCollision, ipc::VertexVertexTangentialCollision
Public Functions
- VertexVertexCandidate(long vertex0_id, long vertex1_id)¶
- inline virtual int num_vertices() const override¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<long, 4> vertex_ids(
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const override¶ Get the indices of the vertices.
-
virtual double compute_distance(
const VectorMax12d& positions) const override¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
-
virtual VectorMax12d compute_distance_gradient(
const 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)
-
virtual MatrixMax12d compute_distance_hessian(
const 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)
-
virtual bool ccd(const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- const VectorMax12d &vertices_t0¶
[in] Stencil vertices at the start of the time step.
- const 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 VertexVertexCandidate& other) const¶
- bool operator!=(const VertexVertexCandidate& other) const¶
- bool operator<(const VertexVertexCandidate& other) const¶
Compare EdgeVertexCandidates for sorting.
Edge-Vertex Candidate¶
-
class EdgeVertexCandidate
: public ipc::ContinuousCollisionCandidate¶ Inheritence diagram for ipc::EdgeVertexCandidate:
Collaboration diagram for ipc::EdgeVertexCandidate:
Subclassed by ipc::EdgeVertexNormalCollision, ipc::EdgeVertexTangentialCollision
Public Functions
- inline virtual int num_vertices() const override¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<long, 4> vertex_ids(
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const override¶ Get the vertex IDs of the collision stencil.
- Parameters:¶
- Returns:¶
The vertex IDs of the collision stencil. Size is always 4, but elements i > num_vertices() are -1.
-
virtual double compute_distance(
const VectorMax12d& positions) const override¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
-
virtual VectorMax12d compute_distance_gradient(
const 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)
-
virtual MatrixMax12d compute_distance_hessian(
const 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)
-
virtual bool ccd(const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- const VectorMax12d &vertices_t0¶
[in] Stencil vertices at the start of the time step.
- const 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 EdgeVertexCandidate& other) const¶
- bool operator!=(const EdgeVertexCandidate& other) const¶
- bool operator<(const EdgeVertexCandidate& other) const¶
Compare EdgeVertexCandidates for sorting.
Edge-Edge Candidate¶
- class EdgeEdgeCandidate : public ipc::ContinuousCollisionCandidate¶
Inheritence diagram for ipc::EdgeEdgeCandidate:
Collaboration diagram for ipc::EdgeEdgeCandidate:
Subclassed by ipc::EdgeEdgeNormalCollision, ipc::EdgeEdgeTangentialCollision
Public Functions
- inline virtual int num_vertices() const override¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<long, 4> vertex_ids(
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const override¶ Get the vertex IDs of the collision stencil.
- Parameters:¶
- Returns:¶
The vertex IDs of the collision stencil. Size is always 4, but elements i > num_vertices() are -1.
-
virtual double compute_distance(
const VectorMax12d& positions) const override¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
-
virtual VectorMax12d compute_distance_gradient(
const 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)
-
virtual MatrixMax12d compute_distance_hessian(
const 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)
-
virtual bool ccd(const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- const VectorMax12d &vertices_t0¶
[in] Stencil vertices at the start of the time step.
- const 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 EdgeEdgeCandidate& other) const¶
- bool operator!=(const EdgeEdgeCandidate& other) const¶
- bool operator<(const EdgeEdgeCandidate& other) const¶
Compare EdgeEdgeCandidates for sorting.
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
- bool operator==(const EdgeFaceCandidate& other) const¶
- bool operator!=(const EdgeFaceCandidate& other) const¶
- bool operator<(const EdgeFaceCandidate& other) const¶
Compare EdgeFaceCandidate for sorting.
Face-Vertex Candidate¶
-
class FaceVertexCandidate
: public ipc::ContinuousCollisionCandidate¶ Inheritence diagram for ipc::FaceVertexCandidate:
Collaboration diagram for ipc::FaceVertexCandidate:
Subclassed by ipc::FaceVertexNormalCollision, ipc::FaceVertexTangentialCollision
Public Functions
- inline virtual int num_vertices() const override¶
Get the number of vertices in the collision stencil.
-
inline virtual std::array<long, 4> vertex_ids(
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const override¶ Get the vertex IDs of the collision stencil.
- Parameters:¶
- Returns:¶
The vertex IDs of the collision stencil. Size is always 4, but elements i > num_vertices() are -1.
-
virtual double compute_distance(
const VectorMax12d& positions) const override¶ Compute the distance of the stencil.
Note
positions can be computed as stencil.dof(vertices, edges, faces)
-
virtual VectorMax12d compute_distance_gradient(
const 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)
-
virtual MatrixMax12d compute_distance_hessian(
const 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)
-
virtual bool ccd(const VectorMax12d& vertices_t0,
const VectorMax12d& vertices_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD) const override¶ Perform narrow-phase CCD on the candidate.
- Parameters:¶
- const VectorMax12d &vertices_t0¶
[in] Stencil vertices at the start of the time step.
- const 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 FaceVertexCandidate& other) const¶
- bool operator!=(const FaceVertexCandidate& other) const¶
- bool operator<(const FaceVertexCandidate& other) const¶
Compare FaceVertexCandidate for sorting.