Barrier

ipctk.barrier(d: SupportsFloat, dhat: SupportsFloat) float

Function that grows to infinity as d approaches 0 from the right.

\[b(d) = -(d-\hat{d})^2\ln\left(\frac{d}{\hat{d}}\right)\]
Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

ipctk.barrier_first_derivative(d: SupportsFloat, dhat: SupportsFloat) float

Derivative of the barrier function.

\[b'(d) = (\hat{d}-d) \left( 2\ln\left( \frac{d}{\hat{d}} \right) - \frac{\hat{d}}{d} + 1\right)\]
Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The derivative of the barrier wrt d.

ipctk.barrier_second_derivative(d: SupportsFloat, dhat: SupportsFloat) float

Second derivative of the barrier function.

\[b''(d) = \left( \frac{\hat{d}}{d} + 2 \right) \frac{\hat{d}}{d} - 2\ln\left( \frac{d}{\hat{d}} \right) - 3\]
Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The second derivative of the barrier wrt d.

Barrier Force Magnitude

ipctk.barrier_force_magnitude(distance_squared: SupportsFloat, barrier: ipctk.Barrier, dhat: SupportsFloat, barrier_stiffness: SupportsFloat, dmin: SupportsFloat = 0) float

Compute the magnitude of the force due to a barrier.

Parameters:
distance_squared: SupportsFloat

The squared distance between elements.

barrier: ipctk.Barrier

The barrier function.

dhat: SupportsFloat

The activation distance of the barrier.

barrier_stiffness: SupportsFloat

The stiffness of the barrier.

dmin: SupportsFloat = 0

The minimum distance offset to the barrier.

Returns:

The magnitude of the force.

ipctk.barrier_force_magnitude_gradient(distance_squared: SupportsFloat, distance_squared_gradient: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], barrier: ipctk.Barrier, dhat: SupportsFloat, barrier_stiffness: SupportsFloat, dmin: SupportsFloat = 0) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']

Compute the gradient of the magnitude of the force due to a barrier.

Parameters:
distance_squared: SupportsFloat

The squared distance between elements.

distance_squared_gradient: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']

The gradient of the squared distance.

barrier: ipctk.Barrier

The barrier function.

dhat: SupportsFloat

The activation distance of the barrier.

barrier_stiffness: SupportsFloat

The stiffness of the barrier.

dmin: SupportsFloat = 0

The minimum distance offset to the barrier.

Returns:

The gradient of the force.

Adaptive Barrier Stiffness

ipctk.initial_barrier_stiffness(bbox_diagonal: SupportsFloat, barrier: ipctk.Barrier, dhat: SupportsFloat, average_mass: SupportsFloat, grad_energy: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], grad_barrier: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], min_barrier_stiffness_scale: SupportsFloat = 100000000000.0, dmin: SupportsFloat = 0) tuple[float, float]

Compute an inital barrier stiffness using the barrier potential gradient.

Parameters:
bbox_diagonal: SupportsFloat

Length of the diagonal of the bounding box of the scene.

barrier: ipctk.Barrier

Barrier function.

dhat: SupportsFloat

Activation distance of the barrier.

average_mass: SupportsFloat

Average mass of all bodies.

grad_energy: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']

Gradient of the elasticity energy function.

grad_barrier: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']

Gradient of the barrier potential.

min_barrier_stiffness_scale: SupportsFloat = 100000000000.0

Scale used to premultiply the minimum barrier stiffness.

dmin: SupportsFloat = 0

Minimum distance between elements.

Returns:

The initial barrier stiffness. Maximum stiffness of the barrier.

Return type:

Tuple of

ipctk.update_barrier_stiffness(prev_min_distance: SupportsFloat, min_distance: SupportsFloat, max_barrier_stiffness: SupportsFloat, barrier_stiffness: SupportsFloat, bbox_diagonal: SupportsFloat, dhat_epsilon_scale: SupportsFloat = 1e-09, dmin: SupportsFloat = 0) float

Update the barrier stiffness if the distance is decreasing and less than dhat_epsilon_scale * diag.

Parameters:
prev_min_distance: SupportsFloat

Previous minimum distance between elements.

min_distance: SupportsFloat

Current minimum distance between elements.

max_barrier_stiffness: SupportsFloat

Maximum stiffness of the barrier.

barrier_stiffness: SupportsFloat

Current barrier stiffness.

bbox_diagonal: SupportsFloat

Length of the diagonal of the bounding box of the scene.

dhat_epsilon_scale: SupportsFloat = 1e-09

Update if distance is less than this fraction of the diagonal.

dmin: SupportsFloat = 0

Minimum distance between elements.

Returns:

The updated barrier stiffness.

Semi-Implicit Stiffness

ipctk.semi_implicit_stiffness(*args, **kwargs)

