Continuous Collision Detection¶
-
bool ipc::is_step_collision_free(const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices_t0,
const Eigen::MatrixXd& vertices_t1,
const double min_distance = 0.0,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD)¶ 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 double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD¶
The broad phase method to use.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
- 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 double min_distance = 0.0,
const BroadPhaseMethod broad_phase_method
= DEFAULT_BROAD_PHASE_METHOD,
const NarrowPhaseCCD& narrow_phase_ccd
= DEFAULT_NARROW_PHASE_CCD)¶ 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 double min_distance = 0.0¶
The minimum distance allowable between any two elements.
- const BroadPhaseMethod broad_phase_method = DEFAULT_BROAD_PHASE_METHOD¶
The broad phase method to use.
- const NarrowPhaseCCD &narrow_phase_ccd = DEFAULT_NARROW_PHASE_CCD¶
The narrow phase CCD algorithm to use.
- 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.
Narrow Phase CCD¶
- class NarrowPhaseCCD¶
Inheritence diagram for ipc::NarrowPhaseCCD:
Subclassed by ipc::AdditiveCCD, ipc::TightInclusionCCD
Public Functions
- NarrowPhaseCCD() = default¶
- virtual ~NarrowPhaseCCD() = default¶
-
virtual bool 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 = 0¶
-
virtual bool 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 = 0¶
-
virtual bool 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 = 0¶
-
virtual bool 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 = 0¶
Tight Inclusion CCD¶
- class TightInclusionCCD : public ipc::NarrowPhaseCCD¶
Inheritence diagram for ipc::TightInclusionCCD:
Collaboration diagram for ipc::TightInclusionCCD:
Public Functions
-
TightInclusionCCD(const double tolerance = DEFAULT_TOLERANCE,
const long max_iterations = DEFAULT_MAX_ITERATIONS,
const double conservative_rescaling
= DEFAULT_CONSERVATIVE_RESCALING)¶ Construct a new AdditiveCCD object.
- Parameters:¶
- const double conservative_rescaling = DEFAULT_CONSERVATIVE_RESCALING¶
The conservative rescaling of the time of impact.
-
virtual bool 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 override¶ 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.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- Returns:¶
True if a collision was detected, false otherwise.
Public Members
- double tolerance¶
Solver tolerance.
- long max_iterations¶
Maximum number of iterations.
- double conservative_rescaling¶
Conservative rescaling of the time of impact.
Public Static Attributes
- static constexpr double DEFAULT_TOLERANCE = 1e-6¶
The default tolerance used with Tight-Inclusion CCD.
- static constexpr long DEFAULT_MAX_ITERATIONS = 10'000'000l¶
The default maximum number of iterations used with Tight-Inclusion CCD.
- static constexpr double DEFAULT_CONSERVATIVE_RESCALING = 0.8¶
The default conservative rescaling value used to avoid taking steps exactly to impact.
- static constexpr double SMALL_TOI = 1e-6¶
Tolerance for small time of impact which triggers rerunning CCD without a minimum separation.
Private Functions
-
bool point_point_ccd_3D(const Eigen::Vector3d& p0_t0,
const Eigen::Vector3d& p1_t0, const Eigen::Vector3d& p0_t1,
const Eigen::Vector3d& p1_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0) const¶ Computes the time of impact between two points in 3D using continuous collision detection.
- Parameters:¶
- const Eigen::Vector3d &p0_t0¶
[in] The initial position of the first point.
- const Eigen::Vector3d &p1_t0¶
[in] The initial position of the second point.
- const Eigen::Vector3d &p0_t1¶
[in] The final position of the first point.
- const Eigen::Vector3d &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 objects.
- const double tmax = 1.0¶
[in] The maximum time to check for collisions.
- Returns:¶
True if a collision was detected, false otherwise.
-
bool point_edge_ccd_3D(const Eigen::Vector3d& p_t0,
const Eigen::Vector3d& e0_t0, const Eigen::Vector3d& e1_t0,
const Eigen::Vector3d& p_t1, const Eigen::Vector3d& e0_t1,
const Eigen::Vector3d& e1_t1, double& toi,
const double min_distance = 0.0, const double tmax = 1.0) const¶ Computes the time of impact between a point and an edge in 3D using continuous collision detection.
- Parameters:¶
- const Eigen::Vector3d &p_t0¶
[in] The initial position of the point.
- const Eigen::Vector3d &e0_t0¶
[in] The initial position of the first endpoint of the edge.
- const Eigen::Vector3d &e1_t0¶
[in] The initial position of the second endpoint of the edge.
- const Eigen::Vector3d &p_t1¶
[in] The final position of the point.
- const Eigen::Vector3d &e0_t1¶
[in] The final position of the first endpoint of the edge.
- const Eigen::Vector3d &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.
- Returns:¶
True if a collision was detected, false otherwise.
Private Static Functions
-
static bool ccd_strategy(
const std::function<bool(double, bool, double&)>& ccd,
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(double, bool, double&)> &ccd¶
[in] The continuous collision detection function.
- 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.
-
TightInclusionCCD(const double tolerance = DEFAULT_TOLERANCE,
Additive CCD¶
- class AdditiveCCD : public ipc::NarrowPhaseCCD¶
Inheritence diagram for ipc::AdditiveCCD:
Collaboration diagram for ipc::AdditiveCCD:
Additive Continuous Collision Detection (CCD) from [Li et al. 2021].
Public Functions
-
AdditiveCCD(const double conservative_rescaling
= DEFAULT_CONSERVATIVE_RESCALING)¶ Construct a new AdditiveCCD object.
- Parameters:¶
- const double conservative_rescaling = DEFAULT_CONSERVATIVE_RESCALING¶
The conservative rescaling of the time of impact.
-
virtual bool 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 override¶ 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.
- conservative_rescaling
The conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- conservative_rescaling
The conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- conservative_rescaling
The conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.
-
virtual bool 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 override¶ 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.
- conservative_rescaling
The conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.
Public Members
- double conservative_rescaling¶
The conservative rescaling value used to avoid taking steps exactly to impact.
Public Static Attributes
- static constexpr double DEFAULT_CONSERVATIVE_RESCALING = 0.9¶
The default conservative rescaling value used to avoid taking steps exactly to impact.
Value choosen to based on [Li et al. 2021].
Private Static Functions
-
static bool 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_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_CONSERVATIVE_RESCALING¶
The amount to rescale the objects by to ensure conservative advancement.
- Returns:¶
True if a collision was detected, false otherwise.
-
AdditiveCCD(const double conservative_rescaling
Inexact CCD¶
Note
This method is disabled by default. To enable it, set the
IPC_TOOLKIT_WITH_INEXACT_CCD
CMake option to ON
.
-
bool ipc::inexact_point_edge_ccd_2D(const Eigen::Vector2d& p_t0,
const Eigen::Vector2d& e0_t0, const Eigen::Vector2d& e1_t0,
const Eigen::Vector2d& p_t1, const Eigen::Vector2d& e0_t1,
const Eigen::Vector2d& e1_t1, double& toi,
const double conservative_rescaling)¶ Inexact continuous collision detection between a point and an edge in 2D.
- Parameters:¶
- const Eigen::Vector2d &p_t0¶
[in] The initial position of the point.
- const Eigen::Vector2d &e0_t0¶
[in] The initial position of the first endpoint of the edge.
- const Eigen::Vector2d &e1_t0¶
[in] The initial position of the second endpoint of the edge.
- const Eigen::Vector2d &p_t1¶
[in] The final position of the point.
- const Eigen::Vector2d &e0_t1¶
[in] The final position of the first endpoint of the edge.
- const Eigen::Vector2d &e1_t1¶
[in] The final position of the second endpoint of the edge.
- double &toi¶
[out] Output time of impact.
- const double conservative_rescaling¶
[in] The conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.
Nonlinear CCD¶
- class NonlinearTrajectory¶
Inheritence diagram for ipc::NonlinearTrajectory:
A nonlinear trajectory is a function that maps time to a point in space.
Subclassed by ipc::IntervalNonlinearTrajectory
Public Functions
- virtual ~NonlinearTrajectory() = default¶
-
class IntervalNonlinearTrajectory
: public virtual ipc::NonlinearTrajectory¶ Inheritence diagram for ipc::IntervalNonlinearTrajectory:
Collaboration diagram for ipc::IntervalNonlinearTrajectory:
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.
-
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 = TightInclusionCCD::DEFAULT_TOLERANCE,
const long max_iterations
= TightInclusionCCD::DEFAULT_MAX_ITERATIONS,
const double conservative_rescaling
= TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE¶
[in] Tolerance for the linear CCD algorithm
- const long max_iterations = TightInclusionCCD::DEFAULT_MAX_ITERATIONS¶
[in] Maximum number of iterations for the linear CCD algorithm
- const double conservative_rescaling = TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE,
const long max_iterations
= TightInclusionCCD::DEFAULT_MAX_ITERATIONS,
const double conservative_rescaling
= TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE¶
[in] Tolerance for the linear CCD algorithm
- const long max_iterations = TightInclusionCCD::DEFAULT_MAX_ITERATIONS¶
[in] Maximum number of iterations for the linear CCD algorithm
- const double conservative_rescaling = TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE,
const long max_iterations
= TightInclusionCCD::DEFAULT_MAX_ITERATIONS,
const double conservative_rescaling
= TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE¶
[in] Tolerance for the linear CCD algorithm
- const long max_iterations = TightInclusionCCD::DEFAULT_MAX_ITERATIONS¶
[in] Maximum number of iterations for the linear CCD algorithm
- const double conservative_rescaling = TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE,
const long max_iterations
= TightInclusionCCD::DEFAULT_MAX_ITERATIONS,
const double conservative_rescaling
= TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_TOLERANCE¶
[in] Tolerance for the linear CCD algorithm
- const long max_iterations = TightInclusionCCD::DEFAULT_MAX_ITERATIONS¶
[in] Maximum number of iterations for the linear CCD algorithm
- const double conservative_rescaling = TightInclusionCCD::DEFAULT_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
= TightInclusionCCD::DEFAULT_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 = TightInclusionCCD::DEFAULT_CONSERVATIVE_RESCALING¶
[in] Conservative rescaling of the time of impact.
- Returns:¶
Miscellaneous¶
-
bool ipc::point_static_plane_ccd(const VectorMax3d& p_t0,
const VectorMax3d& p_t1, const VectorMax3d& plane_origin,
const VectorMax3d& plane_normal, double& toi,
const double conservative_rescaling
= TightInclusionCCD::DEFAULT_CONSERVATIVE_RESCALING)¶ Computes the time of impact between a point and a static plane in 3D using continuous collision detection.
- Parameters:¶
- const VectorMax3d &p_t0¶
[in] The initial position of the point.
- const VectorMax3d &p_t1¶
[in] The final position of the point.
- const VectorMax3d &plane_origin¶
[in] The origin of the plane.
- const VectorMax3d &plane_normal¶
[in] The normal of the plane.
- double &toi¶
[out] Output time of impact.
- const double conservative_rescaling = TightInclusionCCD::DEFAULT_CONSERVATIVE_RESCALING¶
[in] Conservative rescaling of the time of impact.
- Returns:¶
True if a collision was detected, false otherwise.