Normal Collisions¶
Normal Collisions¶
- class NormalCollisions;¶
Public Types¶
- using value_type = NormalCollision;¶
The type of the collisions.
Public Functions¶
- NormalCollisions() = default;¶
Initialize the set of collisions used to compute the barrier potential.
The collision mesh.
Vertices of the collision mesh.
The activation distance of the barrier.
Minimum distance.
- broad_phase_method
Broad-phase method to use.
-
void build(const Candidates& candidates, const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices, const double dhat,
const double dmin = 0);¶ Initialize the set of collisions used to compute the barrier potential.
- Parameters:¶
- const Candidates &candidates¶
Distance candidates from which the collision set is built.
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices¶
Vertices of the collision mesh.
- const double dhat¶
The activation distance of the barrier.
- const double dmin = 0¶
Minimum distance.
-
double compute_minimum_distance(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices) const;¶ Computes the minimum distance between any non-adjacent elements.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices¶
Vertices of the collision mesh.
- Returns:¶
The minimum distance between any non-adjacent elements.
- size_t size() const;¶
Get the number of collisions.
- bool empty() const;¶
Get if the collision set are empty.
- void clear();¶
Clear the collision set.
- NormalCollision& operator[](size_t i);¶
Get a reference to collision at index i.
- const NormalCollision& operator[](size_t i) const;¶
Get a const reference to collision at index i.
- inline bool use_area_weighting() const;¶
Get if the collision set should use area weighting.
Note
If not empty, this is the current value not necessarily the value used to build the collisions.
- Returns:¶
If the collision set should use area weighting.
- void set_use_area_weighting(const bool use_area_weighting);¶
Set if the collision set should use area weighting.
Warning
This must be set before the collisions are built.
- inline bool use_improved_max_approximator() const;¶
Get if the collision set should use the improved max approximator.
Note
If not empty, this is the current value not necessarily the value used to build the collisions.
- Returns:¶
If the collision set should use the improved max approximator.
-
void set_use_improved_max_approximator(
const bool use_improved_max_approximator);¶ Set if the collision set should use the improved max approximator.
Warning
This must be set before the collisions are built.
- inline bool enable_shape_derivatives() const;¶
Get if the collision set are using the convergent formulation.
Note
If not empty, this is the current value not necessarily the value used to build the collisions.
- Returns:¶
If the collision set are using the convergent formulation.
-
void set_enable_shape_derivatives(
const bool enable_shape_derivatives);¶ Set if the collision set should enable shape derivative computation.
Warning
This must be set before the collisions are built.
-
std::string to_string(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices) const;¶
Public Members¶
- std::vector<VertexVertexNormalCollision> vv_collisions;¶
Vertex-vertex normal collisions.
- std::vector<EdgeVertexNormalCollision> ev_collisions;¶
Edge-vertex normal collisions.
- std::vector<EdgeEdgeNormalCollision> ee_collisions;¶
Edge-edge normal collisions.
- std::vector<FaceVertexNormalCollision> fv_collisions;¶
Face-vertex normal collisions.
- std::vector<PlaneVertexNormalCollision> pv_collisions;¶
Plane-vertex normal collisions.
Normal Collision¶
- class NormalCollision : public virtual ipc::CollisionStencil;¶
Inheritence diagram for ipc::NormalCollision:
Collaboration diagram for ipc::NormalCollision:
Subclassed by ipc::EdgeEdgeNormalCollision, ipc::EdgeVertexNormalCollision, ipc::FaceVertexNormalCollision, ipc::PlaneVertexNormalCollision, ipc::VertexVertexNormalCollision
Public Functions¶
- NormalCollision() = default;¶
-
NormalCollision(const double weight,
const Eigen::SparseVector<double>& weight_gradient);¶
- virtual ~NormalCollision() = default;¶
- inline virtual bool is_mollified() const;¶
Does the distance potentially have to be mollified?
-
inline virtual double mollifier_threshold(
const VectorMax12d& rest_positions) const;¶ Compute the mollifier threshold for the distance.
- virtual double mollifier(const VectorMax12d& positions) const;¶
Compute the mollifier for the distance.
-
virtual double mollifier(
const VectorMax12d& positions, double eps_x) const;¶ Compute the mollifier for the distance.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions) const;¶ Compute the gradient of the mollifier for the distance wrt the positions.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions, double eps_x) const;¶ Compute the gradient of the mollifier for the distance wrt the positions.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions) const;¶ Compute the Hessian of the mollifier for the distance wrt the positions.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions, double eps_x) const;¶ Compute the Hessian of the mollifier for the distance wrt the positions.
-
virtual Vector12d mollifier_gradient_wrt_x(
const VectorMax12d& rest_positions,
const VectorMax12d& positions) const;¶ Compute the gradient of the mollifier for the distance w.r.t.
rest positions.
-
virtual Matrix12d mollifier_gradient_jacobian_wrt_x(
const VectorMax12d& rest_positions,
const VectorMax12d& positions) const;¶ Compute the jacobian of the distance mollifier’s gradient w.r.t.
rest positions.
Vertex-Vertex Normal Collision¶
-
class VertexVertexNormalCollision
: public ipc::VertexVertexCandidate,
public ipc::NormalCollision;¶ Inheritence diagram for ipc::VertexVertexNormalCollision:
Collaboration diagram for ipc::VertexVertexNormalCollision:
Public Functions¶
-
inline VertexVertexNormalCollision(
const VertexVertexCandidate& candidate);¶
-
inline VertexVertexNormalCollision(const long _vertex0_id,
const long _vertex1_id, const double _weight,
const Eigen::SparseVector<double>& _weight_gradient);¶
- VertexVertexCandidate(long vertex0_id, long vertex1_id);¶
-
inline VertexVertexNormalCollision(
Edge-Vertex Normal Collision¶
-
class EdgeVertexNormalCollision : public ipc::EdgeVertexCandidate,
public ipc::NormalCollision;¶ Inheritence diagram for ipc::EdgeVertexNormalCollision:
Collaboration diagram for ipc::EdgeVertexNormalCollision:
Public Functions¶
-
inline EdgeVertexNormalCollision(
const EdgeVertexCandidate& candidate);¶
-
inline EdgeVertexNormalCollision(const long _edge_id,
const long _vertex_id, const double _weight,
const Eigen::SparseVector<double>& _weight_gradient);¶
- inline virtual PointEdgeDistanceType known_dtype() const override;¶
-
inline EdgeVertexNormalCollision(
Edge-Edge Normal Collision¶
-
class EdgeEdgeNormalCollision : public ipc::EdgeEdgeCandidate,
public ipc::NormalCollision;¶ Inheritence diagram for ipc::EdgeEdgeNormalCollision:
Collaboration diagram for ipc::EdgeEdgeNormalCollision:
Public Functions¶
-
EdgeEdgeNormalCollision(const long edge0_id, const long edge1_id,
const double eps_x,
const EdgeEdgeDistanceType dtype = EdgeEdgeDistanceType::AUTO);¶
-
EdgeEdgeNormalCollision(const EdgeEdgeCandidate& candidate,
const double eps_x,
const EdgeEdgeDistanceType dtype = EdgeEdgeDistanceType::AUTO);¶
-
EdgeEdgeNormalCollision(const long edge0_id, const long edge1_id,
const double eps_x, const double weight,
const Eigen::SparseVector<double>& weight_gradient,
const EdgeEdgeDistanceType dtype = EdgeEdgeDistanceType::AUTO);¶
- inline virtual bool is_mollified() const override;¶
Does the distance potentially have to be mollified?
-
virtual double mollifier_threshold(
const VectorMax12d& rest_positions) const override;¶ Compute the mollifier threshold for the distance.
-
virtual double mollifier(
const VectorMax12d& positions) const override;¶ Compute the mollifier for the distance.
-
virtual double mollifier(
const VectorMax12d& positions, double eps_x) const override;¶ Compute the mollifier for the distance.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions) const override;¶ Compute the gradient of the mollifier for the distance w.r.t.
positions.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions, double eps_x) const override;¶ Compute the gradient of the mollifier for the distance wrt the positions.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions) const override;¶ Compute the Hessian of the mollifier for the distance w.r.t.
positions.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions, double eps_x) const override;¶ Compute the Hessian of the mollifier for the distance wrt the positions.
-
virtual Vector12d mollifier_gradient_wrt_x(
const VectorMax12d& rest_positions,
const VectorMax12d& positions) const override;¶ Compute the gradient of the mollifier for the distance w.r.t.
rest positions.
-
virtual Matrix12d mollifier_gradient_jacobian_wrt_x(
const VectorMax12d& rest_positions,
const VectorMax12d& positions) const override;¶ Compute the jacobian of the distance mollifier’s gradient w.r.t.
rest positions.
- inline virtual EdgeEdgeDistanceType known_dtype() const override;¶
- bool operator==(const EdgeEdgeNormalCollision& other) const;¶
- bool operator!=(const EdgeEdgeNormalCollision& other) const;¶
- bool operator<(const EdgeEdgeNormalCollision& other) const;¶
Public Members¶
- double eps_x;¶
Mollifier activation threshold.
See also
- EdgeEdgeDistanceType dtype;¶
Cached distance type.
Some EE collisions are mollified EV or VV collisions.
-
EdgeEdgeNormalCollision(const long edge0_id, const long edge1_id,
Face-Vertex Normal Collision¶
-
class FaceVertexNormalCollision : public ipc::FaceVertexCandidate,
public ipc::NormalCollision;¶ Inheritence diagram for ipc::FaceVertexNormalCollision:
Collaboration diagram for ipc::FaceVertexNormalCollision:
Public Functions¶
-
inline FaceVertexNormalCollision(
const FaceVertexCandidate& candidate);¶
-
inline FaceVertexNormalCollision(const long _face_id,
const long _vertex_id, const double _weight,
const Eigen::SparseVector<double>& _weight_gradient);¶
-
inline virtual PointTriangleDistanceType
known_dtype() const override;¶
-
inline FaceVertexNormalCollision(
Plane-Vertex Normal Collision¶
- class PlaneVertexNormalCollision : public ipc::NormalCollision;¶
Inheritence diagram for ipc::PlaneVertexNormalCollision:
Collaboration diagram for ipc::PlaneVertexNormalCollision:
Public Functions¶
-
PlaneVertexNormalCollision(const VectorMax3d& plane_origin,
const VectorMax3d& plane_normal, const long vertex_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 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& point) const override;¶ Compute the distance between the point and plane.
-
virtual VectorMax12d compute_distance_gradient(
const VectorMax12d& point) const override;¶ Compute the gradient of the distance w.r.t.
the point’s positions.
-
virtual MatrixMax12d compute_distance_hessian(
const VectorMax12d& point) const override;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual VectorMax4d compute_coefficients(
const VectorMax12d& positions) const override;¶ Compute the coefficients of the stencil.
-
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 VectorMax4d compute_coefficients(
const Eigen::MatrixXd& vertices, const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const;¶ Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
virtual VectorMax4d compute_coefficients(
const VectorMax12d& positions) const
= 0; Compute the coefficients of the stencil s.t.
d(x) = ‖∑ cᵢ xᵢ‖².
-
inline double compute_distance(const Eigen::MatrixXd& vertices,
const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const;¶ Compute the distance of the stencil.
-
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)
-
inline VectorMax12d compute_distance_gradient(
const Eigen::MatrixXd& vertices, const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
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)
-
inline MatrixMax12d compute_distance_hessian(
const Eigen::MatrixXd& vertices, const Eigen::MatrixXi& edges,
const Eigen::MatrixXi& faces) const;¶ Compute the distance Hessian of the stencil w.r.t.
the stencil’s vertex positions.
-
virtual MatrixMax12d compute_distance_hessian(
const 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)
-
PlaneVertexNormalCollision(const VectorMax3d& plane_origin,