TangentialPotential

class TangentialPotential
   
 : public ipc::Potential<TangentialCollisions>;

A tangential dissipative potential.

Subclassed by ipc::FrictionPotential, ipc::TangentialAdhesionPotential

Public Types

enum class DiffWRT : uint8_t;

Variable to differentiate the friction force with respect to.

Values:

enumerator REST_POSITIONS;

Differentiate w.r.t. rest positions.

enumerator LAGGED_DISPLACEMENTS;

Differentiate w.r.t. lagged displacements.

enumerator VELOCITIES;

Differentiate w.r.t. current velocities.

Public Functions

inline virtual ~TangentialPotential();
Eigen::VectorXd force(const TangentialCollisionscollisions,
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
    
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
    
Eigen::ConstRef<Eigen::MatrixXd> velocities,
    
const NormalPotentialnormal_potential,
    
const double normal_stiffnessconst double dmin = 0,
    
const bool no_mu = false) const;

Compute the friction force from the given velocities.

Parameters:
const TangentialCollisions &collisions

The set of collisions.

const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> rest_positions

Rest positions of the vertices (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements

Previous displacements of the vertices (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> velocities

Current displacements of the vertices (rowwise).

const NormalPotential &normal_potential

Normal potential (used for normal force magnitude).

const double normal_stiffness

Normal stiffness (used for normal force magnitude).

const double dmin = 0

Minimum distance (used for normal force magnitude).

const bool no_mu = false

whether to not multiply by mu

Returns:

The friction force.

Eigen::SparseMatrix<double> force_jacobian(
    
const TangentialCollisionscollisions,
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
    
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
    
Eigen::ConstRef<Eigen::MatrixXd> velocities,
    
const NormalPotentialnormal_potential,
    
const double normal_stiffnessconst DiffWRT wrt,
    
const double dmin = 0) const;

Compute the Jacobian of the friction force wrt the velocities.

Parameters:
const TangentialCollisions &collisions

The set of collisions.

const CollisionMesh &mesh

The collision mesh.

Eigen::ConstRef<Eigen::MatrixXd> rest_positions

Rest positions of the vertices (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements

Previous displacements of the vertices (rowwise).

Eigen::ConstRef<Eigen::MatrixXd> velocities

Current displacements of the vertices (rowwise).

const NormalPotential &normal_potential

Normal potential (used for normal force magnitude).

const double normal_stiffness

Normal stiffness (used for normal force magnitude).

const DiffWRT wrt

The variable to take the derivative with respect to.

const double dmin = 0

Minimum distance (used for normal force magnitude).

Returns:

The Jacobian of the friction force wrt the velocities.

Eigen::VectorXd smooth_contact_force(
    
const TangentialCollisionscollisions,
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
    
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
    
Eigen::ConstRef<Eigen::MatrixXd> velocities,
    
const double dmin = 0const bool no_mu = false) const;
Eigen::SparseMatrix<double> smooth_contact_force_jacobian(
    
const TangentialCollisionscollisions,
    
const CollisionMeshmesh,
    
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
    
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
    
Eigen::ConstRef<Eigen::MatrixXd> velocities,
    
const SmoothContactParametersparamsconst DiffWRT wrt,
    
const double dmin = 0const bool no_mu = false) const;
double operator()(const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMax12d> velocities) const override;

Compute the potential for a single collision.

Parameters:
const TangentialCollision &collision

The collision

Eigen::ConstRef<VectorMax12d> velocities

The collision stencil’s velocities.

Returns:

The potential.

VectorMax12d gradient(const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMax12d> velocities) const override;

Compute the gradient of the potential for a single collision.

Parameters:
const TangentialCollision &collision

The collision

Eigen::ConstRef<VectorMax12d> velocities

The collision stencil’s velocities.

Returns:

The gradient of the potential.

MatrixMax12d hessian(const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMax12d> velocities,
    
const PSDProjectionMethod project_hessian_to_psd
    
= PSDProjectionMethod::NONE
) const override;

Compute the hessian of the potential for a single collision.

Parameters:
const TangentialCollision &collision

The collision

Eigen::ConstRef<VectorMax12d> velocities

The collision stencil’s velocities.

const PSDProjectionMethod project_hessian_to_psd = PSDProjectionMethod::NONE

Whether to project the hessian to the positive semi-definite cone.

Returns:

The hessian of the potential.

VectorMax12d force(const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMax12d> rest_positions,
    
Eigen::ConstRef<VectorMax12d> lagged_displacements,
    
Eigen::ConstRef<VectorMax12d> velocities,
    
const NormalPotentialnormal_potential,
    
const double normal_stiffnessconst double dmin = 0,
    
const bool no_mu = false) const;

Compute the friction force.

Parameters:
const TangentialCollision &collision

The collision

Eigen::ConstRef<VectorMax12d> rest_positions

Rest positions of the vertices (rowwise).

Eigen::ConstRef<VectorMax12d> lagged_displacements

Previous displacements of the vertices (rowwise).

Eigen::ConstRef<VectorMax12d> velocities

Current displacements of the vertices (rowwise).

const NormalPotential &normal_potential

Normal potential (used for normal force magnitude).

const double normal_stiffness

Normal stiffness (used for normal force magnitude).

const double dmin = 0

Minimum distance (used for normal force magnitude).

const bool no_mu = false

Whether to not multiply by mu

Returns:

Friction force

MatrixMax12d force_jacobian(const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMax12d> rest_positions,
    
Eigen::ConstRef<VectorMax12d> lagged_displacements,
    
Eigen::ConstRef<VectorMax12d> velocities,
    
const NormalPotentialnormal_potential,
    
const double normal_stiffnessconst DiffWRT wrt,
    
const double dmin = 0) const;

Compute the friction force Jacobian.

Parameters:
const TangentialCollision &collision

The collision

Eigen::ConstRef<VectorMax12d> rest_positions

Rest positions of the vertices (rowwise).

Eigen::ConstRef<VectorMax12d> lagged_displacements

Previous displacements of the vertices (rowwise).

Eigen::ConstRef<VectorMax12d> velocities

Current displacements of the vertices (rowwise).

const NormalPotential &normal_potential

Normal potential (used for normal force magnitude).

const double normal_stiffness

Noraml stiffness (used for normal force magnitude).

const DiffWRT wrt

Variable to differentiate the friction force with respect to.

const double dmin = 0

Minimum distance (used for normal force magnitude).

Returns:

Friction force Jacobian

VectorMaxNd smooth_contact_force(
    
const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMaxNd> rest_positions,
    
Eigen::ConstRef<VectorMaxNd> lagged_displacements,
    
Eigen::ConstRef<VectorMaxNd> velocities,
    
const bool no_mu = false,
    
const bool no_contact_force_multiplier = false) const;
Eigen::MatrixXd smooth_contact_force_jacobian(
    
const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMaxNd> rest_positions,
    
Eigen::ConstRef<VectorMaxNd> lagged_displacements,
    
Eigen::ConstRef<VectorMaxNd> velocitiesconst DiffWRT wrt,
    
const bool no_mu) const;
MatrixMaxNd smooth_contact_force_jacobian_unit(
    
const TangentialCollisioncollision,
    
Eigen::ConstRef<VectorMaxNd> lagged_positions,
    
Eigen::ConstRef<VectorMaxNd> velocitiesconst DiffWRT wrt,
    
const bool no_mu = false) const;

Compute the friction force Jacobian assuming the contact force magnitude being 1.

Parameters:
const TangentialCollision &collision

The collision

rest_positions

Rest positions of the vertices (rowwise).

lagged_displacements

Previous displacements of the vertices (rowwise).

Eigen::ConstRef<VectorMaxNd> velocities

Current displacements of the vertices (rowwise).

barrier_potential

Barrier potential (used for normal force magnitude).

barrier_stiffness

Barrier stiffness (used for normal force magnitude).

const DiffWRT wrt

Variable to differentiate the friction force with respect to.

dmin

Minimum distance (used for normal force magnitude).

Returns:

Friction force Jacobian

Protected Functions

virtual double mu_f0(
    
const double xconst double mu_sconst double mu_k) const
   
 = 0;

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

Parameters:
const double x

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

Returns:

The value of the integral at x.

virtual double mu_f1_over_x(
    
const double xconst double mu_sconst double mu_k) const
   
 = 0;

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

Parameters:
const double x

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

Returns:

The value of the product at x.

virtual double mu_f2_x_minus_mu_f1_over_x3(
    
const double xconst double mu_sconst double mu_k) const
   
 = 0;

Compute the value of [(d/dx (μ(y) f₁(y))) x - μ(y) f₁(y)] / x³, where f₁ is the first derivative of the smooth mollifier.

Parameters:
const double x

The tangential relative speed.

const double mu_s

Coefficient of static friction.

const double mu_k

Coefficient of kinetic friction.

Returns:

The value of the derivative at x.

virtual bool is_dynamic(const double speed) const = 0;

Check if the speed is dynamic.

Parameters:
const double speed

The tangential relative speed.

Returns:

True if the speed is dynamic, false otherwise.

Private Types

using Super = Potential<TangentialCollisions>;
using VectorMaxNd = Super::VectorMaxNd;
using MatrixMaxNd = Super::MatrixMaxNd;