Collision Constraints¶
Collision Constraints¶
- class CollisionConstraints;¶
Public Functions¶
- inline CollisionConstraints();¶
-
void build(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices, const double dhat,
const double dmin = 0,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD);¶ Initialize the set of constraints used to compute the barrier potential.
- Parameters:¶
- 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.
- const BroadPhaseMethod broad_phase_method = DEFAULT_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 constraints used to compute the barrier potential.
- Parameters:¶
- const Candidates &candidates¶
Distance candidates from which the constraint 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_potential(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices, const double dhat) const;¶ Compute the barrier potential for a given constraint set.
-
Eigen::VectorXd compute_potential_gradient(
const CollisionMesh& mesh, const Eigen::MatrixXd& vertices,
const double dhat) const;¶ Compute the gradient of the barrier potential.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices¶
Vertices of the collision mesh.
- const double dhat¶
The activation distance of the barrier.
- Returns:¶
The gradient of all barrier potentials (not scaled by the barrier stiffness). This will have a size of |vertices|.
-
Eigen::SparseMatrix<double> compute_potential_hessian(
const CollisionMesh& mesh, const Eigen::MatrixXd& vertices,
const double dhat,
const bool project_hessian_to_psd = false) const;¶ Compute the hessian of the barrier potential.
- Parameters:¶
- 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 bool project_hessian_to_psd = false¶
Make sure the hessian is positive semi-definite.
- Returns:¶
The hessian of all barrier potentials (not scaled by the barrier stiffness). This will have a size of |vertices|x|vertices|.
-
Eigen::SparseMatrix<double> compute_shape_derivative(
const CollisionMesh& mesh, const Eigen::MatrixXd& vertices,
const double dhat) const;¶ Compute the barrier shape derivative.
- Parameters:¶
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices¶
Vertices of the collision mesh.
- const double dhat¶
The activation distance of the barrier.
- Throws:¶
std::runtime_error
– If the collision constraints were not built with shape derivatives enabled.- Returns:¶
The derivative of the force with respect to X, the rest vertices.
-
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 collision constraints.
- bool empty() const;¶
Get if the collision constraints are empty.
- void clear();¶
Clear the collision constraints.
- CollisionConstraint& operator[](size_t idx);¶
Get a reference to constriant idx.
- const CollisionConstraint& operator[](size_t idx) const;¶
Get a const reference to constriant idx.
- inline bool use_convergent_formulation() const;¶
Get if the collision constraints should use the convergent formulation.
Note
If not empty, this is the current value not necessarily the value used to build the constraints.
- Returns:¶
If the collision constraints should use the convergent formulation.
-
void set_use_convergent_formulation(
const bool use_convergent_formulation);¶ Set if the collision constraints should use the convergent formulation.
Warning
This must be set before the constraints are built.
- inline bool are_shape_derivatives_enabled() const;¶
Get if the collision constraints are using the convergent formulation.
Note
If not empty, this is the current value not necessarily the value used to build the constraints.
- Returns:¶
If the collision constraints are using the convergent formulation.
-
void set_are_shape_derivatives_enabled(
const bool are_shape_derivatives_enabled);¶ Set if the collision constraints should enable shape derivative computation.
Warning
This must be set before the constraints are built.
Public Members¶
- std::vector<VertexVertexConstraint> vv_constraints;¶
- std::vector<EdgeVertexConstraint> ev_constraints;¶
- std::vector<EdgeEdgeConstraint> ee_constraints;¶
- std::vector<FaceVertexConstraint> fv_constraints;¶
- std::vector<PlaneVertexConstraint> pv_constraints;¶
Collision Constraint¶
- class CollisionConstraint : public virtual ipc::CollisionStencil;¶
Inheritance diagram for ipc::CollisionConstraint:
Collaboration diagram for ipc::CollisionConstraint:
Subclassed by ipc::EdgeEdgeConstraint, ipc::EdgeVertexConstraint, ipc::FaceVertexConstraint, ipc::PlaneVertexConstraint, ipc::VertexVertexConstraint
Public Functions¶
- inline virtual ~CollisionConstraint();¶
-
virtual double compute_potential(const Eigen::MatrixXd& vertices,
const Eigen::MatrixXi& edges, const Eigen::MatrixXi& faces,
const double dhat) const;¶
Vertex-Vertex Collision Constraint¶
-
class VertexVertexConstraint : public ipc::VertexVertexCandidate,
public ipc::CollisionConstraint;¶ Inheritance diagram for ipc::VertexVertexConstraint:
Collaboration diagram for ipc::VertexVertexConstraint:
Public Functions¶
-
inline VertexVertexConstraint(
const VertexVertexCandidate& candidate);¶
- VertexVertexCandidate(long vertex0_id, long vertex1_id);¶
-
inline VertexVertexConstraint(
Edge-Vertex Collision Constraint¶
-
class EdgeVertexConstraint : public ipc::EdgeVertexCandidate,
public ipc::CollisionConstraint;¶ Inheritance diagram for ipc::EdgeVertexConstraint:
Collaboration diagram for ipc::EdgeVertexConstraint:
Public Functions¶
- inline EdgeVertexConstraint(const EdgeVertexCandidate& candidate);¶
Protected Functions¶
- inline virtual PointEdgeDistanceType known_dtype() const override;¶
Edge-Edge Collision Constraint¶
-
class EdgeEdgeConstraint : public ipc::EdgeEdgeCandidate,
public ipc::CollisionConstraint;¶ Inheritance diagram for ipc::EdgeEdgeConstraint:
Collaboration diagram for ipc::EdgeEdgeConstraint:
Public Functions¶
-
EdgeEdgeConstraint(
const EdgeEdgeCandidate& candidate, double eps_x);¶
-
virtual double compute_potential(const Eigen::MatrixXd& vertices,
const Eigen::MatrixXi& edges, const Eigen::MatrixXi& faces,
const double dhat) const override;¶
-
EdgeEdgeConstraint(
Face-Vertex Collision Constraint¶
-
class FaceVertexConstraint : public ipc::FaceVertexCandidate,
public ipc::CollisionConstraint;¶ Inheritance diagram for ipc::FaceVertexConstraint:
Collaboration diagram for ipc::FaceVertexConstraint:
Public Functions¶
- inline FaceVertexConstraint(const FaceVertexCandidate& candidate);¶
Protected Functions¶
-
inline virtual PointTriangleDistanceType
known_dtype() const override;¶
Plane-Vertex Collision Constraint¶
- class PlaneVertexConstraint : public ipc::CollisionConstraint;¶
Inheritance diagram for ipc::PlaneVertexConstraint:
Collaboration diagram for ipc::PlaneVertexConstraint:
Public Functions¶
-
PlaneVertexConstraint(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 contact stencil.
Protected Functions¶
-
virtual double compute_distance(
const VectorMax12d& point) const override;¶ Compute the distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
const VectorMax12d& point) const override;¶ Compute the distance gradient of the stencil w.r.t.
the stencil’s vertex positions.
-
PlaneVertexConstraint(const VectorMax3d& plane_origin,