Barrier

double ipc::barrier(const double dconst double dhat);

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

double ipc::barrier_first_derivative(
    
const double dconst double dhat);

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The derivative of the barrier wrt d.

double ipc::barrier_second_derivative(
    
const double dconst double dhat);

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The second derivative of the barrier wrt d.

Barrier Force Magnitude

double ipc::barrier_force_magnitude(const double distance_squared,
    
const Barrierbarrierconst double dhat,
    
const double barrier_stiffnessconst double dmin = 0);

Compute the magnitude of the force due to a barrier.

Parameters:
const double distance_squared

The squared distance between elements.

const Barrier &barrier

The barrier function.

const double dhat

The activation distance of the barrier.

const double barrier_stiffness

The stiffness of the barrier.

const double dmin = 0

The minimum distance offset to the barrier.

Returns:

The magnitude of the force.

VectorMax12d ipc::barrier_force_magnitude_gradient(
    
const double distance_squared,
    
const VectorMax12ddistance_squared_gradient,
    
const Barrierbarrierconst double dhat,
    
const double barrier_stiffnessconst double dmin = 0);

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

Parameters:
const double distance_squared

The squared distance between elements.

const VectorMax12d &distance_squared_gradient

The gradient of the squared distance.

const Barrier &barrier

The barrier function.

const double dhat

The activation distance of the barrier.

const double barrier_stiffness

The stiffness of the barrier.

const double dmin = 0

The minimum distance offset to the barrier.

Returns:

The gradient of the force.

Adaptive Barrier Stiffness

double ipc::initial_barrier_stiffness(const double bbox_diagonal,
    
const Barrierbarrierconst double dhat,
    
const double average_massconst Eigen::VectorXdgrad_energy,
    
const Eigen::VectorXdgrad_barrier,
    
doublemax_barrier_stiffness,
    
const double min_barrier_stiffness_scale = 1e11,
    
const double dmin = 0);

Compute an inital barrier stiffness using the barrier potential gradient.

Parameters:
const double bbox_diagonal

[in] Length of the diagonal of the bounding box of the scene.

const Barrier &barrier

[in] Barrier function.

const double dhat

[in] Activation distance of the barrier.

const double average_mass

[in] Average mass of all bodies.

const Eigen::VectorXd &grad_energy

[in] Gradient of the elasticity energy function.

const Eigen::VectorXd &grad_barrier

[in] Gradient of the barrier potential.

double &max_barrier_stiffness

[out] Maximum stiffness of the barrier.

const double min_barrier_stiffness_scale = 1e11

[in] Scale used to premultiply the minimum barrier stiffness.

const double dmin = 0

[in] Minimum distance between elements.

Returns:

The initial barrier stiffness.

double ipc::update_barrier_stiffness(const double prev_min_distance,
    
const double min_distanceconst double max_barrier_stiffness,
    
const double barrier_stiffnessconst double bbox_diagonal,
    
const double dhat_epsilon_scale = 1e-9const double dmin = 0);

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

Parameters:
const double prev_min_distance

[in] Previous minimum distance between elements.

const double min_distance

[in] Current minimum distance between elements.

const double max_barrier_stiffness

[in] Maximum stiffness of the barrier.

const double barrier_stiffness

[in] Current barrier stiffness.

const double bbox_diagonal

[in] Length of the diagonal of the bounding box of the scene.

const double dhat_epsilon_scale = 1e-9

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

const double dmin = 0

[in] Minimum distance between elements.

Returns:

The updated barrier stiffness.

Semi-Implicit Stiffness

template <typename StencilsT>
Eigen::VectorXd ipc::semi_implicit_stiffness(
    
const CollisionMeshmeshconst Eigen::MatrixXdvertices,
    
const StencilsTcollisions,
    
const Eigen::VectorXdvertex_masses,
    
const Eigen::SparseMatrix<double>hessconst double dmin);

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

See [Ando 2024] for details.

Parameters:
const CollisionMesh &mesh

Collision mesh.

const Eigen::MatrixXd &vertices

Vertex positions.

const StencilsT &collisions

Normal collisions or collision candidates.

const Eigen::VectorXd &vertex_masses

Lumped vertex masses.

const Eigen::SparseMatrix<double> &hess

Hessian of the elasticity energy function.

const double dmin

Minimum distance between elements.

Returns:

The semi-implicit stiffness’s.

Barrier Class

class Barrier;

Inheritence diagram for ipc::Barrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "1" [label="ipc::Barrier" tooltip="ipc::Barrier" fillcolor="#BFBFBF"] "2" [label="ipc::ClampedLogBarrier" tooltip="ipc::ClampedLogBarrier"] "3" [label="ipc::ClampedLogSqBarrier" tooltip="ipc::ClampedLogSqBarrier"] "4" [label="ipc::CubicBarrier" tooltip="ipc::CubicBarrier"] "2" -> "1" [dir=forward tooltip="public-inheritance"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "1" [dir=forward tooltip="public-inheritance"] }

