Candidates

Candidates

class Candidates

Public Functions

Candidates() = default
void build(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices
,
    const
 double inflation_radius = 0
,
    const
 BroadPhaseMethod broad_phase_method
    
= DEFAULT_BROAD_PHASE_METHOD
)

Initialize the set of discrete collision detection candidates.

Parameters:
const CollisionMesh &mesh

The surface of the collision mesh.

const Eigen::MatrixXd &vertices

Surface vertex positions (rowwise).

const double inflation_radius = 0

Amount to inflate the bounding boxes.

const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD

Broad phase method to use.

void build(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices_t0
,
    const
 Eigen::MatrixXd& vertices_t1
,
    const
 double inflation_radius = 0
,
    const
 BroadPhaseMethod broad_phase_method
    
= DEFAULT_BROAD_PHASE_METHOD
)

Initialize the set of continuous collision detection candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The surface of the collision mesh.

const Eigen::MatrixXd &vertices_t0

Surface vertex starting positions (rowwise).

const Eigen::MatrixXd &vertices_t1

Surface vertex ending positions (rowwise).

const double inflation_radius = 0

Amount to inflate the bounding boxes.

const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD

Broad phase method to use.

size_t size() const
bool empty() const
void clear()
ContinuousCollisionCandidate& operator[](size_t i)
const ContinuousCollisionCandidate& operator[](size_t i) const
bool is_step_collision_free(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices_t0
,
    const
 Eigen::MatrixXd& vertices_t1
,
    const
 double min_distance = 0.0
,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
)
 const

Determine if the step is collision free from the set of candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &vertices_t0

Surface vertex starting positions (rowwise).

const Eigen::MatrixXd &vertices_t1

Surface vertex ending positions (rowwise).

const double min_distance = 0.0

The minimum distance allowable between any two elements.

const double tolerance = DEFAULT_CCD_TOLERANCE

The tolerance for the CCD algorithm.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

The maximum number of iterations for the CCD algorithm.

Returns:

True if any collisions occur.

double compute_collision_free_stepsize(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices_t0
,
    const
 Eigen::MatrixXd& vertices_t1
,
    const
 double min_distance = 0.0
,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
)
 const

Computes a maximal step size that is collision free using the set of collision candidates.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &vertices_t0

Surface vertex starting positions (rowwise). Assumed to be intersection free.

const Eigen::MatrixXd &vertices_t1

Surface vertex ending positions (rowwise).

const double min_distance = 0.0

The minimum distance allowable between any two elements.

const double tolerance = DEFAULT_CCD_TOLERANCE

The tolerance for the CCD algorithm.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

The maximum number of iterations for the CCD algorithm.

Returns:

A step-size \(\in [0, 1]\) that is collision free. A value of 1.0 if a full step and 0.0 is no step.

double compute_noncandidate_conservative_stepsize(
    const
 CollisionMesh& mesh
, const Eigen::MatrixXd& displacements,
    const
 double dhat
)
 const

Computes a conservative bound on the largest-feasible step size for surface primitives not in collision.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &displacements

Surface vertex displacements (rowwise).

const double dhat

Barrier activation distance.

double compute_cfl_stepsize(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices_t0
,
    const
 Eigen::MatrixXd& vertices_t1
, const double dhat,
    const
 BroadPhaseMethod broad_phase_method
    
= DEFAULT_BROAD_PHASE_METHOD
,
    const
 double min_distance = 0.0
,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
)
 const

Computes a CFL-inspired CCD maximum step step size.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &vertices_t0

Surface vertex starting positions (rowwise).

const Eigen::MatrixXd &vertices_t1

Surface vertex ending positions (rowwise).

const double dhat

Barrier activation distance.

const double min_distance = 0.0

The minimum distance allowable between any two elements.

const double tolerance = DEFAULT_CCD_TOLERANCE

The tolerance for the CCD algorithm.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

The maximum number of iterations for the CCD algorithm.

bool save_obj(const std::string& filename,
    const
 Eigen::MatrixXd& vertices
, const Eigen::MatrixXi& edges,
    const
 Eigen::MatrixXi& faces
)
 const

Public Members

std::vector<VertexVertexCandidate> vv_candidates
std::vector<EdgeVertexCandidate> ev_candidates
std::vector<EdgeEdgeCandidate> ee_candidates
std::vector<FaceVertexCandidate> fv_candidates

Collision Stencil

class CollisionStencil

A stencil representing a collision between at most four vertices.

Subclassed by ipc::Collision, ipc::ContinuousCollisionCandidate, ipc::FrictionCollision

Public Functions

virtual ~CollisionStencil() = default
virtual int num_vertices() const = 0

Get the number of vertices in the collision stencil.

