Normal Collisions

Normal Collisions

class NormalCollisions;

Public Types

using value_type = NormalCollision;

The type of the collisions.

Public Functions

NormalCollisions() = default;
void build(const CollisionMeshmesh,
    
const Eigen::MatrixXdverticesconst double dhat,
    
const double dmin = 0,
    
const std::shared_ptr<BroadPhase> broad_phase
   
 = make_default_broad_phase()
);

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.

broad_phase_method

Broad-phase method to use.

void build(const Candidatescandidatesconst CollisionMeshmesh,
    
const Eigen::MatrixXdverticesconst 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 CollisionMeshmesh,
    
const Eigen::MatrixXdvertices) 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.

NormalCollisionoperator[](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 NormalCollisionoperator[](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 CollisionMeshmesh,
    
const Eigen::MatrixXdvertices) 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.

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:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "3" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision"] "4" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision"] "5" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision"] "1" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision" fillcolor="#BFBFBF"] "6" [label="ipc::PlaneVertexNormalCollision" tooltip="ipc::PlaneVertexNormalCollision"] "7" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] "5" -> "1" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "6" -> "1" [dir=forward tooltip="public-inheritance"] "7" -> "1" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::NormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

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 VectorMax12drest_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 VectorMax12dpositions) const;

Compute the mollifier for the distance.

Parameters:
const VectorMax12d &positions

The stencil’s vertex positions.

Returns:

The mollifier value.

virtual double mollifier(
    
const VectorMax12dpositionsdouble 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 VectorMax12dpositions) 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 VectorMax12dpositionsdouble 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 VectorMax12dpositions) 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 VectorMax12dpositionsdouble 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 VectorMax12drest_positions,
    
const VectorMax12dpositions) 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 VectorMax12drest_positions,
    
const VectorMax12dpositions) 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:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" [label="ipc::VertexVertexCandidate" tooltip="ipc::VertexVertexCandidate"] "1" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision" fillcolor="#BFBFBF"] "4" -> "3" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::VertexVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" [label="ipc::VertexVertexCandidate" tooltip="ipc::VertexVertexCandidate"] "1" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision" fillcolor="#BFBFBF"] "4" -> "3" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] }

Public Functions

inline VertexVertexNormalCollision(
    
const VertexVertexCandidatecandidate);
inline VertexVertexNormalCollision(const long _vertex0_id,
    
const long _vertex1_idconst double _weight,
    
const Eigen::SparseVector<double>_weight_gradient);
VertexVertexCandidate(long vertex0_idlong vertex1_id);

Friends

template <typename H>
inline friend H AbslHashValue(
    
H hconst VertexVertexNormalCollisionvv);

Edge-Vertex Normal Collision

class EdgeVertexNormalCollision : public ipc::EdgeVertexCandidate,
                                 
 public ipc::NormalCollision;

Inheritence diagram for ipc::EdgeVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate"] "1" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::EdgeVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate"] "1" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Public Functions

inline EdgeVertexNormalCollision(
    
const EdgeVertexCandidatecandidate);
inline EdgeVertexNormalCollision(const long _edge_id,
    
const long _vertex_idconst double _weight,
    
const Eigen::SparseVector<double>_weight_gradient);
inline virtual PointEdgeDistanceType known_dtype() const override;
EdgeVertexCandidate(long edge_idlong vertex_id);

Friends

template <typename H>
inline friend H AbslHashValue(
    
H hconst EdgeVertexNormalCollisionev);

Edge-Edge Normal Collision

class EdgeEdgeNormalCollision : public ipc::EdgeEdgeCandidate,
                               
 public ipc::NormalCollision;

Inheritence diagram for ipc::EdgeEdgeNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate"] "1" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::EdgeEdgeNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate"] "1" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Public Functions

EdgeEdgeNormalCollision(const long edge0_idconst long edge1_id,
    
const double eps_x,
    
const EdgeEdgeDistanceType dtype = EdgeEdgeDistanceType::AUTO);
EdgeEdgeNormalCollision(const EdgeEdgeCandidatecandidate,
    
const double eps_x,
    
const EdgeEdgeDistanceType dtype = EdgeEdgeDistanceType::AUTO);
EdgeEdgeNormalCollision(const long edge0_idconst long edge1_id,
    
const double eps_xconst 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 VectorMax12drest_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 VectorMax12dpositions) 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 VectorMax12dpositionsdouble 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 VectorMax12dpositions) 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 VectorMax12dpositionsdouble 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 VectorMax12dpositions) 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 VectorMax12dpositionsdouble 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 VectorMax12drest_positions,
    
const VectorMax12dpositions) 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 VectorMax12drest_positions,
    
