Normal Collisions
Normal Collisions
-
class NormalCollisions
Public Types
-
using value_type = NormalCollision
The type of the collisions.
Public Functions
-
NormalCollisions() = default
-
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 collisions 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 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.
- Parameters:
- size_t i
The index of the collision.
- Returns:
A reference to the collision.
-
const NormalCollision& operator[](size_t i) const
Get a const reference to collision at index i.
- Parameters:
- size_t i
The index of the collision.
- Returns:
A const reference to the collision.
-
bool is_vertex_vertex(size_t i) const
Get if the collision at i is a vertex-vertex collision.
- Parameters:
- size_t i
The index of the collision.
- Returns:
If the collision at i is a vertex-vertex collision.
-
bool is_edge_vertex(size_t i) const
Get if the collision at i is an edge-vertex collision.
- Parameters:
- size_t i
The index of the collision.
- Returns:
If the collision at i is an edge-vertex collision.
-
bool is_edge_edge(size_t i) const
Get if the collision at i is an edge-edge collision.
- Parameters:
- size_t i
The index of the collision.
- Returns:
If the collision at i is an edge-edge collision.
-
bool is_face_vertex(size_t i) const
Get if the collision at i is an face-vertex collision.
- Parameters:
- size_t i
The index of the collision.
- Returns:
If the collision at i is an face-vertex collision.
-
bool is_plane_vertex(size_t i) const
Get if the collision at i is an plane-vertex collision.
- Parameters:
- size_t i
The index of the collision.
- Returns:
If the collision at i is an plane-vertex collision.
-
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.
- Parameters:
- const bool use_area_weighting
If the collision set should use area weighting.
-
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.
- Parameters:
- const bool use_improved_max_approximator
If the collision set should use the improved max approximator.
-
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.
- Parameters:
- const bool enable_shape_derivatives
If the collision set should enable shape derivative computation.
-
std::string to_string(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices) const
Protected Attributes
-
bool m_use_area_weighting = false
-
bool m_use_improved_max_approximator = false
-
bool m_enable_shape_derivatives = false
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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- Returns:
The mollifier threshold.
-
virtual double mollifier(const VectorMax12d& positions) const
Compute the mollifier for the distance.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier value.
-
virtual double mollifier(
const VectorMax12d& positions, double eps_x) const
Compute the mollifier for the distance.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s threshold.
- Returns:
The mollifier value.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions) const
Compute the gradient of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier gradient.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions, double eps_x) const
Compute the gradient of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s threshold.
- Returns:
The mollifier gradient.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions) const
Compute the Hessian of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier Hessian.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions, double eps_x) const
Compute the Hessian of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s threshold.
- Returns:
The mollifier Hessian.
-
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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier gradient 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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The jacobian of the mollifier’s gradient w.r.t. rest positions.
Public Members
-
double dmin = 0
The minimum separation distance.
-
double weight = 1
The term’s weight (e.g., collision area)
-
Eigen::SparseVector<double> weight_gradient
The gradient of the term’s weight wrt the rest positions.
Vertex-Vertex Normal Collision
-
class VertexVertexNormalCollision
: public ipc::VertexVertexCandidate,
public ipc::NormalCollision
Inheritence diagram for ipc::VertexVertexNormalCollision:
Collaboration diagram for ipc::VertexVertexNormalCollision:
Edge-Vertex Normal Collision
-
class EdgeVertexNormalCollision : public ipc::EdgeVertexCandidate,
public ipc::NormalCollision
Inheritence diagram for ipc::EdgeVertexNormalCollision:
Collaboration diagram for ipc::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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- Returns:
The mollifier threshold.
-
virtual double mollifier(
const VectorMax12d& positions) const override
Compute the mollifier for the distance.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier value.
-
virtual double mollifier(
const VectorMax12d& positions, double eps_x) const override
Compute the mollifier for the distance.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s tolerance.
- Returns:
The mollifier value.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions) const override
Compute the gradient of the mollifier for the distance w.r.t.
positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier gradient.
-
virtual VectorMax12d mollifier_gradient(
const VectorMax12d& positions, double eps_x) const override
Compute the gradient of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s tolerance.
- Returns:
The mollifier gradient.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions) const override
Compute the Hessian of the mollifier for the distance w.r.t.
positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier Hessian.
-
virtual MatrixMax12d mollifier_hessian(
const VectorMax12d& positions, double eps_x) const override
Compute the Hessian of the mollifier for the distance wrt the positions.
- Parameters:
- const VectorMax12d &positions
The stencil’s vertex positions.
- double eps_x
The mollifier’s tolerance.
- Returns:
The mollifier Hessian.
-
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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The mollifier gradient 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.
- Parameters:
- const VectorMax12d &rest_positions
The stencil’s rest vertex positions.
- const VectorMax12d &positions
The stencil’s vertex positions.
- Returns:
The jacobian of the 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.
-
EdgeEdgeDistanceType dtype
Cached distance type.
Some EE collisions are mollified EV or VV collisions.
Face-Vertex Normal Collision
-
class FaceVertexNormalCollision : public ipc::FaceVertexCandidate,
public ipc::NormalCollision
Inheritence diagram for ipc::FaceVertexNormalCollision:
Collaboration diagram for ipc::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:
- const Eigen::MatrixXi &edges
Collision mesh edges
- const Eigen::MatrixXi &faces
Collision mesh faces
- 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.
- Parameters:
- const VectorMax12d &point
Point’s position.
- Returns:
Distance of the stencil.
-
virtual VectorMax12d compute_distance_gradient(
const VectorMax12d& point) const override
Compute the gradient of the distance w.r.t.
the point’s positions.
- Parameters:
- const VectorMax12d &point
Point’s position.
- Returns:
Distance gradient 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.
- Parameters:
- const VectorMax12d &point
Point’s position.
- Returns:
Distance Hessian w.r.t. the point’s positions.
Public Members
-
VectorMax3d plane_origin
The plane’s origin.
-
VectorMax3d plane_normal
The plane’s normal.
-
long vertex_id
The vertex’s id.
Last update:
Nov 16, 2024