Friction

See also

Advanced Friction explains the friction model, the static/kinetic transition, and anisotropic usage. Full derivation and plots are in notebooks/anisotropic_friction_math.ipynb.

Smooth Mollifier

double ipc::smooth_friction_f0(const double yconst double eps_v);

Smooth friction mollifier function.

\[ f_0(y)= \begin{cases} -\frac{y^3}{3\epsilon_v^2} + \frac{y^2}{\epsilon_v} + \frac{\epsilon_v}{3}, & |y| < \epsilon_v \newline y, & |y| \geq \epsilon_v \end{cases} \]

Parameters:
const double y

The tangential relative speed.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the mollifier function at y.

double ipc::smooth_friction_f1(const double yconst double eps_v);

The first derivative of the smooth friction mollifier.

\[ f_1(y) = f_0'(y) = \begin{cases} -\frac{y^2}{\epsilon_v^2}+\frac{2 y}{\epsilon_v}, & |y| < \epsilon_v \newline 1, & |y| \geq \epsilon_v \end{cases} \]

Parameters:
const double y

The tangential relative speed.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the derivative of the smooth friction mollifier at y.

double ipc::smooth_friction_f2(const double yconst double eps_v);

The second derivative of the smooth friction mollifier.

\[ f_2(y) = f_0''(y) = \begin{cases} -\frac{2 y}{\epsilon_v^2}+\frac{2}{\epsilon_v}, & |y| < \epsilon_v \newline 0, & |y| \geq \epsilon_v \end{cases} \]

Parameters:
const double y

The tangential relative speed.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the second derivative of the smooth friction mollifier at y.

double ipc::smooth_friction_f1_over_x(
    
const double yconst double eps_v);

Compute the derivative of the smooth friction mollifier divided by y ( \(\frac{f_0'(y)}{y}\)).

\[ \frac{f_1(y)}{y} = \begin{cases} -\frac{y}{\epsilon_v^2}+\frac{2}{\epsilon_v}, & |y| < \epsilon_v \newline \frac{1}{y}, & |y| \geq \epsilon_v \end{cases} \]

Note

The x in the function name refers to the parameter y.

Parameters:
const double y

The tangential relative speed.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the derivative of smooth_friction_f0 divided by y.

double ipc::smooth_friction_f2_x_minus_f1_over_x3(
    
const double yconst double eps_v);

The derivative of f1 times y minus f1 all divided by y cubed.

\[ \frac{f_1'(y) y - f_1(y)}{y^3} = \begin{cases} -\frac{1}{y \epsilon_v^2}, & |y| < \epsilon_v \newline -\frac{1}{y^3}, & |y| \geq \epsilon_v \end{cases} \]

Note

The x in the function name refers to the parameter y.

Parameters:
const double y

The tangential relative speed.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The derivative of f1 times y minus f1 all divided by y cubed.

Smooth \(\mu\)

double ipc::smooth_mu(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Smooth coefficient from static to kinetic friction.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the μ at y.

double ipc::smooth_mu_derivative(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Compute the derivative of the smooth coefficient from static to kinetic friction.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the derivative at y.

double ipc::smooth_mu_f0(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Compute the value of the ∫ μ(y) f₁(y) dy, where f₁ is the first derivative of the smooth friction mollifier.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the integral at y.

double ipc::smooth_mu_f1(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Compute the value of the μ(y) f₁(y), where f₁ is the first derivative of the smooth friction mollifier.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the product at y.

double ipc::smooth_mu_f2(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Compute the value of d/dy (μ(y) f₁(y)), where f₁ is the first derivative of the smooth friction mollifier.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the derivative at y.

double ipc::smooth_mu_f1_over_x(const double yconst double mu_s,
    
const double mu_kconst double eps_v);

Compute the value of the μ(y) f₁(y) / y, where f₁ is the first derivative of the smooth friction mollifier.

Note

The x in the function name refers to the parameter y.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the product at y.

double ipc::smooth_mu_f2_x_minus_mu_f1_over_x3(const double y,
    
const double mu_sconst double mu_kconst double eps_v);

Compute the value of the [(d/dy μ(y) f₁(y)) ⋅ y - μ(y) f₁(y)] / y³, where f₁ and f₂ are the first and second derivatives of the smooth friction mollifier.

Note

The x in the function name refers to the parameter y.

Parameters:
const double y

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

const double eps_v

Velocity threshold below which static friction force is applied.

Returns:

The value of the expression at y.

Anisotropic Friction Helpers

Effective friction follows an elliptical \(L^2\) projection (matchstick cone):

\[\mu_{\text{eff}} = \sqrt{(\mu_0 t_0)^2 + (\mu_1 t_1)^2}\]

with \(t = \tau / \|\tau\|\).

  • Use anisotropic_mu_eff_from_tau_aniso when you have \(\tau_{\text{aniso}}\) and need \(\mu_s\), \(\mu_k\) for the smooth transition;

  • use anisotropic_mu_eff_f when you have the unit direction.

  • Zero mu_s_aniso and mu_k_aniso falls back to scalar \(\mu_s\), \(\mu_k\).

For tangential collision evaluation, effective directional coefficients are used as lagged scalars (stored on each collision), i.e., no in-evaluation \(\partial \mu_{\text{eff}} / \partial \tau\) terms are taken in the force or Jacobian/Hessian paths. Call TangentialCollisions::update_lagged_anisotropic_friction_coefficients(...) after build(...) and whenever lagged geometry/velocities used for slip are updated (typically each Newton iteration).

See Erleben et al. [2019] for the Matchstick model; code: erleben/matchstick.

std::pair<double, double> ipc::anisotropic_mu_eff_f(
    
Eigen::ConstRef<Eigen::Vector2d> tau_dir,
    
Eigen::ConstRef<Eigen::Vector2d> mu_s_aniso,
    
Eigen::ConstRef<Eigen::Vector2d> mu_k_aniso);

Elliptical L2 (matchstick cone) anisotropic friction.

Call anisotropic_x_from_tau_aniso, then anisotropic_mu_eff_f.

Compute effective friction coefficients for elliptical anisotropy (L2 projection):

\(\mu_{\text{eff}} = f(x) = \sqrt{(\mu_0 t_0)^2 + (\mu_1 t_1)^2}\) at direction \(x = \tau_{\text{dir}}\)

.

For anisotropic friction, the friction coefficient depends on the direction of tangential velocity. This function computes the effective friction coefficients along a given direction using the elliptical L2 projection model:

\(\mu_{\text{eff}} = \sqrt{(\mu_0 t_0)^2 + (\mu_1 t_1)^2}\), where \(t = \tau / \|\tau\|\) is the unit direction vector.

See also

Erleben et al., CGF 2019, DOI 10.1111/cgf.13885; https://github.com/erleben/matchstick

Note

The \(f\) in the name refers to the effective-μ formula; \(x\) is the unit direction.

Note

If mu_s_aniso and mu_k_aniso are zero vectors, this function returns (0, 0) as the direct ellipse-formula result. Isotropic behavior is handled by anisotropic_mu_eff_from_tau_aniso.

Parameters:
Eigen::ConstRef<Eigen::Vector2d> tau_dir

Unit 2D direction of tangential velocity (tau / ||tau||). Must be a unit vector.

Eigen::ConstRef<Eigen::Vector2d> mu_s_aniso

Static friction ellipse axes (2D vector). Each component represents the friction coefficient along the corresponding tangent basis direction.

Eigen::ConstRef<Eigen::Vector2d> mu_k_aniso

Kinetic friction ellipse axes (2D vector). Each component represents the friction coefficient along the corresponding tangent basis direction.

Returns:

A pair containing (mu_s_eff, mu_k_eff), the effective static and kinetic friction coefficients along the direction tau_dir.

Eigen::Vector2d ipc::anisotropic_x_from_tau_aniso(
    
Eigen::ConstRef<Eigen::Vector2d> tau_aniso);

Compute unit direction \(x = \tau_{\text{aniso}} / \|\tau_{\text{aniso}}\|\) from tau_aniso, handling edge cases.

Parameters:
Eigen::ConstRef<Eigen::Vector2d> tau_aniso

Anisotropically-scaled tangential velocity (2D vector).

Returns:

Unit direction vector \(x\). Returns (1, 0) if ||tau_aniso|| ≈ 0.

std::pair<double, double> ipc::anisotropic_mu_eff_from_tau_aniso(
    
Eigen::ConstRef<Eigen::Vector2d> tau_aniso,
    
Eigen::ConstRef<Eigen::Vector2d> mu_s_aniso,
    
Eigen::ConstRef<Eigen::Vector2d> mu_k_aniso,
    
const double mu_s_isotropicconst double mu_k_isotropic,
    
const bool no_mu = false);

Compute effective friction coefficients from tau_aniso, handling both anisotropic and isotropic cases.

This function encapsulates the logic for determining whether to use anisotropic or isotropic friction coefficients. If anisotropic friction is enabled (mu_s_aniso.squaredNorm() > 0), it computes the effective mu based on the direction of tau_aniso. Otherwise, it returns the isotropic coefficients.

Parameters:
Eigen::ConstRef<Eigen::Vector2d> tau_aniso

Anisotropically-scaled tangential velocity (2D vector).

Eigen::ConstRef<Eigen::Vector2d> mu_s_aniso

Static friction ellipse axes (2D vector). Zero vector indicates isotropic friction.

Eigen::ConstRef<Eigen::Vector2d> mu_k_aniso

Kinetic friction ellipse axes (2D vector). Zero vector indicates isotropic friction.

const double mu_s_isotropic

Isotropic static friction coefficient (used when anisotropic is disabled).

const double mu_k_isotropic

Isotropic kinetic friction coefficient (used when anisotropic is disabled).

const bool no_mu = false

If true, returns (1.0, 1.0) regardless of input coefficients.

Returns:

A pair containing (mu_s, mu_k) to use in friction calculations.