Continuous Collision Detection

bool ipc::is_step_collision_free(const CollisionMesh& mesh,
    const
 Eigen::MatrixXd& vertices_t0
,
    const
 Eigen::MatrixXd& vertices_t1
,
    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
)

Determine if the step is collision free.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &vertices_t0

Surface vertex vertices at start as rows of a matrix.

const Eigen::MatrixXd &vertices_t1

Surface vertex vertices at end as rows of a matrix.

const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD

The broad phase method to use.

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 ipc::compute_collision_free_stepsize(
    const
 CollisionMesh& mesh
, const Eigen::MatrixXd& vertices_t0,
    const
 Eigen::MatrixXd& vertices_t1
,
    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
)

Computes a maximal step size that is collision free.

Note

Assumes the trajectory is linear.

Parameters:
const CollisionMesh &mesh

The collision mesh.

const Eigen::MatrixXd &vertices_t0

Vertex vertices at start as rows of a matrix. Assumes vertices_t0 is intersection free.

const Eigen::MatrixXd &vertices_t1

Surface vertex vertices at end as rows of a matrix.

const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD

The broad phase method to use.

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.

static constexpr double ipc::DEFAULT_CCD_TOLERANCE = 1e-6

The default tolerance used with Tight-Inclusion CCD.

static constexpr long ipc::DEFAULT_CCD_MAX_ITERATIONS = 10'000'000l

The default maximum number of iterations used with Tight-Inclusion CCD.

static constexpr double ipc::DEFAULT_CCD_CONSERVATIVE_RESCALING
    
= 0.8

The default conservative rescaling value used to avoid taking steps exactly to impact.

Individual CCD Functions

bool ipc::point_point_ccd(const VectorMax3d& p0_t0,
    const
 VectorMax3d& p1_t0
, const VectorMax3d& p0_t1,
    const
 VectorMax3d& p1_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
)

Computes the time of impact between two points using continuous collision detection.

Parameters:
const VectorMax3d &p0_t0

[in] The initial position of the first point.

const VectorMax3d &p1_t0

[in] The initial position of the second point.

const VectorMax3d &p0_t1

[in] The final position of the first point.

const VectorMax3d &p1_t1

[in] The final position of the second point.

double &toi

[out] The time of impact between the two points.

const double min_distance = 0.0

[in] The minimum distance between the points.

const double tmax = 1.0

[in] The maximum time to check for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] The error tolerance for the time of impact.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] The maximum number of iterations to perform.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::point_edge_ccd(const VectorMax3d& p_t0,
    const
 VectorMax3d& e0_t0
, const VectorMax3d& e1_t0,
    const
 VectorMax3d& p_t1
, const VectorMax3d& e0_t1,
    const
 VectorMax3d& e1_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
)

Computes the time of impact between a point and an edge using continuous collision detection.

Parameters:
const VectorMax3d &p_t0

[in] The initial position of the point.

const VectorMax3d &e0_t0

[in] The initial position of the first endpoint of the edge.

const VectorMax3d &e1_t0

[in] The initial position of the second endpoint of the edge.

const VectorMax3d &p_t1

[in] The final position of the point.

const VectorMax3d &e0_t1

[in] The final position of the first endpoint of the edge.

const VectorMax3d &e1_t1

[in] The final position of the second endpoint of the edge.

double &toi

[out] The time of impact between the point and the edge.

const double min_distance = 0.0

[in] The minimum distance between the objects.

const double tmax = 1.0

[in] The maximum time to check for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] The error tolerance for the time of impact.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] The maximum number of iterations to perform.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::edge_edge_ccd(const Eigen::Vector3d& ea0_t0,
    const
 Eigen::Vector3d& ea1_t0
, const Eigen::Vector3d& eb0_t0,
    const
 Eigen::Vector3d& eb1_t0
, const Eigen::Vector3d& ea0_t1,
    const
 Eigen::Vector3d& ea1_t1
, const Eigen::Vector3d& eb0_t1,
    const
 Eigen::Vector3d& eb1_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
)