Base class for barrier functions.

Subclassed by ipc::ClampedLogBarrier, ipc::ClampedLogSqBarrier, ipc::CubicBarrier

Public Functions

Barrier() = default;
virtual ~Barrier() = default;
virtual double operator()(const double dconst double dhat) const
   
 = 0;

Evaluate the barrier function.

Parameters:
const double d

Distance.

const double dhat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

virtual double first_derivative(
    
const double dconst double dhat) const
   
 = 0;

Evaluate the first derivative of the barrier function wrt d.

Parameters:
const double d

Distance.

const double dhat

Activation distance of the barrier.

Returns:

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

virtual double second_derivative(
    
const double dconst double dhat) const
   
 = 0;

Evaluate the second derivative of the barrier function wrt d.

Parameters:
const double d

Distance.

const double dhat

Activation distance of the barrier.

Returns:

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

virtual double units(const double dhat) const = 0;

Get the units of the barrier function.

Parameters:
const double dhat

The activation distance of the barrier.

Returns:

Clamped Log Barrier

class ClampedLogBarrier : public ipc::Barrier;

Inheritence diagram for ipc::ClampedLogBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::ClampedLogBarrier" tooltip="ipc::ClampedLogBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::ClampedLogBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::ClampedLogBarrier" tooltip="ipc::ClampedLogBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

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

Public Functions

ClampedLogBarrier() = default;
inline virtual double operator()(
    
const double dconst double dhat) const override;

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

inline virtual double first_derivative(
    
const double dconst double dhat) const override;

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The derivative of the barrier wrt d.

inline virtual double second_derivative(
    
const double dconst double dhat) const override;

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:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The second derivative of the barrier wrt d.

inline virtual double units(const double dhat) const override;

Get the units of the barrier function.

Parameters:
const double dhat

The activation distance of the barrier.

Returns:

The units of the barrier function.

Normalized Clamped Log Barrier

Warning

doxygenclass: Cannot find class “ipc::NormalizedClampedLogBarrier” in doxygen xml output for project “IPC Toolkit” from directory: ../build/doxyoutput/xml

Clamped Log Squared Barrier

class ClampedLogSqBarrier : public ipc::Barrier;

Inheritence diagram for ipc::ClampedLogSqBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::ClampedLogSqBarrier" tooltip="ipc::ClampedLogSqBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::ClampedLogSqBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::ClampedLogSqBarrier" tooltip="ipc::ClampedLogSqBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

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

Public Functions

ClampedLogSqBarrier() = default;
virtual double operator()(
    
const double dconst double dhat) const override;

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

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

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

virtual double first_derivative(
    
const double dconst double dhat) const override;

Derivative of the barrier function.

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

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The derivative of the barrier wrt d.

virtual double second_derivative(
    
const double dconst double dhat) const override;

Second derivative of the barrier function.

\[ b''(d) = 2 \left(\ln^2\left(\frac{d}{\hat{d}}\right) - \left( \ln\left(\frac{d}{\hat{d}}\right) - 1\right) \frac{(\hat{d} - d)^2}{d^2} - 4 \ln\left(\frac{d}{\hat{d}}\right) \frac{\hat{d} - d}{d}\right) \]

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The second derivative of the barrier wrt d.

inline virtual double units(const double dhat) const override;

Get the units of the barrier function.

Parameters:
const double dhat

The activation distance of the barrier.

Returns:

The units of the barrier function.

Cubic Barrier

class CubicBarrier : public ipc::Barrier;

Inheritence diagram for ipc::CubicBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::CubicBarrier" tooltip="ipc::CubicBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::CubicBarrier:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Barrier" tooltip="ipc::Barrier"] "1" [label="ipc::CubicBarrier" tooltip="ipc::CubicBarrier" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Cubic barrier function from [Ando 2024].

Public Functions

CubicBarrier() = default;
virtual double operator()(
    
const double dconst double dhat) const override;

Weak barrier function.

\[ b(d) = -\frac{2}{3\hat{d}} (d - \hat{d})^3 \]

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The value of the barrier function at d.

virtual double first_derivative(
    
const double dconst double dhat) const override;

Derivative of the barrier function.

\[ b'(d) = -2 (d - \hat{d})^2 \]

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The derivative of the barrier wrt d.

virtual double second_derivative(
    
const double dconst double dhat) const override;

Second derivative of the barrier function.

\[ b''(d) = -4 (d - \hat{d}) \]

Parameters:
const double d

The distance.

const double dhat

Activation distance of the barrier.

Returns:

The second derivative of the barrier wrt d.

inline virtual double units(const double dhat) const override;

Get the units of the barrier function.

Parameters:
const double dhat

The activation distance of the barrier.

Returns:

The units of the barrier function.


Last update: Feb 18, 2025