TangentialPotential¶
-
class TangentialPotential
: public ipc::Potential<TangentialCollisions>;¶ A tangential dissipative potential.
Subclassed by ipc::FrictionPotential, ipc::TangentialAdhesionPotential
Public Types¶
Public Functions¶
- inline virtual ~TangentialPotential();¶
-
Eigen::VectorXd force(const TangentialCollisions& collisions,
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
Eigen::ConstRef<Eigen::MatrixXd> velocities,
const NormalPotential& normal_potential,
const double normal_stiffness, const 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 TangentialCollisions& collisions,
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> rest_positions,
Eigen::ConstRef<Eigen::MatrixXd> lagged_displacements,
Eigen::ConstRef<Eigen::MatrixXd> velocities,
const NormalPotential& normal_potential,
const double normal_stiffness, const 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.
-
double operator()(const TangentialCollision& collision,
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 TangentialCollision& collision,
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 TangentialCollision& collision,
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 TangentialCollision& collision,
Eigen::ConstRef<VectorMax12d> rest_positions,
Eigen::ConstRef<VectorMax12d> lagged_displacements,
Eigen::ConstRef<VectorMax12d> velocities,
const NormalPotential& normal_potential,
const double normal_stiffness, const 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 TangentialCollision& collision,
Eigen::ConstRef<VectorMax12d> rest_positions,
Eigen::ConstRef<VectorMax12d> lagged_displacements,
Eigen::ConstRef<VectorMax12d> velocities,
const NormalPotential& normal_potential,
const double normal_stiffness, const 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
-
Eigen::VectorXd gradient(const TCollisions& collisions,
const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> X) const;¶ Compute the gradient of the potential.
- Parameters:¶
- const TCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> X¶
Degrees of freedom of the collision mesh (e.g., vertices or velocities).
- Returns:¶
The gradient of the potential w.r.t. X. This will have a size of |X|.
-
virtual VectorMax12d gradient(const TCollision& collision,
Eigen::ConstRef<VectorMax12d> x) const
= 0;¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const TCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> x¶
The collision stencil’s degrees of freedom.
- Returns:¶
The gradient of the potential.
-
Eigen::SparseMatrix<double> hessian(const TCollisions& collisions,
const CollisionMesh& mesh, Eigen::ConstRef<Eigen::MatrixXd> X,
const PSDProjectionMethod project_hessian_to_psd
= PSDProjectionMethod::NONE) const;¶ Compute the hessian of the potential.
- Parameters:¶
- const TCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> X¶
Degrees of freedom of the collision mesh (e.g., vertices or velocities).
- const PSDProjectionMethod project_hessian_to_psd = PSDProjectionMethod::NONE¶
Make sure the hessian is positive semi-definite.
- Returns:¶
The Hessian of the potential w.r.t. X. This will have a size of |X|×|X|.
-
virtual MatrixMax12d hessian(const TCollision& collision,
Eigen::ConstRef<VectorMax12d> x,
const PSDProjectionMethod project_hessian_to_psd
= PSDProjectionMethod::NONE) const
= 0;¶ Compute the hessian of the potential for a single collision.
- Parameters:¶
- const TCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> x¶
The collision stencil’s degrees of freedom.
- 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.
Protected Functions¶
Private Types¶
- using Super = Potential;¶