Candidates¶
Candidates¶
- class Candidates;¶
Public Functions¶
- Candidates() = default;¶
Initialize the set of discrete collision detection candidates.
The surface of the collision mesh.
Surface vertex positions (rowwise).
Amount to inflate the bounding boxes.
- broad_phase_method
Broad phase method to use.
Initialize the set of continuous collision detection candidates.
Note
Assumes the trajectory is linear.
The surface of the collision mesh.
Surface vertex starting positions (rowwise).
Surface vertex ending positions (rowwise).
Amount to inflate the bounding boxes.
- broad_phase_method
Broad phase method to use.
- size_t size() const;¶
- bool empty() const;¶
- void clear();¶
- CollisionStencil& operator[](size_t i);¶
- const CollisionStencil& operator[](size_t i) const;¶
-
bool is_step_collision_free(const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<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.
- 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 CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices_t0,
Eigen::ConstRef<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.
- 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 CollisionMesh& mesh,
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.
Computes a CFL-inspired CCD maximum step step size.
The collision mesh.
Surface vertex starting positions (rowwise).
Surface vertex ending positions (rowwise).
Barrier activation distance.
The minimum distance allowable between any two elements.
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;¶
Inheritance diagram for ipc::CollisionStencil:
A stencil representing a collision between at most four vertices.
Subclassed by ipc::EdgeEdgeCandidate, ipc::EdgeVertexCandidate, ipc::FaceVertexCandidate, ipc::NormalCollision, ipc::TangentialCollision, ipc::VertexVertexCandidate
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(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<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.
-
inline std::array<VectorMax3d, 4> 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.
- Template Parameters:¶
- 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.
-
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.
-
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.
-
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.
-
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.
-
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) = ‖∑ cᵢ xᵢ‖².
-
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)
-
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)
-
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)
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<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:¶
- 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::ostream& write_ccd_query(std::ostream& out,
Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<VectorMax12d> vertices_t1) const;¶ Write the CCD query to a stream.
Vertex-Vertex Candidate¶
- class VertexVertexCandidate : public virtual ipc::CollisionStencil;¶
Inheritance 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(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<Eigen::MatrixXi> faces) const override;¶ Get the indices of the vertices.
-
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)
-
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)
-
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)
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<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:¶
- 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 VertexVertexCandidate& other) const;¶
- bool operator!=(const VertexVertexCandidate& other) const;¶
- bool operator<(const VertexVertexCandidate& other) 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) = ‖∑ cᵢ xᵢ‖².
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
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.
-
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)
-
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.
-
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)
-
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.
-
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)
Edge-Vertex Candidate¶
- class EdgeVertexCandidate : public virtual ipc::CollisionStencil;¶
Inheritance 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(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<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(
Eigen::ConstRef<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(
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)
-
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)
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<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:¶
- 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 EdgeVertexCandidate& other) const;¶
- bool operator!=(const EdgeVertexCandidate& other) const;¶
- bool operator<(const EdgeVertexCandidate& other) 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) = ‖∑ cᵢ xᵢ‖².
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
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.
-
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)
-
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.
-
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)
-
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.
-
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)
Edge-Edge Candidate¶
- class EdgeEdgeCandidate : public virtual ipc::CollisionStencil;¶
Inheritance 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(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<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(
Eigen::ConstRef<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(
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)
-
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)
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<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:¶
- 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 EdgeEdgeCandidate& other) const;¶
- bool operator!=(const EdgeEdgeCandidate& other) const;¶
- bool operator<(const EdgeEdgeCandidate& other) 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) = ‖∑ cᵢ xᵢ‖².
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
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.
-
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)
-
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.
-
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)
-
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.
-
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)
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 virtual ipc::CollisionStencil;¶
Inheritance 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(
Eigen::ConstRef<Eigen::MatrixXi> edges,
Eigen::ConstRef<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(
Eigen::ConstRef<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(
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)
-
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)
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual bool ccd(Eigen::ConstRef<VectorMax12d> vertices_t0,
Eigen::ConstRef<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:¶
- 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 FaceVertexCandidate& other) const;¶
- bool operator!=(const FaceVertexCandidate& other) const;¶
- bool operator<(const FaceVertexCandidate& other) 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) = ‖∑ cᵢ xᵢ‖².
-
virtual VectorMax4d compute_coefficients(
Eigen::ConstRef<VectorMax12d> positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
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.
-
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)
-
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.
-
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)
-
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.
-
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)