inline int dim(const int ndof) const

Get the dimension of the collision stencil.

Parameters:
const int ndof

Number of degrees of freedom in the stencil.

Returns:

The dimension of the collision stencil.

virtual std::array<long, 4> vertex_ids(const Eigen::MatrixXi& edges,
    const
 Eigen::MatrixXi& faces
)
 const = 0

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.

template <typename T>
inline
 std::array<VectorMax3<T>, 4> vertices(
    const
 MatrixX<T>& vertices
, const Eigen::MatrixXi& edges,
    const
 Eigen::MatrixXi& faces
)
 const

Get the vertex attributes of the collision stencil.

Template Parameters:
typename T

Type of the attributes

Parameters:
const MatrixX<T> &vertices

Vertex attributes

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

The vertex positions of the collision stencil. Size is always 4, but elements i > num_vertices() are NaN.

template <typename T>
inline
 VectorMax12<T> dof(const MatrixX<T>& X,
    const
 Eigen::MatrixXi& edges
,
    const
 Eigen::MatrixXi& faces
)
 const

Select this stencil’s DOF from the full matrix of DOF.

Template Parameters:
typename T

Type of the DOF

Parameters:
const MatrixX<T> &X

Full matrix of DOF (rowwise).

const Eigen::MatrixXi &edges

Collision mesh edges

const Eigen::MatrixXi &faces

Collision mesh faces

Returns:

This stencil’s DOF.

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)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Distance of the stencil.

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)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

Distance gradient 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)

Parameters:
const VectorMax12d &positions

Stencil’s vertex positions.

Returns:

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

Continuous Collision Candidate

class ContinuousCollisionCandidate
    
: public virtual ipc::CollisionStencil

Virtual class for candidates that support CCD.

Subclassed by ipc::EdgeEdgeCandidate, ipc::EdgeVertexCandidate, ipc::FaceVertexCandidate, ipc::VertexVertexCandidate

Public Functions

virtual ~ContinuousCollisionCandidate() = default
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
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)
 const = 0

Perform narrow-phase CCD on the candidate.

Parameters:
const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

double &toi

Computed time of impact (normalized).

const double min_distance = 0.0

Minimum separation distance between primitives.

const double tmax = 1.0

Maximum time (normalized) to look for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] CCD tolerance used by Tight-Inclusion CCD.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum iterations used by Tight-Inclusion CCD.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling value used to avoid taking steps exactly to impact.

Returns:

If the candidate had a collision over the time interval.

std::ostream& write_ccd_query(std::ostream& out,
    const
 VectorMax12d& vertices_t0
,
    const
 VectorMax12d& vertices_t1
)
 const

Write the CCD query to a stream.

Parameters:
std::ostream &out

Stream to write to.

const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

Returns:

The stream.

Vertex-Vertex Candidate

class VertexVertexCandidate
    
: public ipc::ContinuousCollisionCandidate

Subclassed by ipc::VertexVertexCollision, ipc::VertexVertexFrictionCollision

Public Functions

VertexVertexCandidate(long vertex0_id, long vertex1_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 indices of the vertices.

Parameters:
const Eigen::MatrixXi &edges

edge matrix of mesh

const Eigen::MatrixXi &faces

face matrix of mesh

Returns:

List of vertex indices

virtual double compute_distance(
    const
 VectorMax12d& positions
)
 const override

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.

virtual VectorMax12d compute_distance_gradient(
    const
 VectorMax12d& positions
)
 const override

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.

virtual MatrixMax12d compute_distance_hessian(
    const
 VectorMax12d& positions
)
 const override

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.

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
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)
 const override

Perform narrow-phase CCD on the candidate.

Parameters:
const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

double &toi

Computed time of impact (normalized).

const double min_distance = 0.0

Minimum separation distance between primitives.

const double tmax = 1.0

Maximum time (normalized) to look for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] CCD tolerance used by Tight-Inclusion CCD.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum iterations used by Tight-Inclusion CCD.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling value used to avoid taking steps exactly to impact.

Returns:

If the candidate had a collision over the time interval.

bool operator==(const VertexVertexCandidate& other) const
bool operator!=(const VertexVertexCandidate& other) const
bool operator<(const VertexVertexCandidate& other) const

Compare EdgeVertexCandidates for sorting.

Public Members

long vertex0_id

ID of the first vertex.

long vertex1_id

ID of the second vertex.

Friends

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

Edge-Vertex Candidate

class EdgeVertexCandidate
    
: public ipc::ContinuousCollisionCandidate

Subclassed by ipc::EdgeVertexCollision, ipc::EdgeVertexFrictionCollision

Public Functions

EdgeVertexCandidate(long edge_id, 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& positions
)
 const override

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.

