PolyFEM
Loading...
Searching...
No Matches
SmoothContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "ContactForm.hpp"
4
6#include <polysolve/nonlinear/PostStepData.hpp>
7
8#include <ipc/smooth_contact/smooth_collisions.hpp>
9#include <ipc/smooth_contact/smooth_contact_potential.hpp>
10
11namespace polyfem::solver
12{
14 {
16
17 public:
18 SmoothContactForm(const ipc::CollisionMesh &collision_mesh,
19 const double dhat,
20 const double avg_mass,
21 const double alpha_t,
22 const double alpha_n,
23 const bool use_adaptive_dhat,
24 const double min_distance_ratio,
26 const bool is_time_dependent,
27 const bool enable_shape_derivatives,
28 const ipc::BroadPhaseMethod broad_phase_method,
29 const double ccd_tolerance,
30 const int ccd_max_iterations);
31
32 virtual std::string name() const override { return "smooth-contact"; }
33
34 void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override;
35
39 void post_step(const polysolve::nonlinear::PostStepData &data) override;
40
41 bool using_adaptive_dhat() const { return use_adaptive_dhat; }
42 const ipc::SmoothContactParameters &get_params() const { return params; }
43
44 const ipc::SmoothCollisions &collision_set() const { return collision_set_; }
45
46 protected:
50 double value_unweighted(const Eigen::VectorXd &x) const override;
51
55 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
56
60 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
61
65 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
66
67 double barrier_support_size() const override { return dhat_; }
68
69 void update_collision_set(const Eigen::MatrixXd &displaced_surface) override;
70
71 private:
72 ipc::SmoothContactParameters params;
74
76 ipc::SmoothCollisions collision_set_;
77
79 ipc::SmoothContactPotential barrier_potential_;
80 };
81} // namespace polyfem::solver
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.
const ipc::SmoothContactParameters & get_params() const
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
ipc::SmoothContactParameters params
ipc::SmoothCollisions collision_set_
Cached constraint set for the current solution.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24