const VectorMax12dpositions) 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 EdgeEdgeNormalCollisionother) const;
bool operator!=(const EdgeEdgeNormalCollisionother) const;
bool operator<(const EdgeEdgeNormalCollisionother) const;

Public Members

double eps_x;

Mollifier activation threshold.

EdgeEdgeDistanceType dtype;

Cached distance type.

Some EE collisions are mollified EV or VV collisions.

Friends

template <typename H>
inline friend H AbslHashValue(
    
H hconst EdgeEdgeNormalCollisionee);

Face-Vertex Normal Collision

class FaceVertexNormalCollision : public ipc::FaceVertexCandidate,
                                 
 public ipc::NormalCollision;

Inheritence diagram for ipc::FaceVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate"] "1" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::FaceVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate"] "1" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision" fillcolor="#BFBFBF"] "4" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "4" [dir=forward tooltip="public-inheritance"] "4" -> "3" [dir=forward tooltip="public-inheritance"] }

Public Functions

inline FaceVertexNormalCollision(
    
const FaceVertexCandidatecandidate);
inline FaceVertexNormalCollision(const long _face_id,
    
const long _vertex_idconst double _weight,
    
const Eigen::SparseVector<double>_weight_gradient);
inline virtual PointTriangleDistanceType
known_dtype() const override;
FaceVertexCandidate(long face_idlong vertex_id);

Friends

template <typename H>
inline friend H AbslHashValue(
    
H hconst FaceVertexNormalCollisionfv);

Plane-Vertex Normal Collision

class PlaneVertexNormalCollision : public ipc::NormalCollision;

Inheritence diagram for ipc::PlaneVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "1" [label="ipc::PlaneVertexNormalCollision" tooltip="ipc::PlaneVertexNormalCollision" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::PlaneVertexNormalCollision:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "1" [label="ipc::PlaneVertexNormalCollision" tooltip="ipc::PlaneVertexNormalCollision" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Public Functions

PlaneVertexNormalCollision(const VectorMax3dplane_origin,
    
const VectorMax3dplane_normalconst 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::MatrixXiedges,
    
const Eigen::MatrixXifaces) 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 VectorMax12dpoint) 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 VectorMax12dpoint) 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 VectorMax12dpoint) 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.

virtual VectorMax4d compute_coefficients(
    
const VectorMax12dpositions) const override;

Compute the coefficients of the stencil.

Parameters:
const VectorMax12d &positions

Vertex positions.

Returns:

Coefficients of the stencil.

virtual bool ccd(const VectorMax12dvertices_t0,
    
const VectorMax12dvertices_t1doubletoi,
    
const double min_distance = 0.0const double tmax = 1.0,
    
const NarrowPhaseCCDnarrow_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::MatrixXdverticesconst Eigen::MatrixXiedges,
    
const Eigen::MatrixXifaces) const;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
const Eigen::MatrixXd &vertices

Collision mesh vertices

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

Coefficients of the stencil.

virtual VectorMax4d compute_coefficients(
    
const VectorMax12dpositions) const
   
 = 0;

Compute the coefficients of the stencil s.t.

d(x) = ‖∑ cᵢ xᵢ‖².

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Coefficients of the stencil.

inline double compute_distance(const Eigen::MatrixXdvertices,
    
const Eigen::MatrixXiedges,
    
const Eigen::MatrixXifaces) const;

Compute the distance of the stencil.

Parameters:
const Eigen::MatrixXd &vertices

Collision mesh vertices

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

Distance of the stencil.

virtual double compute_distance(const VectorMax12dpositions) const
   
 = 0;

Compute the distance of the stencil.

Note

positions can be computed as stencil.dof(vertices, edges, faces)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

inline VectorMax12d compute_distance_gradient(
    
const Eigen::MatrixXdverticesconst Eigen::MatrixXiedges,
    
const Eigen::MatrixXifaces) const;

Compute the distance gradient of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
const Eigen::MatrixXd &vertices

Collision mesh vertices

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

virtual VectorMax12d compute_distance_gradient(
    
const VectorMax12dpositions) 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)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Distance gradient of the stencil w.r.t. the stencil’s vertex positions.

inline MatrixMax12d compute_distance_hessian(
    
const Eigen::MatrixXdverticesconst Eigen::MatrixXiedges,
    
const Eigen::MatrixXifaces) const;

Compute the distance Hessian of the stencil w.r.t.

the stencil’s vertex positions.

Parameters:
const Eigen::MatrixXd &vertices

Collision mesh vertices

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex positions.

virtual MatrixMax12d compute_distance_hessian(
    
const VectorMax12dpositions) 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)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Distance Hessian of the stencil w.r.t. the stencil’s vertex 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: Feb 18, 2025