NormalPotential¶
- class NormalPotential : public ipc::Potential<NormalCollisions>;¶
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,
Eigen::ConstRef<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.
- Eigen::ConstRef<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.
-
Eigen::VectorXd gauss_newton_hessian_diagonal(
const NormalCollisions& collisions, const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices) const;¶ Compute the diagonal of the cumulative Gauss-Newton Hessian of the potential.
Uses the distance-vector formulation to efficiently compute the diagonal of the Gauss-Newton Hessian without forming full local 12×12 Hessian matrices. This is useful as a Jacobi preconditioner for iterative solvers.
Note
This is a Gauss-Newton approximation (drops derivatives of closest-point coefficients), not the exact Hessian.
- Parameters:¶
- const NormalCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Vertices of the collision mesh.
- Returns:¶
The diagonal of the Gauss-Newton Hessian as a vector of size vertices.size().
-
double gauss_newton_hessian_quadratic_form(
const NormalCollisions& collisions, const CollisionMesh& mesh,
Eigen::ConstRef<Eigen::MatrixXd> vertices,
Eigen::ConstRef<Eigen::VectorXd> p) const;¶ Compute the product pᵀHp for the cumulative Gauss-Newton Hessian of the potential.
Uses the distance-vector formulation to efficiently compute the quadratic form pᵀHp without forming full local 12×12 Hessian matrices nor the global sparse Hessian. This is useful for nonlinear conjugate gradient methods.
Note
This is a Gauss-Newton approximation (drops derivatives of closest-point coefficients), not the exact Hessian.
- Parameters:¶
- const NormalCollisions &collisions¶
The set of collisions.
- const CollisionMesh &mesh¶
The collision mesh.
- Eigen::ConstRef<Eigen::MatrixXd> vertices¶
Vertices of the collision mesh.
- Eigen::ConstRef<Eigen::VectorXd> p¶
The direction vector of size vertices.size().
- Returns:¶
The scalar value pᵀHp.
-
double operator()(const NormalCollision& collision,
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> positions¶
The collision stencil’s positions.
- Returns:¶
The potential.
-
VectorMax12d gradient(const NormalCollision& collision,
Eigen::ConstRef<VectorMax12d> positions) const override;¶ Compute the gradient of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> positions¶
The collision stencil’s positions.
- Returns:¶
The gradient of the potential.
-
MatrixMax12d hessian(const NormalCollision& collision,
Eigen::ConstRef<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.
- Eigen::ConstRef<VectorMax12d> positions¶
The collision stencil’s positions.
- 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 gauss_newton_hessian_diagonal(
const NormalCollision& collision,
Eigen::ConstRef<VectorMax12d> positions) const;¶ Compute the diagonal of the Gauss-Newton Hessian of the potential for a single collision using the distance-vector formulation (Eqs.
10-12).
diag(∂²b/∂x²) ≈ w · [4·f”(d²) · diag((∂dv/∂x·dv)(∂dv/∂x·dv)ᵀ)
2·f’(d²) · diag((∂dv/∂x)(∂dv/∂x)ᵀ)]
Note
This is a Gauss-Newton approximation, not the exact Hessian diagonal. It drops derivatives of the closest-point coefficients.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> positions¶
The collision stencil’s positions.
- Returns:¶
The diagonal of the Gauss-Newton Hessian as a vector of size ndof.
-
double gauss_newton_hessian_quadratic_form(
const NormalCollision& collision,
Eigen::ConstRef<VectorMax12d> positions,
Eigen::ConstRef<VectorMax12d> p) const;¶ Compute pᵀHp for a single collision using the Gauss-Newton Hessian via the distance-vector formulation (Eqs.
10, 13-14).
pᵀHp ≈ w · [4·f”(d²) · (pᵀ·∂dv/∂x·dv)² + 2·f’(d²) · ‖pᵀ·∂dv/∂x‖²]
Note
This is a Gauss-Newton approximation, not the exact Hessian quadratic form. It drops derivatives of the closest-point coefficients.
- Parameters:¶
- const NormalCollision &collision¶
The collision.
- Eigen::ConstRef<VectorMax12d> positions¶
The collision stencil’s positions.
- Eigen::ConstRef<VectorMax12d> p¶
The local direction vector (size ndof, i.e., the stencil’s DOF slice of the global p).
- Returns:¶
The scalar value pᵀHp (approximate).
-
void shape_derivative(const NormalCollision& collision,
const std::array<index_t, 4>& vertex_ids,
Eigen::ConstRef<VectorMax12d> rest_positions,
Eigen::ConstRef<VectorMax12d> positions,
std::vector<Eigen::Triplet<double>>& out,
const int n_total_verts = -1) const;¶ Compute the shape derivative of the potential for a single collision.
- Parameters:¶
- const NormalCollision &collision¶
[in] The collision.
- const std::array<index_t, 4> &vertex_ids¶
[in] The collision stencil’s vertex ids.
- Eigen::ConstRef<VectorMax12d> rest_positions¶
[in] The collision stencil’s rest positions.
- Eigen::ConstRef<VectorMax12d> positions¶
[in] The collision stencil’s positions.
- std::vector<Eigen::Triplet<double>> &out¶
[inout] Store the triplets of the shape derivative here.
- const int n_total_verts = -1¶
[in] The total number of vertices in the mesh, used for computing global indices in the triplets. See also
local_hessian_to_global_triplets.
-
virtual double force_magnitude(
const double distance_squared, const double dmin) const
= 0;¶ Compute the force magnitude for a collision.
-
virtual VectorMax12d force_magnitude_gradient(
const double distance_squared,
Eigen::ConstRef<VectorMax12d> distance_squared_gradient,
const double dmin) 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>;¶