virtual VectorMax12d compute_distance_gradient(
    const
 VectorMax12d& positions
)
 const override

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.

virtual MatrixMax12d compute_distance_hessian(
    const
 VectorMax12d& positions
)
 const override

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.

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
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)
 const override

Perform narrow-phase CCD on the candidate.

Parameters:
const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

double &toi

Computed time of impact (normalized).

const double min_distance = 0.0

Minimum separation distance between primitives.

const double tmax = 1.0

Maximum time (normalized) to look for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] CCD tolerance used by Tight-Inclusion CCD.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum iterations used by Tight-Inclusion CCD.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling value used to avoid taking steps exactly to impact.

Returns:

If the candidate had a collision over the time interval.

inline virtual PointEdgeDistanceType known_dtype() const
bool operator==(const EdgeVertexCandidate& other) const
bool operator!=(const EdgeVertexCandidate& other) const
bool operator<(const EdgeVertexCandidate& other) const

Compare EdgeVertexCandidates for sorting.

Public Members

long edge_id

ID of the edge.

long vertex_id

ID of the vertex.

Friends

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

Edge-Edge Candidate

class EdgeEdgeCandidate : public ipc::ContinuousCollisionCandidate

Subclassed by ipc::EdgeEdgeCollision, ipc::EdgeEdgeFrictionCollision

Public Functions

EdgeEdgeCandidate(long edge0_id, long edge1_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& positions
)
 const override

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.

virtual VectorMax12d compute_distance_gradient(
    const
 VectorMax12d& positions
)
 const override

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.

virtual MatrixMax12d compute_distance_hessian(
    const
 VectorMax12d& positions
)
 const override

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.

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
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)
 const override

Perform narrow-phase CCD on the candidate.

Parameters:
const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

double &toi

Computed time of impact (normalized).

const double min_distance = 0.0

Minimum separation distance between primitives.

const double tmax = 1.0

Maximum time (normalized) to look for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] CCD tolerance used by Tight-Inclusion CCD.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum iterations used by Tight-Inclusion CCD.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling value used to avoid taking steps exactly to impact.

Returns:

If the candidate had a collision over the time interval.

inline virtual EdgeEdgeDistanceType known_dtype() const
bool operator==(const EdgeEdgeCandidate& other) const
bool operator!=(const EdgeEdgeCandidate& other) const
bool operator<(const EdgeEdgeCandidate& other) const

Compare EdgeEdgeCandidates for sorting.

Public Members

long edge0_id

ID of the first edge.

long edge1_id

ID of the second edge.

Friends

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

Edge-Face Candidate

class EdgeFaceCandidate

Candidate for intersection between edge and face.

Not included in Candidates because it is not a collision candidate.

Public Functions

EdgeFaceCandidate(long edge_id, long face_id)
bool operator==(const EdgeFaceCandidate& other) const
bool operator!=(const EdgeFaceCandidate& other) const
bool operator<(const EdgeFaceCandidate& other) const

Compare EdgeFaceCandidate for sorting.

Public Members

long edge_id

ID of the edge.

long face_id

ID of the face.

Friends

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

Face-Vertex Candidate

class FaceVertexCandidate
    
: public ipc::ContinuousCollisionCandidate

Subclassed by ipc::FaceVertexCollision, ipc::FaceVertexFrictionCollision

Public Functions

FaceVertexCandidate(long face_id, 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& positions
)
 const override

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.

virtual VectorMax12d compute_distance_gradient(
    const
 VectorMax12d& positions
)
 const override

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.

virtual MatrixMax12d compute_distance_hessian(
    const
 VectorMax12d& positions
)
 const override

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.

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
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)
 const override

Perform narrow-phase CCD on the candidate.

Parameters:
const VectorMax12d &vertices_t0

Stencil vertices at the start of the time step.

const VectorMax12d &vertices_t1

Stencil vertices at the end of the time step.

double &toi

Computed time of impact (normalized).

const double min_distance = 0.0

Minimum separation distance between primitives.

const double tmax = 1.0

Maximum time (normalized) to look for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] CCD tolerance used by Tight-Inclusion CCD.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum iterations used by Tight-Inclusion CCD.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling value used to avoid taking steps exactly to impact.

Returns:

If the candidate had a collision over the time interval.

inline virtual PointTriangleDistanceType known_dtype() const
bool operator==(const FaceVertexCandidate& other) const
bool operator!=(const FaceVertexCandidate& other) const
bool operator<(const FaceVertexCandidate& other) const

Compare FaceVertexCandidate for sorting.

Public Members

long face_id

ID of the face.

long vertex_id

ID of the vertex.

Friends

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

Last update: Apr 03, 2024