Main Functions¶
- ipctk.construct_friction_constraint_set(*args, **kwargs)¶
Overloaded function.
construct_friction_constraint_set(mesh: ipc::CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], contact_constraint_set: ipc::Constraints, dhat: float, barrier_stiffness: float, mu: float) -> ipctk.FrictionConstraints
construct_friction_constraint_set(mesh: ipc::CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], contact_constraint_set: ipc::Constraints, dhat: float, barrier_stiffness: float, mus: numpy.ndarray[numpy.float64[m, 1]]) -> ipctk.FrictionConstraints
construct_friction_constraint_set(mesh: ipc::CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], contact_constraint_set: ipc::Constraints, dhat: float, barrier_stiffness: float, mus: numpy.ndarray[numpy.float64[m, 1]], blend_mu: std::function<double (double, double)>) -> ipctk.FrictionConstraints
- ipctk.compute_barrier_potential(mesh: ipctk.CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], constraint_set: ipctk.Constraints, dhat: float) float ¶
Compute the barrier potential for a given constraint set.
- Parameters:¶
- mesh: ipctk.CollisionMesh¶
The collision mesh.
- V: numpy.ndarray[numpy.float64[m, n]]¶
Vertices of the collision mesh.
- constraint_set: ipctk.Constraints¶
The set of constraints.
- dhat: float¶
The activation distance of the barrier.
- Returns:¶
The sum of all barrier potentials (not scaled by the barrier stiffness).
- ipctk.compute_barrier_potential_gradient(mesh: ipctk.CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], constraint_set: ipctk.Constraints, dhat: float) numpy.ndarray[numpy.float64[m, 1]] ¶
Compute the gradient of the barrier potential.
- Parameters:¶
- mesh: ipctk.CollisionMesh¶
The collision mesh.
- V: numpy.ndarray[numpy.float64[m, n]]¶
Vertices of the collision mesh.
- constraint_set: ipctk.Constraints¶
The set of constraints.
- dhat: float¶
The activation distance of the barrier.
- Returns:¶
The gradient of all barrier potentials (not scaled by the barrier stiffness). This will have a size of V.size.
-
ipctk.compute_barrier_potential_hessian(mesh: ipctk.CollisionMesh, V: numpy.ndarray[numpy.float64[m, n]], constraint_set: ipctk.Constraints, dhat: float, project_hessian_to_psd: bool =
True
) scipy.sparse.csc_matrix[numpy.float64] ¶ Compute the hessian of the barrier potential.
- Parameters:¶
- mesh: ipctk.CollisionMesh¶
The collision mesh.
- V: numpy.ndarray[numpy.float64[m, n]]¶
Vertices of the collision mesh.
- constraint_set: ipctk.Constraints¶
The set of constraints.
- dhat: float¶
The activation distance of the barrier.
- project_hessian_to_psd: bool =
True
¶ Make sure the hessian is positive semi-definite.
- Returns:¶
The hessian of all barrier potentials (not scaled by the barrier stiffness). This will have a shape of `(V.size, V.size).
- ipctk.compute_friction_potential(mesh: ipc::CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], friction_constraint_set: ipctk.FrictionConstraints, epsv_times_h: float) float ¶
Compute the friction potential between to positions.
- Parameters:¶
- V0
Vertex positions at start of time-step (rowwise)
- V1
Current vertex positions (rowwise)
- E
Edge vertex indicies
- F
Face vertex indicies (empty in 2D)
- friction_constraint_set
- epsv_times_h
- ipctk.compute_friction_potential_gradient(mesh: ipc::CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], friction_constraint_set: ipctk.FrictionConstraints, epsv_times_h: float) numpy.ndarray[numpy.float64[m, 1]] ¶
- ipctk.compute_friction_potential_hessian(mesh: ipc::CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], friction_constraint_set: ipctk.FrictionConstraints, epsv_times_h: float, project_hessian_to_psd: bool = True) scipy.sparse.csc_matrix[numpy.float64] ¶
- ipctk.is_step_collision_free(*args, **kwargs)¶
Overloaded function.
is_step_collision_free(mesh: ipctk.CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], method: ipctk.BroadPhaseMethod = <BroadPhaseMethod.HASH_GRID: 1>, min_distance: float = 0.0, tolerance: float = 1e-06, max_iterations: int = 10000000) -> bool
Determine if the step is collision free.
- Note:
Assumes the trajectory is linear.
- Parameters:
mesh: The collision mesh. V0: Surface vertex positions at start as rows of a matrix. V1: Surface vertex positions at end as rows of a matrix.
- Returns:
True if <b>any</b> collisions occur.
is_step_collision_free(candidates: ipctk.Candidates, mesh: ipctk.CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], min_distance: float = 0.0, tolerance: float = 1e-06, max_iterations: int = 10000000) -> bool
Determine if the step is collision free from a set of candidates.
- Note:
Assumes the trajectory is linear.
- Parameters:
candidates: Set of candidates to check for collisions. mesh: The collision mesh. V0: Surface vertex positions at start as rows of a matrix. V1: Surface vertex positions at end as rows of a matrix.
- Returns:
True if <b>any</b> collisions occur.
- ipctk.compute_collision_free_stepsize(*args, **kwargs)¶
Overloaded function.
compute_collision_free_stepsize(mesh: ipctk.CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], method: ipctk.BroadPhaseMethod = <BroadPhaseMethod.HASH_GRID: 1>, min_distance: float = 0.0, tolerance: float = 1e-06, max_iterations: int = 10000000) -> float
Computes a maximal step size that is collision free.
- Note:
Assumes the trajectory is linear.
- Parameters:
mesh: The collision mesh. V0: Vertex positions at start as rows of a matrix. Assumes V0 is intersection free. V1: Surface vertex positions at end as rows of a matrix.
- 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.
compute_collision_free_stepsize(candidates: ipctk.Candidates, mesh: ipctk.CollisionMesh, V0: numpy.ndarray[numpy.float64[m, n]], V1: numpy.ndarray[numpy.float64[m, n]], min_distance: float = 0.0, tolerance: float = 1e-06, max_iterations: int = 10000000) -> float
Computes a maximal step size that is collision free using a set of collision candidates.
- Note:
Assumes the trajectory is linear.
- Parameters:
candidates: Set of candidates to check for collisions. mesh: The collision mesh. V0: Vertex positions at start as rows of a matrix. Assumes V0 is intersection free. V1: Surface vertex positions at end as rows of a matrix.
- 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.