Potentials

Generic Potential

template <class TCollisions> class Potential

Base class for potentials.

Template Parameters:
class TCollisions

The type of the collisions.

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:

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

Collaboration diagram for ipc::NormalPotential:

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

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.

Parameters:
const double distance_squared

The squared distance between elements.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The barrier stiffness.

Returns:

The force magnitude.

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.

Parameters:
const double distance_squared

The squared distance between elements.

const VectorMax12d &distance_squared_gradient

The gradient of the squared distance.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The stiffness of the barrier.

Returns:

The gradient of the force.

Protected Functions

virtual double operator()(
    const
 double distance_squared
, const double dmin = 0)
 const = 0

Compute the unmollified distance-based potential for a collisions.

Parameters:
const double distance_squared

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 gradient(
    const
 double distance_squared
, const double dmin = 0)
 const = 0

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

Parameters:
const double distance_squared

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 hessian(
    const
 double distance_squared
, const double dmin = 0)
 const = 0

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

Parameters:
const double distance_squared

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<NormalCollisions>

Barrier Potential

class BarrierPotential : public ipc::NormalPotential

Inheritence 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< NormalCollisions >" tooltip="ipc::Potential< NormalCollisions >"] "1" [label="ipc::BarrierPotential" tooltip="ipc::BarrierPotential" fillcolor="#BFBFBF"] "2" [label="ipc::NormalPotential" tooltip="ipc::NormalPotential"] "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< NormalCollisions >" tooltip="ipc::Potential< NormalCollisions >"] "1" [label="ipc::BarrierPotential" tooltip="ipc::BarrierPotential" fillcolor="#BFBFBF"] "2" [label="ipc::NormalPotential" tooltip="ipc::NormalPotential"] "1" -> "2" [dir=forward tooltip="public-inheritance"] "2" -> "3" [dir=forward tooltip="public-inheritance"] }

The barrier collision potential.

Public Functions

explicit BarrierPotential(
    const
 double dhat
, const 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 dhat
, const bool use_physical_barrier = false)

Construct a barrier potential.

Parameters:
const std::shared_ptr<Barrier> barrier

The barrier function.

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 Barrier& barrier() 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.

virtual double force_magnitude(const double distance_squared,
    const
 double dmin
,
    const
 double barrier_stiffness
)
 const override

Compute the force magnitude for a collision.

Parameters:
const double distance_squared

The squared distance between elements.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The barrier stiffness.

Returns:

The force magnitude.

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.

Parameters:
const double distance_squared

The squared distance between elements.

const VectorMax12d &distance_squared_gradient

The gradient of the squared distance.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The stiffness of the barrier.

Returns:

The gradient of the force.

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.

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.

Parameters:
const double distance_squared

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.

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.

Parameters:
const double distance_squared

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.

Protected Functions

virtual double operator()(const double distance_squared,
    const
 double dmin = 0
)
 const override

Compute the barrier potential for a collision.

Parameters:
const double distance_squared

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 gradient(const double distance_squared,
    const
 double dmin = 0
)
 const override

Compute the gradient of the barrier potential for a collision.

Parameters:
const double distance_squared

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 hessian(const double distance_squared,
    const
 double dmin = 0
)
 const override

Compute the hessian of the barrier potential for a collision.

Parameters:
const double distance_squared

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.

Private Types

using Super = NormalPotential

Normal Adhesion Potential

class NormalAdhesionPotential : public ipc::NormalPotential

Inheritence diagram for ipc::NormalAdhesionPotential:

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

Collaboration diagram for ipc::NormalAdhesionPotential:

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

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.

Parameters:
const double distance_squared

The squared distance between elements.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The barrier stiffness.

Returns:

The force magnitude.

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.

Parameters:
const double distance_squared

The squared distance between elements.

const VectorMax12d &distance_squared_gradient

The gradient of the squared distance.

const double dmin

The minimum distance offset to the barrier.

const double barrier_stiffness

The stiffness of the barrier.

Returns:

The gradient of the force.

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.

Parameters:
const double distance_squared

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.

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.

Parameters:
const double distance_squared

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.

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.

Parameters:
const double distance_squared

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 gradient(const double distance_squared,
    const
 double dmin = 0
)
 const override

Compute the gradient of the barrier potential for a collision.

Parameters:
const double distance_squared

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 hessian(const double distance_squared,
    const
 double dmin = 0
)
 const override

Compute the hessian of the barrier potential for a collision.

Parameters:
const double distance_squared

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.

Private Types

using Super = NormalPotential

Private Functions

std::array<double, 3> normal_adhesion_potential_args(
    const
 double dmin
)
 const

Tangential Potentials

class TangentialPotential
    
: public ipc::Potential<TangentialCollisions>

Inheritence diagram for ipc::TangentialPotential:

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

Collaboration diagram for ipc::TangentialPotential:

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

A tangential dissipative potential.

Subclassed by ipc::FrictionPotential, ipc::TangentialAdhesionPotential

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

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

virtual double f0(const double x) const = 0
virtual double f1_over_x(const double x) const = 0
virtual double f2_x_minus_f1_over_x3(const double x) const = 0
virtual bool is_dynamic(const double speed) const = 0

Private Types

using Super = Potential

Friction Potential

class FrictionPotential : public ipc::TangentialPotential

Inheritence diagram for ipc::FrictionPotential:

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

The friction dissipative potential.

Public Functions

explicit FrictionPotential(const double eps_v)

Construct a friction potential.

Parameters:
const double eps_v

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

inline double eps_v() const

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

inline void set_eps_v(const double eps_v)

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

Parameters:
const double eps_v

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

Protected Functions

virtual double f0(const double x) const override
virtual double f1_over_x(const double x) const override
virtual double f2_x_minus_f1_over_x3(const double x) const override
inline virtual bool is_dynamic(const double speed) const override

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:

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

Collaboration diagram for ipc::TangentialAdhesionPotential:

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

The tangential adhesion potential.

Public Functions

explicit TangentialAdhesionPotential(const double eps_a)

Construct a tangential adhesion potential.

Parameters:
const double eps_a

The tangential adhesion mollifier parameter \(\epsilon_a\).

inline double eps_a() const

Get the tangential adhesion mollifier parameter \(\epsilon_a\).

inline void set_eps_a(const double eps_a)

Set the tangential adhesion mollifier parameter \(\epsilon_v\).

Parameters:
const double eps_a

The tangential adhesion mollifier parameter \(\epsilon_v\).

Protected Functions

virtual double f0(const double x) const override
virtual double f1_over_x(const double x) const override
virtual double f2_x_minus_f1_over_x3(const double x) const override
inline virtual bool is_dynamic(const double speed) const override

Protected Attributes

double m_eps_a

The tangential adhesion mollifier parameter \(\epsilon_a\).

Private Types

using Super = TangentialPotential

Last update: Nov 16, 2024