Computes the time of impact between two edges in 3D using continuous collision detection.

Parameters:
const Eigen::Vector3d &ea0_t0

[in] The initial position of the first endpoint of the first edge.

const Eigen::Vector3d &ea1_t0

[in] The initial position of the second endpoint of the first edge.

const Eigen::Vector3d &eb0_t0

[in] The initial position of the first endpoint of the second edge.

const Eigen::Vector3d &eb1_t0

[in] The initial position of the second endpoint of the second edge.

const Eigen::Vector3d &ea0_t1

[in] The final position of the first endpoint of the first edge.

const Eigen::Vector3d &ea1_t1

[in] The final position of the second endpoint of the first edge.

const Eigen::Vector3d &eb0_t1

[in] The final position of the first endpoint of the second edge.

const Eigen::Vector3d &eb1_t1

[in] The final position of the second endpoint of the second edge.

double &toi

[out] The time of impact between the two edges.

const double min_distance = 0.0

[in] The minimum distance between the objects.

const double tmax = 1.0

[in] The maximum time to check for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] The error tolerance for the time of impact.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] The maximum number of iterations to perform.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::point_triangle_ccd(const Eigen::Vector3d& p_t0,
    const
 Eigen::Vector3d& t0_t0
, const Eigen::Vector3d& t1_t0,
    const
 Eigen::Vector3d& t2_t0
, const Eigen::Vector3d& p_t1,
    const
 Eigen::Vector3d& t0_t1
, const Eigen::Vector3d& t1_t1,
    const
 Eigen::Vector3d& t2_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
)

Computes the time of impact between a point and a triangle in 3D using continuous collision detection.

Parameters:
const Eigen::Vector3d &p_t0

[in] The initial position of the point.

const Eigen::Vector3d &t0_t0

[in] The initial position of the first vertex of the triangle.

const Eigen::Vector3d &t1_t0

[in] The initial position of the second vertex of the triangle.

const Eigen::Vector3d &t2_t0

[in] The initial position of the third vertex of the triangle.

const Eigen::Vector3d &p_t1

[in] The final position of the point.

const Eigen::Vector3d &t0_t1

[in] The final position of the first vertex of the triangle.

const Eigen::Vector3d &t1_t1

[in] The final position of the second vertex of the triangle.

const Eigen::Vector3d &t2_t1

[in] The final position of the third vertex of the triangle.

double &toi

[out] The time of impact between the point and the triangle.

const double min_distance = 0.0

[in] The minimum distance between the objects.

const double tmax = 1.0

[in] The maximum time to check for collisions.

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] The error tolerance for the time of impact.

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] The maximum number of iterations to perform.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

Generic Interface

bool ipc::ccd_strategy(
    const
 std::function<bool(long, double, bool, double&)>& ccd
,
    const
 long max_iterations
, const double min_distance,
    const
 double initial_distance
,
    const
 double conservative_rescaling
, double& toi)

Perform the CCD strategy outlined by Li et al.

[2020].

Parameters:
const std::function<bool(long, double, bool, double&)> &ccd

[in] The continuous collision detection function.

const long max_iterations

[in] The maximum number of iterations to perform.

const double min_distance

[in] The minimum distance between the objects.

const double initial_distance

[in] The initial distance between the objects.

const double conservative_rescaling

[in] The conservative rescaling of the time of impact.

double &toi

[out] Output time of impact.

Returns:

True if a collision was detected, false otherwise.

Additive CCD

