Potentials

Generic Potential

template <class TCollisionsclass Potential;

Base class for potentials.

Template Parameters:
class TCollisions

The type of the collisions.

Public Functions

Potential() = default;
virtual ~Potential() = default;
double operator()(const TCollisionscollisions,
    
const CollisionMeshmeshconst Eigen::MatrixXdX) 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 TCollisionscollisions,
    
const CollisionMeshmeshconst Eigen::MatrixXdX) 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 TCollisionscollisions,
    
const CollisionMeshmeshconst Eigen::MatrixXdX,
    
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 TCollisioncollisionconst VectorMax12dx) 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 TCollisioncollisionconst VectorMax12dx) 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 TCollisioncollision,
    
const VectorMax12dx,
    
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;

Distance-based Potential

class DistanceBasedPotential : public ipc::Potential<Collisions>;

Inheritance diagram for ipc::DistanceBasedPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Potential< Collisions >" tooltip="ipc::Potential< Collisions >"] "3" [label="ipc::BarrierPotential" tooltip="ipc::BarrierPotential"] "1" [label="ipc::DistanceBasedPotential" tooltip="ipc::DistanceBasedPotential" fillcolor="#BFBFBF"] "3" -> "1" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::DistanceBasedPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Potential< Collisions >" tooltip="ipc::Potential< Collisions >"] "1" [label="ipc::DistanceBasedPotential" tooltip="ipc::DistanceBasedPotential" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Base class for distance-based potentials.

Subclassed by ipc::BarrierPotential

Public Functions

DistanceBasedPotential() = default;
virtual ~DistanceBasedPotential() = default;
Eigen::SparseMatrix<double> shape_derivative(
    
const Collisionscollisionsconst CollisionMeshmesh,
    
const Eigen::MatrixXdvertices) const;

Compute the shape derivative of the potential.

Parameters:
const Collisions &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 Collisioncollision,
    
const VectorMax12dpositions) const override;

Compute the potential for a single collision.

Parameters:
const Collision &collision

The collision.

const VectorMax12d &positions

The collision stencil’s positions.

Returns:

The potential.

VectorMax12d gradient(const Collisioncollision,
    
const VectorMax12dpositions) const override;

Compute the gradient of the potential for a single collision.

Parameters:
const Collision &collision

The collision.

const VectorMax12d &positions

The collision stencil’s positions.

Returns:

The gradient of the potential.

MatrixMax12d hessian(const Collisioncollision,
    
const VectorMax12dpositions,
    
const PSDProjectionMethod project_hessian_to_psd
    
= PSDProjectionMethod::NONE
) const override;

Compute the hessian of the potential for a single collision.

Parameters:
const Collision &collision

The collision.

const VectorMax12d &positions

The collision stencil’s positions.

Returns:

The hessian of the potential.

void shape_derivative(const Collisioncollision,
    
const std::array<long, 4>& vertex_ids,
    
const VectorMax12drest_positions,
    
const VectorMax12dpositions,
    
std::vector<Eigen::Triplet<double>>& out) const;

Compute the shape derivative of the potential for a single collision.

Parameters:
const Collision &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.

Protected Functions

virtual double distance_based_potential(
    
const double distance_sqrconst double dmin = 0) const
   
 = 0;

Compute the unmollified distance-based potential for a collisions.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The unmollified distance-based potential.

virtual double distance_based_potential_gradient(
    
const double distance_sqrconst double dmin = 0) const
   
 = 0;

Compute the gradient of the unmollified distance-based potential for a collision.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The gradient of the unmollified distance-based potential.

virtual double distance_based_potential_hessian(
    
const double distance_sqrconst double dmin = 0) const
   
 = 0;

Compute the hessian of the unmollified distance-based potential for a collision.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The hessian of the unmollified distance-based potential.

Private Types

using Super = Potential<Collisions>;

Barrier Potential

class BarrierPotential : public ipc::DistanceBasedPotential;

Inheritance diagram for ipc::BarrierPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::Potential< Collisions >" tooltip="ipc::Potential< Collisions >"] "1" [label="ipc::BarrierPotential" tooltip="ipc::BarrierPotential" fillcolor="#BFBFBF"] "2" [label="ipc::DistanceBasedPotential" tooltip="ipc::DistanceBasedPotential"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::BarrierPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="ipc::Potential< Collisions >" tooltip="ipc::Potential< Collisions >"] "1" [label="ipc::BarrierPotential" tooltip="ipc::BarrierPotential" fillcolor="#BFBFBF"] "2" [label="ipc::DistanceBasedPotential" tooltip="ipc::DistanceBasedPotential"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] }

The barrier collision potential.

Public Functions

explicit BarrierPotential(
    
const double dhatconst bool use_physical_barrier = false);

Construct a barrier potential.

Parameters:
const double dhat

The activation distance of the barrier.

BarrierPotential(const std::shared_ptr<Barrier> barrier,
    
const double dhatconst bool use_physical_barrier = false);

Construct a barrier potential.

Parameters:
const double dhat

The activation distance of the barrier.

inline double dhat() const;

Get the activation distance of the barrier.

inline void set_dhat(const double dhat);

Set the activation distance of the barrier.

Parameters:
const double dhat

The activation distance of the barrier.

inline const Barrierbarrier() const;

Get the barrier function used to compute the potential.

inline void set_barrier(const std::shared_ptr<Barrier> barrier);

Set the barrier function used to compute the potential.

Parameters:
const std::shared_ptr<Barrier> barrier

The barrier function used to compute the potential.

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.

Parameters:
bool use_physical_barrier

Whether to use the physical barrier.

Protected Functions

