Adhesion¶
Adhesion simulates the sticking interaction between surfaces, such as glue bonding or material contact.
We provide functionality to compute normal adhesion (perpendicular forces) and tangential adhesion (parallel forces) using the model of Fang et al. [2024].
Normal Adhesion¶
The normal adhesion potential models the attraction force based on the distance between two surfaces.
For a given distance \(d\):
and
The parameters are:
\(\hat{d}_p\) (
dhat_p
) is the threshold distance for adhesion (in units of distance) where the largest adhesion force is applied\(\hat{d}_a\) (
dhat_a
) is the adhesion activation distance (in units of distance), representing the maximum range where adhesion forces are active\(Y\) (
Y
) is the adhesion stiffness (in units of stress, such as Young’s modulus), controlling the intensity of adhesion forces\(\epsilon_c\) (
eps_c
) is the adhesion coefficient (unitless) that defines the critical strain at which adhesion forces decrease
We can build a normal adhesion potential object and compute the adhesion potential for a given set of normal collisions.
const double dhat_p = 1e-3;
const double dhat_a = 2 * dhat_p;
const double Y = 1e3;
const double eps_c = 0.5;
const ipc::NormalAdhesionPotential A_n(dhat_p, dhat_a, Y, eps_c)
double adhesion_potential = A_n(normal_collisions, collision_mesh, vertices);
dhat_p = 1e-3
dhat_a = 2 * dhat_p
Y = 1e3
eps_c = 0.5
A_n = ipctk.NormalAdhesionPotential(dhat_p, dhat_a, Y, eps_c)
adhesion_potential = A_n(normal_collisions, collision_mesh, vertices)
Derivatives¶
We can also compute the first and second derivatives of the normal adhesion potential with respect to the vertex positions.
Eigen::VectorXd adhesion_potential_grad =
A_n.gradient(normal_collisions, collision_mesh, vertices);
Eigen::SparseMatrix<double> adhesion_potential_hess =
A_n.hessian(normal_collisions, collision_mesh, vertices);
adhesion_potential_grad = A_n.gradient(normal_collisions, collision_mesh, vertices)
adhesion_potential_hess = A_n.hessian(normal_collisions, collision_mesh, vertices)
Tangential Adhesion¶
The tangential adhesion potential models resistance to sliding (parallel to surfaces).
It is structured similar to the friction model with a smooth transition between sticking and sliding and lagged normal forces and tangential bases.
where \(C\) is the lagged collisions, \(\lambda_k^a\) is the normal adhesion force magnitude for the \(k\)-th collision, \(\mathbf{T}_k\) is the tangential basis for the \(k\)-th collision, and \(f_0^a\) is the smooth tangential adhesion function used to approximate the non-smooth transition from sticking to sliding.
For relative displacement magnitude \(y\):
where \(\epsilon_a\) (eps_a
) is the adhesion threshold (in units of displacement) used to smoothly transition from sticking to sliding.
We can build a tangential adhesion potential object and compute the adhesion potential for a given set of tangential collisions.
const double eps_a = 0.01;
const ipc::TangentialAdhesionPotential A_t(eps_a);
double adhesion_potential = A_t(tangential_collisions, collision_mesh, displacement);
eps_a = 0.01
A_t = ipctk.TangentialAdhesionPotential(eps_a)
adhesion_potential = A_t(tangential_collisions, collision_mesh, displacement);
Derivatives¶
We can also compute the first and second derivatives of the tangential adhesion potential with respect to the displacement.
Eigen::VectorXd adhesion_potential_grad =
A_t.gradient(tangential_collisions, collision_mesh, displacement);
Eigen::SparseMatrix<double> adhesion_potential_hess =
A_t.hessian(tangential_collisions, collision_mesh, displacement);
adhesion_potential_grad = A_t.gradient(tangential_collisions, collision_mesh, displacement)
adhesion_potential_hess = A_t.hessian(tangential_collisions, collision_mesh, displacement)