Overloaded function.

  1. semi_implicit_stiffness(stencil: ipc::CollisionStencil, vertices: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”], mass: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”], local_hess: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, n]”, “flags.f_contiguous”], dmin: typing.SupportsFloat) -> float

    Compute the semi-implicit stiffness for a single collision.

    See [Ando 2024] for details.

    Parameters:

    stencil: Collision stencil. vertex_ids: Vertex indices of the collision. vertices: Vertex positions. mass: Vertex masses. local_hess: Local hessian of the elasticity energy function. dmin: Minimum distance between elements.

    Returns:

    The semi-implicit stiffness.

  2. semi_implicit_stiffness(mesh: ipc::CollisionMesh, vertices: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, n]”, “flags.f_contiguous”], collisions: ipc::NormalCollisions, vertex_masses: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”], hess: scipy.sparse.csc_matrix[numpy.float64], dmin: typing.SupportsFloat) -> typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”]

    Compute the semi-implicit stiffness’s for all collisions.

    See [Ando 2024] for details.

    Parameters:

    mesh: Collision mesh. vertices: Vertex positions. collisions: Normal collisions. vertex_masses: Lumped vertex masses. hess: Hessian of the elasticity energy function. dmin: Minimum distance between elements.

    Returns:

    The semi-implicit stiffness’s.

  3. semi_implicit_stiffness(mesh: ipc::CollisionMesh, vertices: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, n]”, “flags.f_contiguous”], collisions: ipc::Candidates, vertex_masses: typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”], hess: scipy.sparse.csc_matrix[numpy.float64], dmin: typing.SupportsFloat) -> typing.Annotated[numpy.typing.NDArray[numpy.float64], “[m, 1]”]

    Compute the semi-implicit stiffness’s for all collisions.

    See [Ando 2024] for details.

    Parameters:

    mesh: Collision mesh. vertices: Vertex positions. collisions: Collisions candidates. vertex_masses: Lumped vertex masses. hess: Hessian of the elasticity energy function. dmin: Minimum distance between elements.

    Returns:

    The semi-implicit stiffness’s.

Barrier Class

class ipctk.Barrier

Bases: pybind11_object

__call__(self, d: SupportsFloat, dhat: SupportsFloat) float

Evaluate the barrier function.

Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()
first_derivative(self, d: SupportsFloat, dhat: SupportsFloat) float

Evaluate the first derivative of the barrier function wrt d.

Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The value of the first derivative of the barrier function at d.

second_derivative(self, d: SupportsFloat, dhat: SupportsFloat) float

Evaluate the second derivative of the barrier function wrt d.

Parameters:
d: SupportsFloat

The distance.

dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The value of the second derivative of the barrier function at d.

units(self, dhat: SupportsFloat) float

Get the units of the barrier function.

Parameters:
dhat: SupportsFloat

Activation distance of the barrier.

Returns:

The units of the barrier function.

Clamped Log Barrier

class ipctk.ClampedLogBarrier

Bases: Barrier

Smoothly clamped log barrier functions from [Li et al. 2020].

\[b(d) = -(d-\hat{d})^2\ln\left(\frac{d}{\hat{d}}\right)\]
__annotations__ = {}
__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()

Normalized Clamped Log Barrier

class ipctk.NormalizedClampedLogBarrier

Bases: ClampedLogBarrier

Normalized barrier function from [Li et al. 2023].

\[b(d) = -\left(\frac{d}{\hat{d}}-1\right)^2\ln\left(\frac{d}{\hat{d}}\right)\]
__annotations__ = {}
__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()

Clamped Log Squared Barrier

class ipctk.ClampedLogSqBarrier

Bases: Barrier

Clamped log barrier with a quadratic log term from [Huang et al. 2024].

\[b(d) = (d-\hat{d})^2\ln^2\left(\frac{d}{\hat{d}}\right)\]
__annotations__ = {}
__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()

Cubic Barrier

class ipctk.CubicBarrier

Bases: Barrier

Cubic barrier function from [Ando 2024].

\[b(d) = -\frac{2}{3\hat{d}} (d - \hat{d})^3\]
__annotations__ = {}
__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()

Two-Stage Barrier

class ipctk.TwoStageBarrier

Bases: Barrier

Two-stage activation barrier from [Chen et al. 2025].

\[\begin{split}b(d) = \begin{cases} -\frac{\hat{d}^2}{4} \left(\ln\left(\frac{2d}{\hat{d}}\right) - \tfrac{1}{2}\right) & d < \frac{\hat{d}}{2}\\ \tfrac{1}{2} (\hat{d} - d)^2 & d < \hat{d}\\ 0 & d \ge \hat{d} \end{cases}\end{split}\]
__annotations__ = {}
__init__(self)
__module__ = 'ipctk'
_pybind11_conduit_v1_()