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

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"] "7" [label="ipc::VertexVertexNormalCollision" tooltip="ipc::VertexVertexNormalCollision"] "1" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision" fillcolor="#BFBFBF"] "5" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "3" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision"] "4" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision"] "6" [label="ipc::PlaneVertexNormalCollision" tooltip="ipc::PlaneVertexNormalCollision"] "7" -> "1" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "5" -> "1" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] "6" -> "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"] "1" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision" fillcolor="#BFBFBF"] "2" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "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
 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:

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

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)

Friends

template <typename H>
inline
 friend H AbslHashValue(
    H
 h
, const VertexVertexNormalCollision& vv)

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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate"] "1" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision" fillcolor="#BFBFBF"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "2" [label="ipc::EdgeVertexCandidate" tooltip="ipc::EdgeVertexCandidate"] "1" [label="ipc::EdgeVertexNormalCollision" tooltip="ipc::EdgeVertexNormalCollision" fillcolor="#BFBFBF"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [dir=forward tooltip="public-inheritance"] }

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
EdgeVertexCandidate(long edge_id, long vertex_id)

Friends

template <typename H>
inline
 friend H AbslHashValue(
    H
 h
, const EdgeVertexNormalCollision& ev)

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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision" fillcolor="#BFBFBF"] "2" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [dir=forward tooltip="public-inheritance"] "2" -> "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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "1" [label="ipc::EdgeEdgeNormalCollision" tooltip="ipc::EdgeEdgeNormalCollision" fillcolor="#BFBFBF"] "2" [label="ipc::EdgeEdgeCandidate" tooltip="ipc::EdgeEdgeCandidate"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] }

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.

Friends

template <typename H>
inline
 friend H AbslHashValue(
    H
 h
, const EdgeEdgeNormalCollision& ee)

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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "1" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision" fillcolor="#BFBFBF"] "2" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [dir=forward tooltip="public-inheritance"] "2" -> "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"] "5" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::ContinuousCollisionCandidate" tooltip="ipc::ContinuousCollisionCandidate"] "1" [label="ipc::FaceVertexNormalCollision" tooltip="ipc::FaceVertexNormalCollision" fillcolor="#BFBFBF"] "2" [label="ipc::FaceVertexCandidate" tooltip="ipc::FaceVertexCandidate"] "4" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "5" -> "4" [dir=forward tooltip="public-inheritance"] "3" -> "4" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "1" -> "5" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] }

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
FaceVertexCandidate(long face_id, long vertex_id)

Friends

template <typename H>
inline
 friend H AbslHashValue(
    H
 h
, const FaceVertexNormalCollision& fv)

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"] "2" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "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"] "2" [label="ipc::NormalCollision" tooltip="ipc::NormalCollision"] "3" [label="ipc::CollisionStencil" tooltip="ipc::CollisionStencil"] "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 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