Potentials¶
Generic Potential¶
- template <class TCollisions> class Potential¶
Base class for potentials.
Public Functions
- Potential() = default¶
- virtual ~Potential() = default¶
-
double operator()(const TCollisions& collisions,
const CollisionMesh& mesh, const Eigen::MatrixXd& X) const¶ Compute the potential for a set of collisions.
- Parameters:¶
- const TCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &X¶
Degrees of freedom of the collision mesh (e.g., vertices or velocities).
- Returns:¶
The potential for a set of collisions.
-
Eigen::VectorXd gradient(const TCollisions& collisions,
const CollisionMesh& mesh, const Eigen::MatrixXd& X) const¶ Compute the gradient of the potential.
- Parameters:¶
- const TCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- const 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|.
-
Eigen::SparseMatrix<double> hessian(const TCollisions& collisions,
const CollisionMesh& mesh, const 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.
- const 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 double operator()(
const TCollision& collision, const VectorMax12d& x) const = 0¶ Compute the potential for a single collision.
- Parameters:¶
- const TCollision &collision¶
The collision.
- const VectorMax12d &x¶
The collision stencil’s degrees of freedom.
- Returns:¶
The potential.
-
virtual VectorMax12d gradient(
const TCollision& collision, const VectorMax12d& x) const = 0¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const TCollision &collision¶
The collision.
- const VectorMax12d &x¶
The collision stencil’s degrees of freedom.
- Returns:¶
The gradient of the potential.
-
virtual MatrixMax12d hessian(const TCollision& collision,
const 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.
- const VectorMax12d &x¶
The collision stencil’s degrees of freedom.
- Returns:¶
The hessian of the potential.
Protected Types
- using TCollision = typename TCollisions::value_type¶
Normal Potentials¶
- class NormalPotential : public ipc::Potential<NormalCollisions>¶
Inheritence diagram for ipc::NormalPotential:
Collaboration diagram for ipc::NormalPotential:
Base class for distance-based potentials.
Subclassed by ipc::BarrierPotential, ipc::NormalAdhesionPotential
Public Functions
- NormalPotential() = default¶
- virtual ~NormalPotential() = default¶
-
Eigen::SparseMatrix<double> shape_derivative(
const NormalCollisions& collisions, const CollisionMesh& mesh,
const Eigen::MatrixXd& vertices) const¶ Compute the shape derivative of the potential.
- Parameters:¶
- const NormalCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- const Eigen::MatrixXd &vertices¶
Vertices of the collision mesh.
- Throws:¶
std::runtime_error
– If the collision collisions were not built with shape derivatives enabled.- Returns:¶
The derivative of the force with respect to X, the rest vertices.
-
double operator()(const NormalCollision& collision,
const VectorMax12d& positions) const override¶ Compute the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The potential.
-
VectorMax12d gradient(const NormalCollision& collision,
const VectorMax12d& positions) const override¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The gradient of the potential.
-
MatrixMax12d hessian(const NormalCollision& collision,
const VectorMax12d& positions,
const PSDProjectionMethod project_hessian_to_psd
= PSDProjectionMethod::NONE) const override¶ Compute the hessian of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The hessian of the potential.
-
void shape_derivative(const NormalCollision& collision,
const std::array<long, 4>& vertex_ids,
const VectorMax12d& rest_positions,
const VectorMax12d& positions,
std::vector<Eigen::Triplet<double>>& out) const¶ Compute the shape derivative of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
[in] The collision.
- const std::array<long, 4> &vertex_ids¶
[in] The collision stencil’s vertex ids.
- const VectorMax12d &rest_positions¶
[in] The collision stencil’s rest positions.
- const VectorMax12d &positions¶
[in] The collision stencil’s positions.
- std::vector<Eigen::Triplet<double>> &out¶
[inout] Store the triplets of the shape derivative here.
-
virtual double force_magnitude(const double distance_squared,
const double dmin, const double barrier_stiffness) const = 0¶ Compute the force magnitude for a collision.
-
virtual VectorMax12d force_magnitude_gradient(
const double distance_squared,
const VectorMax12d& distance_squared_gradient,
const double dmin, const double barrier_stiffness) const = 0¶ Compute the gradient of the force magnitude for a collision.
Protected Functions
-
virtual double operator()(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the unmollified distance-based potential for a collisions.
-
virtual double gradient(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the gradient of the unmollified distance-based potential for a collision.
-
virtual double hessian(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the hessian of the unmollified distance-based potential for a collision.
Private Types
- using Super = Potential<NormalCollisions>¶
Barrier Potential¶
- class BarrierPotential : public ipc::NormalPotential¶
Inheritence diagram for ipc::BarrierPotential:
Collaboration diagram for ipc::BarrierPotential:
The barrier collision potential.
Public Functions
-
explicit BarrierPotential(
const double dhat, const bool use_physical_barrier = false)¶ Construct a barrier potential.
Construct a barrier potential.
The barrier function.
The activation distance of the barrier.
- inline double dhat() const¶
Get the activation distance of the barrier.
Set the barrier function used to compute the potential.
The barrier function used to compute the potential.
-
virtual double force_magnitude(const double distance_squared,
const double dmin,
const double barrier_stiffness) const override¶ Compute the force magnitude for a collision.
-
virtual VectorMax12d force_magnitude_gradient(
const double distance_squared,
const VectorMax12d& distance_squared_gradient,
const double dmin,
const double barrier_stiffness) const override¶ Compute the gradient of the force magnitude for a collision.
- inline bool use_physical_barrier() const¶
Get whether to use the physical barrier.
Note
When using the convergent formulation we want the barrier to have units of Pa⋅m, so κ gets units of Pa and the barrier function should have units of m. See notebooks/physical_barrier.ipynb for more details.
- inline void set_use_physical_barrier(bool use_physical_barrier)¶
Set use physical barrier flag.
Note
When using the convergent formulation we want the barrier to have units of Pa⋅m, so κ gets units of Pa and the barrier function should have units of m. See notebooks/physical_barrier.ipynb for more details.
-
VectorMax12d gradient(const NormalCollision& collision,
const VectorMax12d& positions) const override¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The gradient of the potential.
-
double gradient(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the gradient of the unmollified distance-based potential for a collision.
-
MatrixMax12d hessian(const NormalCollision& collision,
const VectorMax12d& positions,
const PSDProjectionMethod project_hessian_to_psd
= PSDProjectionMethod::NONE) const override¶ Compute the hessian of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The hessian of the potential.
-
double hessian(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the hessian of the unmollified distance-based potential for a collision.
Protected Functions
-
virtual double operator()(const double distance_squared,
const double dmin = 0) const override¶ Compute the barrier potential for a collision.
-
virtual double gradient(const double distance_squared,
const double dmin = 0) const override Compute the gradient of the barrier potential for a collision.
-
virtual double hessian(const double distance_squared,
const double dmin = 0) const override Compute the hessian of the barrier potential for a collision.
Protected Attributes
-
std::shared_ptr<Barrier> m_barrier
= std::make_shared<ClampedLogBarrier>()¶ The barrier function used to compute the potential.
- double m_dhat¶
The activation distance of the barrier.
- bool m_use_physical_barrier = false¶
Whether to use the physical barrier.
Note
When using the convergent formulation we want the barrier to have units of Pa⋅m, so κ gets units of Pa and the barrier function should have units of m. See notebooks/physical_barrier.ipynb for more details.
Private Types
- using Super = NormalPotential¶
-
explicit BarrierPotential(
Normal Adhesion Potential¶
- class NormalAdhesionPotential : public ipc::NormalPotential¶
Inheritence diagram for ipc::NormalAdhesionPotential:
Collaboration diagram for ipc::NormalAdhesionPotential:
The normal adhesion potential.
Public Functions
-
NormalAdhesionPotential(const double _dhat_p, const double _dhat_a,
const double _Y, const double _eps_c)¶
-
virtual double force_magnitude(const double distance_squared,
const double dmin,
const double barrier_stiffness) const override¶ Compute the force magnitude for a collision.
-
virtual VectorMax12d force_magnitude_gradient(
const double distance_squared,
const VectorMax12d& distance_squared_gradient,
const double dmin,
const double barrier_stiffness) const override¶ Compute the gradient of the force magnitude for a collision.
-
VectorMax12d gradient(const NormalCollision& collision,
const VectorMax12d& positions) const override¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The gradient of the potential.
-
double gradient(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the gradient of the unmollified distance-based potential for a collision.
-
MatrixMax12d hessian(const NormalCollision& collision,
const VectorMax12d& positions,
const PSDProjectionMethod project_hessian_to_psd
= PSDProjectionMethod::NONE) const override¶ Compute the hessian of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- const VectorMax12d &positions¶
The collision stencil’s positions.
- Returns:¶
The hessian of the potential.
-
double hessian(
const double distance_squared, const double dmin = 0) const = 0¶ Compute the hessian of the unmollified distance-based potential for a collision.
Public Members
- double dhat_p¶
The distance of largest adhesion force ( \(\hat{d}_p\)) ( \(0 < \hat{d}_p < \hat{d}_a\)).
- double dhat_a¶
The adhesion activation distance ( \(\hat{d}_a\)).
- double Y¶
The Young’s modulus ( \(Y\)).
- double eps_c¶
The critical strain ( \(\varepsilon_c\)).
Protected Functions
-
virtual double operator()(const double distance_squared,
const double dmin = 0) const override¶ Compute the barrier potential for a collision.
-
virtual double gradient(const double distance_squared,
const double dmin = 0) const override Compute the gradient of the barrier potential for a collision.
-
virtual double hessian(const double distance_squared,
const double dmin = 0) const override Compute the hessian of the barrier potential for a collision.
Private Types
- using Super = NormalPotential¶
-
NormalAdhesionPotential(const double _dhat_p, const double _dhat_a,
Tangential Potentials¶
-
class TangentialPotential
: public ipc::Potential<TangentialCollisions>¶ Inheritence diagram for ipc::TangentialPotential:
Collaboration diagram for ipc::TangentialPotential:
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,
const Eigen::MatrixXd& rest_positions,
const Eigen::MatrixXd& lagged_displacements,
const 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.
- const Eigen::MatrixXd &rest_positions¶
Rest positions of the vertices (rowwise).
- const Eigen::MatrixXd &lagged_displacements¶
Previous displacements of the vertices (rowwise).
- const 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,
const Eigen::MatrixXd& rest_positions,
const Eigen::MatrixXd& lagged_displacements,
const 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.
- const Eigen::MatrixXd &rest_positions¶
Rest positions of the vertices (rowwise).
- const Eigen::MatrixXd &lagged_displacements¶
Previous displacements of the vertices (rowwise).
- const 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,
const VectorMax12d& velocities) const override¶ Compute the potential for a single collision.
- Parameters:¶
- const TangentialCollision &collision¶
The collision
- const VectorMax12d &velocities¶
The collision stencil’s velocities.
- Returns:¶
The potential.
-
VectorMax12d gradient(const TangentialCollision& collision,
const VectorMax12d& velocities) const override¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const TangentialCollision &collision¶
The collision
- const VectorMax12d &velocities¶
The collision stencil’s velocities.
- Returns:¶
The gradient of the potential.
-
MatrixMax12d hessian(const TangentialCollision& collision,
const 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
- const VectorMax12d &velocities¶
The collision stencil’s velocities.
- Returns:¶
The hessian of the potential.
-
VectorMax12d force(const TangentialCollision& collision,
const VectorMax12d& rest_positions,
const VectorMax12d& lagged_displacements,
const 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
- const VectorMax12d &rest_positions¶
Rest positions of the vertices (rowwise).
- const VectorMax12d &lagged_displacements¶
Previous displacements of the vertices (rowwise).
- const 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,
const VectorMax12d& rest_positions,
const VectorMax12d& lagged_displacements,
const 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
- const VectorMax12d &rest_positions¶
Rest positions of the vertices (rowwise).
- const VectorMax12d &lagged_displacements¶
Previous displacements of the vertices (rowwise).
- const 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, const Eigen::MatrixXd& X) const¶ Compute the gradient of the potential.
- Parameters:¶
- const TCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- const 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|.
-
VectorMax12d gradient(
const TCollision& collision, const VectorMax12d& x) const = 0¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const TCollision &collision¶
The collision.
- const 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, const 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.
- const 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|.
-
MatrixMax12d hessian(const TCollision& collision,
const 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.
- const VectorMax12d &x¶
The collision stencil’s degrees of freedom.
- Returns:¶
The hessian of the potential.
Protected Functions
Friction Potential¶
- class FrictionPotential : public ipc::TangentialPotential¶
Inheritence diagram for ipc::FrictionPotential:
Collaboration diagram for ipc::FrictionPotential:
The friction dissipative potential.
Public Functions
- inline double eps_v() const¶
Get the smooth friction mollifier parameter \(\epsilon_v\).
Protected Functions
Protected Attributes
- double m_eps_v¶
The smooth friction mollifier parameter \(\epsilon_v\).
Private Types
- using Super = TangentialPotential¶
Tangential Adhesion Potential¶
- class TangentialAdhesionPotential : public ipc::TangentialPotential¶
Inheritence diagram for ipc::TangentialAdhesionPotential:
Collaboration diagram for ipc::TangentialAdhesionPotential:
The tangential adhesion potential.
Public Functions
- explicit TangentialAdhesionPotential(const double eps_a)¶
Construct a tangential adhesion potential.
- inline double eps_a() const¶
Get the tangential adhesion mollifier parameter \(\epsilon_a\).
Protected Functions
Protected Attributes
- double m_eps_a¶
The tangential adhesion mollifier parameter \(\epsilon_a\).
Private Types
- using Super = TangentialPotential¶