bool ipc::additive_ccd::point_point_ccd(const VectorMax3d& p0_t0,
    const
 VectorMax3d& p1_t0
, const VectorMax3d& p0_t1,
    const
 VectorMax3d& p1_t1
, double& toi,
    const
 double min_distance = 0.0
, const double tmax = 1.0,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Computes the time of impact between two points using continuous collision detection.

Parameters:
const VectorMax3d &p0_t0

The initial position of the first point.

const VectorMax3d &p1_t0

The initial position of the second point.

const VectorMax3d &p0_t1

The final position of the first point.

const VectorMax3d &p1_t1

The final position of the second point.

double &toi

[out] The time of impact between the two points.

const double min_distance = 0.0

The minimum distance between two objects.

const double tmax = 1.0

The maximum time to check for collisions.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::additive_ccd::point_edge_ccd(const VectorMax3d& p_t0,
    const
 VectorMax3d& e0_t0
, const VectorMax3d& e1_t0,
    const
 VectorMax3d& p_t1
, const VectorMax3d& e0_t1,
    const
 VectorMax3d& e1_t1
, double& toi,
    const
 double min_distance = 0.0
, const double tmax = 1.0,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Computes the time of impact between a point and an edge using continuous collision detection.

Parameters:
const VectorMax3d &p_t0

The initial position of the point.

const VectorMax3d &e0_t0

The initial position of the first endpoint of the edge.

const VectorMax3d &e1_t0

The initial position of the second endpoint of the edge.

const VectorMax3d &p_t1

The final position of the point.

const VectorMax3d &e0_t1

The final position of the first endpoint of the edge.

const VectorMax3d &e1_t1

The final position of the second endpoint of the edge.

double &toi

[out] The time of impact between the point and the edge.

const double min_distance = 0.0

The minimum distance between two objects.

const double tmax = 1.0

The maximum time to check for collisions.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::additive_ccd::edge_edge_ccd(const Eigen::Vector3d& ea0_t0,
    const
 Eigen::Vector3d& ea1_t0
, const Eigen::Vector3d& eb0_t0,
    const
 Eigen::Vector3d& eb1_t0
, const Eigen::Vector3d& ea0_t1,
    const
 Eigen::Vector3d& ea1_t1
, const Eigen::Vector3d& eb0_t1,
    const
 Eigen::Vector3d& eb1_t1
, double& toi,
    const
 double min_distance = 0.0
, const double tmax = 1.0,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Computes the time of impact between two edges using continuous collision detection.

Parameters:
const Eigen::Vector3d &ea0_t0

The initial position of the first endpoint of the first edge.

const Eigen::Vector3d &ea1_t0

The initial position of the second endpoint of the first edge.

const Eigen::Vector3d &eb0_t0

The initial position of the first endpoint of the second edge.

const Eigen::Vector3d &eb1_t0

The initial position of the second endpoint of the second edge.

const Eigen::Vector3d &ea0_t1

The final position of the first endpoint of the first edge.

const Eigen::Vector3d &ea1_t1

The final position of the second endpoint of the first edge.

const Eigen::Vector3d &eb0_t1

The final position of the first endpoint of the second edge.

const Eigen::Vector3d &eb1_t1

The final position of the second endpoint of the second edge.

double &toi

[out] The time of impact between the two edges.

const double min_distance = 0.0

The minimum distance between two objects.

const double tmax = 1.0

The maximum time to check for collisions.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

bool ipc::additive_ccd::point_triangle_ccd(
    const
 Eigen::Vector3d& p_t0
, const Eigen::Vector3d& t0_t0,
    const
 Eigen::Vector3d& t1_t0
, const Eigen::Vector3d& t2_t0,
    const
 Eigen::Vector3d& p_t1
, const Eigen::Vector3d& t0_t1,
    const
 Eigen::Vector3d& t1_t1
, const Eigen::Vector3d& t2_t1,
    double
& toi
, const double min_distance = 0.0,
    const
 double tmax = 1.0
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Computes the time of impact between a point and a triangle using continuous collision detection.

Parameters:
const Eigen::Vector3d &p_t0

The initial position of the point.

const Eigen::Vector3d &t0_t0

The initial position of the first vertex of the triangle.

const Eigen::Vector3d &t1_t0

The initial position of the second vertex of the triangle.

const Eigen::Vector3d &t2_t0

The initial position of the third vertex of the triangle.

const Eigen::Vector3d &p_t1

The final position of the point.

const Eigen::Vector3d &t0_t1

The final position of the first vertex of the triangle.

const Eigen::Vector3d &t1_t1

The final position of the second vertex of the triangle.

const Eigen::Vector3d &t2_t1

The final position of the third vertex of the triangle.

double &toi

[out] The time of impact between the point and the triangle.

const double min_distance = 0.0

The minimum distance between two objects.

const double tmax = 1.0

The maximum time to check for collisions.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

The conservative rescaling of the time of impact.

Returns:

True if a collision was detected, false otherwise.

Generic Interface

bool ipc::additive_ccd::additive_ccd(VectorMax12d x,
    const
 VectorMax12d& dx
,
    const
 std::function<double(const VectorMax12d&)>&
        distance_squared
,
    const
 double max_disp_mag
, double& toi,
    const
 double min_distance = 0.0
, const double tmax = 1.0,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Computes the time of impact between two objects using additive continuous collision detection.

Parameters:
const std::function<double(const VectorMax12d&)> &distance_squared

A function that computes the squared distance between the two objects at a given time.

double &toi

[out] The time of impact between the two objects.

const double min_distance = 0.0

The minimum distance between the objects.

const double tmax = 1.0

The maximum time to check for collisions.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

The amount to rescale the objects by to ensure conservative advancement.

Returns:

True if a collision was detected, false otherwise.

Nonlinear CCD

class NonlinearTrajectory

A nonlinear trajectory is a function that maps time to a point in space.

Subclassed by ipc::IntervalNonlinearTrajectory

Public Functions

virtual ~NonlinearTrajectory() = default
virtual VectorMax3d operator()(const double t) const = 0

Compute the point’s position at time t.

virtual double max_distance_from_linear(
    const
 double t0
, const double t1)
 const = 0

Compute the maximum distance from the nonlinear trajectory to a linearized trajectory.

Parameters:
const double t0

[in] Start time of the trajectory

const double t1

[in] End time of the trajectory

class IntervalNonlinearTrajectory
    
: public virtual ipc::NonlinearTrajectory

A nonlinear trajectory with an implementation of the max_distance_from_linear function using interval arithmetic.

Public Functions

virtual ~IntervalNonlinearTrajectory() = default
virtual VectorMax3I operator()(const filib::Interval& t) const = 0

Compute the point’s position over a time interval t.

virtual double max_distance_from_linear(
    const
 double t0
, const double t1)
 const

Compute the maximum distance from the nonlinear trajectory to a linearized trajectory.

Note

This uses interval arithmetic to compute the maximum distance. If you know a tighter bound on the maximum distance, it is recommended to override this function.

Parameters:
const double t0

[in] Start time of the trajectory

const double t1

[in] End time of the trajectory

bool ipc::point_point_nonlinear_ccd(const NonlinearTrajectory& p0,
    const
 NonlinearTrajectory& p1
, double& toi,
    const
 double tmax = 1.0
, const double min_distance = 0,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Perform nonlinear CCD between two points moving along nonlinear trajectories.

Parameters:
const NonlinearTrajectory &p0

[in] First point’s trajectory

const NonlinearTrajectory &p1

[in] Second point’s trajectory

double &toi

[out] Output time of impact

const double tmax = 1.0

[in] Maximum time to check for collision

const double min_distance = 0

[in] Minimum separation distance between the two points

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] Tolerance for the linear CCD algorithm

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum number of iterations for the linear CCD algorithm

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling of the time of impact

Returns:

True if the two points collide, false otherwise.

bool ipc::point_edge_nonlinear_ccd(const NonlinearTrajectory& p,
    const
 NonlinearTrajectory& e0
, const NonlinearTrajectory& e1,
    double
& toi
, const double tmax = 1.0,
    const
 double min_distance = 0
,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Perform nonlinear CCD between a point and a linear edge moving along nonlinear trajectories.

Parameters:
const NonlinearTrajectory &p

[in] Point’s trajectory

const NonlinearTrajectory &e0

[in] Edge’s first endpoint’s trajectory

const NonlinearTrajectory &e1

[in] Edge’s second endpoint’s trajectory

double &toi

[out] Output time of impact

const double tmax = 1.0

[in] Maximum time to check for collision

const double min_distance = 0

[in] Minimum separation distance between the point and the edge

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] Tolerance for the linear CCD algorithm

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum number of iterations for the linear CCD algorithm

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling of the time of impact

Returns:

True if the point and edge collide, false otherwise.

bool ipc::edge_edge_nonlinear_ccd(const NonlinearTrajectory& ea0,
    const
 NonlinearTrajectory& ea1
, const NonlinearTrajectory& eb0,
    const
 NonlinearTrajectory& eb1
, double& toi,
    const
 double tmax = 1.0
, const double min_distance = 0,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Perform nonlinear CCD between two linear edges moving along nonlinear trajectories.

Parameters:
const NonlinearTrajectory &ea0

[in] First edge’s first endpoint’s trajectory

const NonlinearTrajectory &ea1

[in] First edge’s second endpoint’s trajectory

const NonlinearTrajectory &eb0

[in] Second edge’s first endpoint’s trajectory

const NonlinearTrajectory &eb1

[in] Second edge’s second endpoint’s trajectory

double &toi

[out] Output time of impact

const double tmax = 1.0

[in] Maximum time to check for collision

const double min_distance = 0

[in] Minimum separation distance between the two edges

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] Tolerance for the linear CCD algorithm

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum number of iterations for the linear CCD algorithm

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling of the time of impact

Returns:

True if the two edges collide, false otherwise.

bool ipc::point_triangle_nonlinear_ccd(const NonlinearTrajectory& p,
    const
 NonlinearTrajectory& t0
, const NonlinearTrajectory& t1,
    const
 NonlinearTrajectory& t2
, double& toi,
    const
 double tmax = 1.0
, const double min_distance = 0,
    const
 double tolerance = DEFAULT_CCD_TOLERANCE
,
    const
 long max_iterations = DEFAULT_CCD_MAX_ITERATIONS
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Perform nonlinear CCD between a point and a linear triangle moving along nonlinear trajectories.

Parameters:
const NonlinearTrajectory &p

[in] Point’s trajectory

const NonlinearTrajectory &t0

[in] Triangle’s first vertex’s trajectory

const NonlinearTrajectory &t1

[in] Triangle’s second vertex’s trajectory

const NonlinearTrajectory &t2

[in] Triangle’s third vertex’s trajectory

double &toi

[out] Output time of impact

const double tmax = 1.0

[in] Maximum time to check for collision

const double min_distance = 0

[in] Minimum separation distance between the two edges

const double tolerance = DEFAULT_CCD_TOLERANCE

[in] Tolerance for the linear CCD algorithm

const long max_iterations = DEFAULT_CCD_MAX_ITERATIONS

[in] Maximum number of iterations for the linear CCD algorithm

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling of the time of impact

Returns:

True if the point and triangle collide, false otherwise.

Generic Interface

bool ipc::conservative_piecewise_linear_ccd(
    const
 std::function<double(const double)>& distance
,
    const
 std::function<double(const double, const double)>&
        max_distance_from_linear
,
    const
 std::function<bool(const double, const double,
        
const double, const bool, double&)>& linear_ccd
,
    double
& toi
, const double tmax = 1.0,
    const
 double min_distance = 0
,
    const
 double conservative_rescaling
    
= DEFAULT_CCD_CONSERVATIVE_RESCALING
)

Perform conservative piecewise linear CCD of a nonlinear trajectories.

Parameters:
const std::function<double(const double)> &distance

[in] Return the distance for a given time in [0, 1].

const std::function<double(const double, const double)> &max_distance_from_linear

[in] Return the maximum distance from the linearized trajectory for a given time interval.

const std::function<bool(const double, const double, const double, const bool, double&)> &linear_ccd

[in] Perform linear CCD on a given time interval.

double &toi

[out] Output time of impact.

const double tmax = 1.0

[in] Maximum time to check for collision.

const double min_distance = 0

[in] Minimum separation distance between the objects.

const double conservative_rescaling = DEFAULT_CCD_CONSERVATIVE_RESCALING

[in] Conservative rescaling of the time of impact.

Returns:


Last update: Apr 03, 2024