Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SmoothContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "ContactForm.hpp"
4#include <ipc/smooth_contact/smooth_collisions.hpp>
5#include <ipc/smooth_contact/smooth_contact_potential.hpp>
6#include <cmath>
7
8namespace polyfem::solver
9{
11 {
12 public:
13 SmoothContactForm(const ipc::CollisionMesh &collision_mesh,
14 const double dhat,
15 const double avg_mass,
16 const double alpha_t,
17 const double alpha_n,
18 const bool use_adaptive_dhat,
19 const double min_distance_ratio,
21 const bool is_time_dependent,
22 const bool enable_shape_derivatives,
23 const ipc::BroadPhaseMethod broad_phase_method,
24 const double ccd_tolerance,
25 const int ccd_max_iterations);
26
27 virtual std::string name() const override { return "smooth-contact"; }
28
29 void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override;
30
31 void force_shape_derivative(const ipc::SmoothCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term) const;
32
36 void post_step(const polysolve::nonlinear::PostStepData &data) override;
37
38 bool using_adaptive_dhat() const { return use_adaptive_dhat; }
39 const ipc::ParameterType &get_params() const { return params; }
40
41 const ipc::SmoothCollisions &collision_set() const { return collision_set_; }
42
43 protected:
47 double value_unweighted(const Eigen::VectorXd &x) const override;
48
52 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
53
57 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
58
62 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
63
64 double barrier_support_size() const override { return dhat_; }
65
66 void update_collision_set(const Eigen::MatrixXd &displaced_surface) override;
67
68 private:
69 ipc::ParameterType params;
71
73 ipc::SmoothCollisions collision_set_;
74
76 ipc::SmoothContactPotential barrier_potential_;
77 };
78}
int x
Form representing the contact potential and forces.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
const double dhat_
Barrier activation distance.
void update_collision_set(const Eigen::MatrixXd &displaced_surface) override
Update the cached candidate set for the current solution.
const ipc::SmoothCollisions & collision_set() const
ipc::SmoothContactPotential barrier_potential_
Contact potential.
void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override
Update the barrier stiffness based on the current elasticity energy.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied per element.
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
double barrier_support_size() const override
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
virtual std::string name() const override
const ipc::ParameterType & get_params() const
void force_shape_derivative(const ipc::SmoothCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term) const
ipc::SmoothCollisions collision_set_
Cached constraint set for the current solution.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22