virtual double distance_based_potential(const double distance_sqr,
    
const double dmin = 0) const override;

Compute the barrier potential for a collision.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The barrier potential.

virtual double distance_based_potential_gradient(
    
const double distance_sqr,
    
const double dmin = 0) const override;

Compute the gradient of the barrier potential for a collision.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The gradient of the barrier potential.

virtual double distance_based_potential_hessian(
    
const double distance_sqr,
    
const double dmin = 0) const override;

Compute the hessian of the barrier potential for a collision.

Parameters:
const double distance_sqr

The distance (squared) between the two objects.

const double dmin = 0

The minimum distance (unsquared) between the two objects.

Returns:

The hessian of the barrier potential.

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.

Friction Potential

class FrictionPotential : public ipc::Potential<FrictionCollisions>;

Inheritance diagram for ipc::FrictionPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Potential< FrictionCollisions >" tooltip="ipc::Potential< FrictionCollisions >"] "1" [label="ipc::FrictionPotential" tooltip="ipc::FrictionPotential" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for ipc::FrictionPotential:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="ipc::Potential< FrictionCollisions >" tooltip="ipc::Potential< FrictionCollisions >"] "1" [label="ipc::FrictionPotential" tooltip="ipc::FrictionPotential" fillcolor="#BFBFBF"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

The friction dissipative potential.

Public Types

enum class DiffWRT;

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

explicit FrictionPotential(const double epsv);

Construct a friction potential.

Parameters:
const double epsv

The smooth friction mollifier parameter \(\epsilon_v\).

inline double epsv() const;

Get the smooth friction mollifier parameter \(\epsilon_v\).

inline void set_epsv(const double epsv);

Set the smooth friction mollifier parameter \(\epsilon_v\).

Parameters:
const double epsv

The smooth friction mollifier parameter \(\epsilon_v\).

Eigen::VectorXd force(const FrictionCollisionscollisions,
    
const CollisionMeshmesh,
    
const Eigen::MatrixXdrest_positions,
    
const Eigen::MatrixXdlagged_displacements,
    
const Eigen::MatrixXdvelocities,
    
const BarrierPotentialbarrier_potential,
    
const double barrier_stiffnessconst double dmin = 0,
    
const bool no_mu = false) const;

Compute the friction force from the given velocities.

Parameters:
const FrictionCollisions &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 BarrierPotential &barrier_potential

Barrier potential (used for normal force magnitude).

const double barrier_stiffness

Barrier 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 FrictionCollisionscollisionsconst CollisionMeshmesh,
    
const Eigen::MatrixXdrest_positions,
    
const Eigen::MatrixXdlagged_displacements,
    
const Eigen::MatrixXdvelocities,
    
const BarrierPotentialbarrier_potential,
    
const double barrier_stiffnessconst DiffWRT wrt,
    
const double dmin = 0) const;

Compute the Jacobian of the friction force wrt the velocities.

Parameters:
const FrictionCollisions &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 BarrierPotential &barrier_potential

Barrier potential (used for normal force magnitude).

const double barrier_stiffness

Barrier 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 FrictionCollisioncollision,
    
const VectorMax12dvelocities) const override;

Compute the potential for a single collision.

Parameters:
const FrictionCollision &collision

The collision

const VectorMax12d &velocities

The collision stencil’s velocities.

Returns:

The potential.

VectorMax12d gradient(const FrictionCollisioncollision,
    
const VectorMax12dvelocities) const override;

Compute the gradient of the potential for a single collision.

Parameters:
const FrictionCollision &collision

The collision

const VectorMax12d &velocities

The collision stencil’s velocities.

Returns:

The gradient of the potential.

MatrixMax12d hessian(const FrictionCollisioncollision,
    
const VectorMax12dvelocities,
    
const PSDProjectionMethod project_hessian_to_psd
    
= PSDProjectionMethod::NONE
) const override;

Compute the hessian of the potential for a single collision.

Parameters:
const FrictionCollision &collision

The collision

const VectorMax12d &velocities

The collision stencil’s velocities.

Returns:

The hessian of the potential.

VectorMax12d force(const FrictionCollisioncollision,
    
const VectorMax12drest_positions,
    
const VectorMax12dlagged_displacements,
    
const VectorMax12dvelocities,
    
const BarrierPotentialbarrier_potential,
    
const double barrier_stiffnessconst double dmin = 0,
    
const bool no_mu = false) const;

Compute the friction force.

Parameters:
const FrictionCollision &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 BarrierPotential &barrier_potential

Barrier potential (used for normal force magnitude).

const double barrier_stiffness

Barrier 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 FrictionCollisioncollision,
    
const VectorMax12drest_positions,
    
const VectorMax12dlagged_displacements,
    
const VectorMax12dvelocities,
    
const BarrierPotentialbarrier_potential,
    
const double barrier_stiffnessconst DiffWRT wrt,
    
const double dmin = 0) const;

Compute the friction force Jacobian.

Parameters:
const FrictionCollision &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 BarrierPotential &barrier_potential

Barrier potential (used for normal force magnitude).

const double barrier_stiffness

Barrier 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 TCollisionscollisions,
    
const CollisionMeshmeshconst Eigen::MatrixXdX) 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|.

virtual VectorMax12d gradient(
    
const TCollisioncollisionconst VectorMax12dx) 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 TCollisionscollisions,
    
const CollisionMeshmeshconst Eigen::MatrixXdX,
    
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 MatrixMax12d hessian(const TCollisioncollision,
    
const VectorMax12dx,
    
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 Attributes

double m_epsv;

The smooth friction mollifier parameter \(\epsilon_v\).

Private Types

using Super = Potential;